SCALE-RM
mod_realinput_netcdf.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_file_h
19  use scale_tracer
21  use scale_hash
22 
23  use scale_prc, only: &
24  myrank => prc_myrank, &
25  prc_abort
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
33  public :: parentatmossetupnetcdf
34  public :: parentatmosopennetcdf
36  public :: parentatmosinputnetcdf
37  public :: parentlandsetupnetcdf
38  public :: parentlandopennetcdf
39  public :: parentlandfinalizenetcdf
40  public :: parentlandinputnetcdf
41  public :: parentoceansetupnetcdf
42  public :: parentoceanopennetcdf
44  public :: parentoceaninputnetcdf
45 
46  !-----------------------------------------------------------------------------
47  !
48  !++ Public parameters & variables
49  !
50  !-----------------------------------------------------------------------------
51  !
52  !++ Private procedure
53  !
54  !-----------------------------------------------------------------------------
55  !
56  !++ Private parameters & variables
57  !
58  type :: vinfo
59  character(len=32) :: name = ""
60  logical :: zstg = .false.
61  logical :: xstg = .false.
62  logical :: ystg = .false.
63  real(RP) :: fact = 1.0_rp
64  real(RP) :: offset = 0.0_rp
65  end type vinfo
66 
67  type(hash_table) :: vars_atmos
68  type(hash_table) :: vars_ocean
69  type(hash_table) :: vars_land
70 
71  real(RP), allocatable :: work3d(:,:,:)
72  real(RP), allocatable :: work2d(:,:)
73 
74 
75  logical :: SCALE_tile_atm
76  logical :: SCALE_tile_lnd
77  logical :: SCALE_tile_ocn
78  integer :: SCALE_DOMID_atm = -1
79  integer :: SCALE_DOMID_lnd = -1
80  integer :: SCALE_DOMID_ocn = -1
81  integer :: nfiles_atm = 0
82  integer :: nfiles_lnd = 0
83  integer :: nfiles_ocn = 0
84  integer :: fid_atm = -1
85  integer :: fid_lnd = -1
86  integer :: fid_ocn = -1
87  integer, allocatable :: fids_atm(:)
88  integer, allocatable :: fids_lnd(:)
89  integer, allocatable :: fids_ocn(:)
90  integer, allocatable :: tile_id_atm(:)
91  integer, allocatable :: tile_id_lnd(:)
92  integer, allocatable :: tile_id_ocn(:)
93 
94  integer, parameter :: vars_max = 100
95 
96 
97  character(len=32) :: zname, zhname
98  character(len=32) :: xname, xhname
99  character(len=32) :: yname, yhname
100  character(len=32) :: tname
101 
102  namelist / netcdf_dims / &
103  zname, &
104  zhname, &
105  xname, &
106  xhname, &
107  yname, &
108  yhname, &
109  tname
110 
111  character(len=32) :: item
112  character(len=32) :: name
113  logical :: zstg, xstg, ystg
114  real(RP) :: fact, offset
115 
116  namelist / netcdf_item / &
117  item, &
118  name, &
119  zstg, &
120  xstg, &
121  ystg, &
122  fact, &
123  offset
124 
125 
126  logical :: first_atm
127  logical :: first_ocn
128  logical :: first_lnd
129 
130 
131  !-----------------------------------------------------------------------------
132 contains
133  !-----------------------------------------------------------------------------
135  subroutine parentatmossetupnetcdf( &
136  dims, &
137  timelen, &
138  mixing_ratio, &
139  update_coord, &
140  mapping_info, &
141  qtrc_flag, &
142  lon_all, &
143  lat_all, &
144  basename_org, &
145  basename_num, &
146  same_mp_type, &
147  pt_dry, &
148  serial, &
149  do_read )
150  use scale_prc, only: &
152  use scale_const, only: &
153  undef => const_undef, &
154  grav => const_grav, &
155  d2r => const_d2r
156  use scale_file, only: &
158  file_get_attribute, &
159  file_get_datainfo
160  use scale_mapprojection, only: &
162  use scale_comm_cartesc, only: &
163  comm_bcast
164  use scale_comm_cartesc_nest, only: &
167  use mod_atmos_phy_mp_vars, only: &
168  qs_mp, &
169  qe_mp
170  implicit none
171  integer, intent(out) :: dims(6)
172  integer, intent(out) :: timelen
173  logical, intent(out) :: mixing_ratio
174  logical, intent(out) :: update_coord
175  type(mappinginfo), intent(out) :: mapping_info
176  logical, intent(out) :: qtrc_flag(qa)
177  real(rp), allocatable, intent(out) :: lon_all(:,:)
178  real(rp), allocatable, intent(out) :: lat_all(:,:)
179 
180  character(len=*), intent(in) :: basename_org
181  character(len=*), intent(in) :: basename_num
182  logical, intent(in) :: same_mp_type
183 
184  logical, intent(inout) :: pt_dry
185  logical, intent(inout) :: serial
186  logical, intent(inout) :: do_read
187 
188  character(len=8) :: file_type ! "SCALE-RM", "WRFARW", "NAMELIST", "AUTO"
189  character(len=FILE_HLONG) :: nm_file ! namelist
190  integer :: scale_parent_prc_num_x
191  integer :: scale_parent_prc_num_y
192  character(len=FILE_HLONG) :: scale_latlon_catalogue
193 
194  namelist / param_mkinit_real_atmos_netcdf / &
195  file_type, &
196  nm_file, &
197  mixing_ratio, &
198  scale_parent_prc_num_x, &
199  scale_parent_prc_num_y, &
200  scale_latlon_catalogue
201 
202  character(len=H_SHORT) :: mapping_name
203  real(dp) :: false_easting
204  real(dp) :: false_northing
205  real(dp) :: longitude_of_central_meridian
206  real(dp) :: longitude_of_projection_origin
207  real(dp) :: latitude_of_projection_origin
208  real(dp) :: straight_vertical_longitude_from_pole
209  real(dp) :: standard_parallel(2)
210  real(dp) :: rotation
211 
212  namelist / netcdf_mapprojection / &
213  mapping_name, &
214  false_easting, &
215  false_northing, &
216  longitude_of_central_meridian, &
217  longitude_of_projection_origin, &
218  latitude_of_projection_origin, &
219  straight_vertical_longitude_from_pole, &
220  standard_parallel, &
221  rotation
222 
223  character(len=32) :: items(vars_max)
224  integer :: nvars
225  type(vinfo), pointer :: var_info
226  class(*), pointer :: v
227 
228  character(len=FILE_HLONG) :: basename
229  character(len=FILE_HLONG) :: fname
230  character(len=32) :: map
231 
232  integer :: nmfid
233  integer :: i, n, iq
234  integer :: ierr
235  logical :: exist, error
236  !---------------------------------------------------------------------------
237 
238  log_newline
239  log_info("ParentAtmosSetupNetCDF",*) 'Real Case/Atmos Setup'
240 
241  file_type = "AUTO"
242  nm_file = ""
243  mixing_ratio = .false.
244  scale_parent_prc_num_x = -1
245  scale_parent_prc_num_y = -1
246  scale_latlon_catalogue = ""
247 
248  !--- read namelist
249  rewind(io_fid_conf)
250  read(io_fid_conf,nml=param_mkinit_real_atmos_netcdf,iostat=ierr)
251  if( ierr > 0 ) then
252  log_error("ParentAtmosSetupNetCDF",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_ATMOS_NetCDF. Check!'
253  call prc_abort
254  endif
255  log_nml(param_mkinit_real_atmos_netcdf)
256 
257  basename = trim(basename_org) // trim(basename_num)
258 
259  fid_atm = -1
260  if ( do_read ) then
261  call check_filetype(fid_atm, file_type, basename, scale_tile_atm, "ParentAtmosOpenNetCDF")
262  end if
263 
264  call comm_bcast( file_type )
265 
266  if ( file_type == "SCALE-RM" ) then
267  call comm_bcast( scale_tile_atm )
268  if ( scale_tile_atm ) then
269  do_read = .true.
270  serial = .false.
271  end if
272  end if
273 
274 
275  if ( do_read ) then
276 
277  mapping_name = ""
278  false_easting = undef
279  false_northing = undef
280  longitude_of_central_meridian = undef
281  longitude_of_projection_origin = undef
282  latitude_of_projection_origin = undef
283  straight_vertical_longitude_from_pole = undef
284  standard_parallel = (/ undef, undef /)
285  rotation = undef
286 
287  vars_atmos = hash_table()
288 
289  select case( file_type )
290  case ( "SCALE-RM" )
291  zname = "z"
292  zhname = "zh"
293  xname = "x"
294  xhname = "xh"
295  yname = "y"
296  yhname = "yh"
297  tname = "time"
298 
299  call vars_atmos%put("lon", vinfo("lon"))
300  call vars_atmos%put("lat", vinfo("lat"))
301 
302  call vars_atmos%put("height", vinfo("height"))
303  call vars_atmos%put("pressure", vinfo("PRES"))
304 
305  call vars_atmos%put("DENS", vinfo("DENS"))
306  call vars_atmos%put("W", vinfo("W"))
307  call vars_atmos%put("MOMZ", vinfo("MOMZ", zstg=.true.))
308  call vars_atmos%put("Umet", vinfo("Umet"))
309  call vars_atmos%put("U", vinfo("U"))
310  call vars_atmos%put("MOMX", vinfo("MOMX", xstg=.true.))
311  call vars_atmos%put("Vmet", vinfo("Vmet"))
312  call vars_atmos%put("V", vinfo("V"))
313  call vars_atmos%put("MOMY", vinfo("MOMY", ystg=.true.))
314 
315  call vars_atmos%put("PT", vinfo("PT"))
316  call vars_atmos%put("T", vinfo("T"))
317  call vars_atmos%put("RHOT", vinfo("RHOT"))
318 
319  call vars_atmos%put("QV", vinfo("QV"))
320  call vars_atmos%put("RH", vinfo("RH"))
321 
322  if ( same_mp_type ) then
323  do iq = qs_mp, qe_mp
324  call vars_atmos%put(tracer_name(iq), vinfo(tracer_name(iq)))
325  end do
326  else
327  call vars_atmos%put("QC", vinfo("QC"))
328  call vars_atmos%put("QR", vinfo("QR"))
329  call vars_atmos%put("QI", vinfo("QI"))
330  call vars_atmos%put("QS", vinfo("QS"))
331  call vars_atmos%put("QG", vinfo("QG"))
332 
333  call vars_atmos%put("NC", vinfo("NC"))
334  call vars_atmos%put("NR", vinfo("NR"))
335  call vars_atmos%put("NI", vinfo("NI"))
336  call vars_atmos%put("NS", vinfo("NS"))
337  call vars_atmos%put("NG", vinfo("NG"))
338  end if
339 
340  call vars_atmos%put("topo", vinfo("topo"))
341  call vars_atmos%put("MSLP", vinfo("MSLP"))
342  call vars_atmos%put("SFC_PRES", vinfo("SFC_PRES"))
343  call vars_atmos%put("U10met", vinfo("U10met"))
344  call vars_atmos%put("U10", vinfo("U10"))
345  call vars_atmos%put("V10met", vinfo("V10met"))
346  call vars_atmos%put("V10", vinfo("V10"))
347  call vars_atmos%put("T2", vinfo("T2"))
348  call vars_atmos%put("Q2", vinfo("Q2"))
349  call vars_atmos%put("RH2", vinfo("RH2"))
350 
351  mixing_ratio = .false.
352  update_coord = .false.
353  pt_dry = .false.
354 
355  if ( prc_ismaster ) then
356  call file_get_attribute( fid_atm, "QV", "grid_mapping", map, existed=exist )
357  if ( exist ) then
358  call file_get_attribute( fid_atm, map, "grid_mapping_name", mapping_name )
359 
360  call file_get_attribute( fid_atm, map, "false_easting", false_easting, existed=exist )
361  call file_get_attribute( fid_atm, map, "false_northing", false_northing, existed=exist )
362  call file_get_attribute( fid_atm, map, "longitude_of_central_meridian", longitude_of_central_meridian, existed=exist )
363  call file_get_attribute( fid_atm, map, "longitude_of_projection_origin", longitude_of_projection_origin, existed=exist )
364  call file_get_attribute( fid_atm, map, "latitude_of_projection_origin", latitude_of_projection_origin, existed=exist )
365  call file_get_attribute( fid_atm, map, "straight_vertical_longitude_from_pole", straight_vertical_longitude_from_pole, existed=exist )
366  call file_get_attribute( fid_atm, map, "standard_parallel", standard_parallel(1:1), existed=exist )
367  if ( .not. exist ) &
368  call file_get_attribute( fid_atm, map, "standard_parallel", standard_parallel(:), existed=exist )
369  call file_get_attribute( fid_atm, map, "rotation", rotation, existed=exist )
370  end if
371  end if
372 
373  call comm_bcast( mapping_name )
374 
375  call comm_bcast( false_easting )
376  call comm_bcast( false_northing )
377  call comm_bcast( longitude_of_central_meridian )
378  call comm_bcast( longitude_of_projection_origin )
379  call comm_bcast( latitude_of_projection_origin )
380  call comm_bcast( straight_vertical_longitude_from_pole )
381  call comm_bcast( 2, standard_parallel )
382  call comm_bcast( rotation )
383 
384  case ( "WRFARW" )
385  zname = "bottom_top"
386  zhname = "bottom_top_stag"
387  xname = "west_east"
388  xhname = "west_east_stag"
389  yname = "south_north"
390  yhname = "south_north_stag"
391  tname = "Time"
392 
393  call vars_atmos%put("lon", vinfo("XLONG"))
394  call vars_atmos%put("lat", vinfo("XLAT"))
395 
396  call vars_atmos%put("hbar", vinfo("PHB", zstg=.true., fact=1.0_rp/grav)) ! geopotential height
397  call vars_atmos%put("hdev", vinfo("PH", zstg=.true., fact=1.0_rp/grav))
398 
399  call vars_atmos%put("pbar", vinfo("PB"))
400  call vars_atmos%put("pdev", vinfo("P"))
401 
402  call file_get_datainfo( fid_atm, "U", existed=exist )
403  if ( exist ) then
404  log_info("ParentAtmosSetupNetCDF",*) 'WRF-ARW FILE-TYPE: WRF History Output'
405  call vars_atmos%put("W", vinfo("W",zstg=.true.))
406  call vars_atmos%put("U", vinfo("U",xstg=.true.))
407  call vars_atmos%put("V", vinfo("V",ystg=.true.))
408  call vars_atmos%put("PT", vinfo("T", offset=300.0_rp))
409  else
410  log_info("ParentAtmosSetupNetCDF",*) 'WRF-ARW FILE-TYPE: WRF Restart'
411  call vars_atmos%put("W", vinfo("W_1"))
412  call vars_atmos%put("U", vinfo("U_1"))
413  call vars_atmos%put("V", vinfo("V_1"))
414  call vars_atmos%put("PT", vinfo("T_1", offset=300.0_rp))
415  endif
416 
417  if ( same_mp_type ) then
418  log_error("ParentAtmosSetupNetCDF",*) 'same_mp_type must be .false. for WRF file'
419  call prc_abort
420  end if
421  call vars_atmos%put("QV", vinfo("QVAPOR"))
422  call vars_atmos%put("QC", vinfo("QCLOUD"))
423  call vars_atmos%put("QR", vinfo("QRAIN"))
424  call vars_atmos%put("QI", vinfo("QICE"))
425  call vars_atmos%put("QS", vinfo("QSNOW"))
426  call vars_atmos%put("QG", vinfo("QGRAUP"))
427  call vars_atmos%put("NC", vinfo("NC"))
428  call vars_atmos%put("NR", vinfo("NR"))
429  call vars_atmos%put("NI", vinfo("NI"))
430  call vars_atmos%put("NS", vinfo("NS"))
431  call vars_atmos%put("NG", vinfo("NG"))
432 
433  mixing_ratio = .true.
434  pt_dry = .true.
435 
436  call vars_atmos%put("topo", vinfo("HGT"))
437  call vars_atmos%put("U10", vinfo("U10"))
438  call vars_atmos%put("V10", vinfo("V10"))
439  call vars_atmos%put("T2", vinfo("T2"))
440  call vars_atmos%put("Q2", vinfo("Q2"))
441  call vars_atmos%put("RH2", vinfo("RH2"))
442  call vars_atmos%put("SFC_PRES", vinfo("PSFC"))
443 
444  call file_get_attribute( fid_atm, "global", "MAP_PROJ", i, existed=exist )
445  if ( exist ) then
446  if ( i == 1 ) then ! Lambert Conformal
447  mapping_name = "lambert_conformal_conic"
448  call file_get_attribute( fid_atm, "global", "TRUELAT1", standard_parallel(1) )
449  call file_get_attribute( fid_atm, "global", "TRUELAT2", standard_parallel(2) )
450  call file_get_attribute( fid_atm, "global", "STAND_LON", longitude_of_central_meridian )
451  else if ( i >= 3 ) then ! No rotate
452  ! do nothing
453  else
454  log_warn("ParentAtmodSetupNetCDF",*) "This map projection type is not supported: ", i
455  log_warn_cont(*) "Specify map projection parameter manually"
456  end if
457  end if
458 
459  update_coord = .true.
460 
461  case ( "NAMELIST" )
462 
463  update_coord = .true.
464 
465  case default
466  log_error("ParentAtmosSetupNetCDF",*) 'FILE_TYPE must be "SCALE-RM", "WRFARW", "NAMELIST", or "AUTO", ', trim(file_type)
467  call prc_abort
468  end select
469 
470  !--- read namelist
471  nmfid = -1
472  if ( nm_file /= "" ) then
473  nmfid = io_get_available_fid()
474  call io_get_fname(fname, nm_file)
475  open(nmfid, file=fname, form="formatted", status="old", action="read", iostat=ierr)
476  if ( ierr /= 0 ) then
477  log_error("ParentAtmosSetupNetCDF",*) 'namelist file is not found! ', trim(fname)
478  call prc_abort
479  end if
480 
481  read(nmfid, nml=netcdf_dims, iostat=ierr)
482  if( ierr > 0 ) then
483  log_error("ParentAtmosSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_DIMS in ', trim(fname), '. Check!'
484  call prc_abort
485  end if
486 
487  rewind(nmfid)
488  read(nmfid, nml=netcdf_mapprojection, iostat=ierr)
489  if( ierr > 0 ) then
490  log_error("ParentAtmosSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_MAPPROJECTION in ', trim(fname), '. Check!'
491  call prc_abort
492  end if
493  ! items
494  rewind(nmfid)
495  nvars = 0
496  do n = 1, vars_max
497  read(nmfid, nml=netcdf_item, iostat=ierr)
498  if( ierr > 0 ) then
499  log_error("ParentAtmosSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_ITEM in ', trim(fname), '. Check!'
500  call prc_abort
501  else if( ierr < 0 ) then
502  exit
503  end if
504  nvars = nvars + 1
505  items(nvars) = item
506  end do
507  if ( nvars > vars_max ) then
508  log_error("ParentAtmosSetupNetCDF",*) "The number of item in the namelist file exceeds the limit! ", nvars
509  call prc_abort
510  end if
511  rewind(nmfid)
512  do n = 1, nvars
513  ! set default
514  if ( vars_atmos%has_key(items(n)) ) then
515  item = items(n)
516  v => vars_atmos%get(item)
517  select type( v )
518  type is (vinfo)
519  var_info => v
520  end select
521  name = var_info%name
522  zstg = var_info%zstg
523  xstg = var_info%xstg
524  ystg = var_info%ystg
525  fact = var_info%fact
526  offset = var_info%offset
527  else
528  item = items(n)
529  name = items(n)
530  zstg = .false.
531  xstg = .false.
532  ystg = .false.
533  fact = 1.0_rp
534  offset = 0.0_rp
535  end if
536  read(nmfid, nml=netcdf_item, iostat=ierr)
537  if ( ierr /= 0 ) exit
538  ! set params
539  call vars_atmos%put(item, vinfo(name=name, zstg=zstg, xstg=xstg, ystg=ystg, fact=fact, offset=offset))
540  end do
541 
542  else if ( file_type == "NAMELIST" ) then
543  log_error("ParentAtmosSetupNetCDF",*) 'NM_FILE is necessary'
544  call prc_abort
545  end if
546 
547  mapping_info%mapping_name = mapping_name
548  if ( false_easting /= undef ) mapping_info%false_easting = false_easting
549  if ( false_northing /= undef ) mapping_info%false_northing = false_northing
550  if ( longitude_of_central_meridian /= undef ) mapping_info%longitude_of_central_meridian = longitude_of_central_meridian
551  if ( longitude_of_projection_origin /= undef ) mapping_info%longitude_of_projection_origin = longitude_of_projection_origin
552  if ( latitude_of_projection_origin /= undef ) mapping_info%latitude_of_projection_origin = latitude_of_projection_origin
553  if ( straight_vertical_longitude_from_pole /= undef ) mapping_info%straight_vertical_longitude_from_pole = straight_vertical_longitude_from_pole
554  if ( standard_parallel(1) /= undef ) mapping_info%standard_parallel(1) = standard_parallel(1)
555  if ( standard_parallel(2) /= undef ) mapping_info%standard_parallel(2) = standard_parallel(2)
556  if ( rotation /= undef ) mapping_info%rotation = rotation
557 
558  end if
559 
560  if ( scale_tile_atm ) then
561 
563  scale_domid_atm, & ! [OUT]
564  basename, & ! [IN]
565  scale_parent_prc_num_x, & ! [IN]
566  scale_parent_prc_num_y, & ! [IN]
567  scale_latlon_catalogue ) ! [IN]
568 
570  scale_domid_atm, & ! [IN]
571  kmax=dims(1), & ! [OUT]
572  imaxg=dims(2), & ! [OUT]
573  jmaxg=dims(3), & ! [OUT]
574  num_tile=nfiles_atm ) ! [OUT]
575 
576  dims(4) = dims(1)
577  dims(5) = dims(2)
578  dims(6) = dims(3)
579 
580  allocate( fids_atm(nfiles_atm) )
581  allocate( tile_id_atm(nfiles_atm) )
582 
584  scale_domid_atm, & ! [IN]
585  tile_id = tile_id_atm ) ! [OUT]
586 
587  call parentatmosopennetcdf( basename_org, basename_num )
588 
589  else if ( do_read ) then
590 
591  call file_get_dimlength( fid_atm, zname, dims(1) )
592  call file_get_dimlength( fid_atm, xname, dims(2) )
593  call file_get_dimlength( fid_atm, yname, dims(3) )
594  call file_get_dimlength( fid_atm, zhname, dims(4) )
595  call file_get_dimlength( fid_atm, xhname, dims(5) )
596  call file_get_dimlength( fid_atm, yhname, dims(6) )
597 
598  end if
599 
600  if ( do_read ) then
601 
602  do iq = 1, qa
603  if ( iq >= qs_mp .and. iq <= qe_mp ) cycle
604  qtrc_flag(iq) = .false.
605  if ( vars_atmos%has_key( tracer_name(iq) ) ) then
606  select type ( v => vars_atmos%get( tracer_name(iq) ) )
607  type is ( vinfo )
608  if ( v%name .ne. "" ) then
609  call file_get_datainfo( fid_atm, v%name, existed = qtrc_flag(iq) )
610  end if
611  end select
612  end if
613  end do
614 
615  call file_get_dimlength( fid_atm, tname, timelen, error=error )
616  if ( error ) timelen = 1
617 
618  allocate( lon_all(dims(2), dims(3)) )
619  allocate( lat_all(dims(2), dims(3)) )
620 
621  call read2d( dims(2), 1, dims(2), dims(3), 1, dims(3), &
622  lon_all(:,:), vars_atmos%get("lon"), &
623  1, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm )
624  lon_all(:,:) = lon_all(:,:) * d2r
625  call read2d( dims(2), 1, dims(2), dims(3), 1, dims(3), &
626  lat_all(:,:), vars_atmos%get("lat"), &
627  1, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm )
628  lat_all(:,:) = lat_all(:,:) * d2r
629 
630  end if
631 
632 
633  first_atm = .true.
634 
635  return
636  end subroutine parentatmossetupnetcdf
637 
638  !-----------------------------------------------------------------------------
640  subroutine parentatmosopennetcdf( &
641  basename_org, basename_num )
642  use scale_file, only: &
643  file_open
644  implicit none
645  character(len=*), intent(in) :: basename_org
646  character(len=*), intent(in) :: basename_num
647 
648  character(len=FILE_HLONG) :: basename
649 
650  integer :: n
651  !---------------------------------------------------------------------------
652 
653  log_newline
654  log_info("ParentAtmosOpenNetCDF",*) 'Real Case/Atmos Open'
655 
656 
657  basename = trim(basename_org) // trim(basename_num)
658 
659  if ( scale_tile_atm ) then
660 
661  do n = 1, nfiles_atm
662  call file_open( &
663  basename, & ! [IN]
664  fids_atm(n), & ! [OUT]
665  aggregate=.false., & ! [IN]
666  allnodes=.false., & ! [IN]
667  rankid=tile_id_atm(n) ) ! [IN]
668  end do
669 
670  fid_atm = fids_atm(1)
671 
672  else
673 
674  call file_open(basename, fid_atm, postfix="")
675 
676  end if
677 
678  return
679  end subroutine parentatmosopennetcdf
680 
681  !-----------------------------------------------------------------------------
683  subroutine parentatmosfinalizenetcdf
684  implicit none
685  !---------------------------------------------------------------------------
686 
687  log_newline
688  log_info("ParentAtmosFinalizeNetCDF",*) 'Real Case/Atmos Finalize'
689 
690  deallocate( work2d )
691  deallocate( work3d )
692  if ( allocated(fids_atm) ) deallocate( fids_atm )
693  if ( allocated(tile_id_atm) ) deallocate( tile_id_atm )
694 
695  call vars_atmos%destroy()
696 
697  scale_domid_atm = -1
698  nfiles_atm = 0
699  fid_atm = -1
700 
701  return
702  end subroutine parentatmosfinalizenetcdf
703 
704  !-----------------------------------------------------------------------------
705  subroutine parentatmosinputnetcdf( &
706  KA_org, KS_org, KE_org, &
707  IA_org, IS_org, IE_org, &
708  JA_org, JS_org, JE_org, &
709  QA, &
710  cz_org, &
711  w_org, u_org, v_org, &
712  pres_org, &
713  dens_org, &
714  temp_org, pt_org, &
715  qtrc_org, &
716  qv_org, rh_org, &
717  qhyd_org, qnum_org, &
718  nopres, nodens, &
719  uvmet, &
720  temp2pt, rh2qv, &
721  qnum_flag, &
722  same_mp_type, &
723  sfc_diagnoses, &
724  update_coord, &
725  dims, &
726  it )
727  use scale_const, only: &
728  undef => const_undef
729  use scale_atmos_hydrometeor, only: &
730  hyd_name, &
731  num_name, &
732  n_hyd
733  use mod_atmos_phy_mp_vars, only: &
734  qs_mp, &
735  qe_mp
736  implicit none
737  integer, intent(in) :: ka_org, ks_org, ke_org
738  integer, intent(in) :: ia_org, is_org, ie_org
739  integer, intent(in) :: ja_org, js_org, je_org
740  integer, intent(in) :: qa
741 
742  real(rp), intent(inout) :: cz_org(ka_org,ia_org,ja_org)
743 
744  real(rp), intent(out) :: w_org(ka_org,ia_org,ja_org)
745  real(rp), intent(out) :: u_org(ka_org,ia_org,ja_org)
746  real(rp), intent(out) :: v_org(ka_org,ia_org,ja_org)
747  real(rp), intent(out) :: pres_org(ka_org,ia_org,ja_org)
748  real(rp), intent(out) :: dens_org(ka_org,ia_org,ja_org)
749  real(rp), intent(out) :: temp_org(ka_org,ia_org,ja_org)
750  real(rp), intent(out) :: pt_org(ka_org,ia_org,ja_org)
751  real(rp), intent(out) :: qtrc_org(ka_org,ia_org,ja_org,qa)
752  real(rp), intent(out) :: qv_org(ka_org,ia_org,ja_org)
753  real(rp), intent(out) :: rh_org(ka_org,ia_org,ja_org)
754  real(rp), intent(out) :: qhyd_org(ka_org,ia_org,ja_org,n_hyd)
755  real(rp), intent(out) :: qnum_org(ka_org,ia_org,ja_org,n_hyd)
756  logical, intent(out) :: nopres
757  logical, intent(out) :: nodens
758  logical, intent(out) :: uvmet
759  logical, intent(out) :: temp2pt
760  logical, intent(out) :: rh2qv
761  logical, intent(out) :: qnum_flag
762 
763  logical, intent(in) :: same_mp_type
764  logical, intent(in) :: sfc_diagnoses
765  logical, intent(in) :: update_coord
766  integer, intent(in) :: dims(6)
767  integer, intent(in) :: it
768 
769 
770  logical :: exist
771  integer :: k, i, j, iq
772  !---------------------------------------------------------------------------
773 
774  if ( .not. allocated( work2d ) ) then
775  allocate( work2d(ia_org,ja_org) )
776  allocate( work3d(ka_org-2,ia_org,ja_org) )
777  end if
778 
779  ! height
780  if ( first_atm .or. update_coord ) then
781  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
782  cz_org(:,:,:), vars_atmos%get("height"), &
783  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
784  exist = exist )
785  if ( .not. exist ) then
786  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
787  cz_org(:,:,:), vars_atmos%get("hbar"), &
788  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
789  exist = exist )
790  if ( exist ) then
791  call read3d( ka_org-2, ks_org, ke_org-2, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
792  work3d(:,:,:), vars_atmos%get("hdev"), &
793  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
794  exist = exist )
795  end if
796  if ( .not. exist ) then
797  log_error("ParentAtmosInputNetCDF",*) '"height" or "hbar"+"hdev" is necessary'
798  call prc_abort
799  end if
800  !$omp parallel do collapse(2)
801  do j = 1, ja_org
802  do i = 1, ia_org
803  do k = 1, ka_org-2
804  cz_org(k+2,i,j) = cz_org(k+2,i,j) + work3d(k,i,j)
805  end do
806  end do
807  end do
808  else
809  end if
810  end if
811 
812  ! tracers other than water contents
813  do iq = 1, qa
814  if ( iq >= qs_mp .and. iq <= qe_mp ) cycle
815  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
816  qtrc_org(:,:,:,iq), vars_atmos%get(tracer_name(iq)), &
817  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
818  exist = exist )
819  if ( .not. exist ) then
820  !$omp parallel do collapse(2)
821  do j = 1, ja_org
822  do i = 1, ia_org
823  do k = 1, ka_org-2
824  qtrc_org(k+2,i,j,iq) = undef
825  end do
826  end do
827  end do
828  end if
829  end do
830 
831  ! water contents
832  if ( same_mp_type ) then
833 
834  do iq = qs_mp, qe_mp
835  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
836  qtrc_org(:,:,:,iq), vars_atmos%get(tracer_name(iq)), &
837  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
838  exist = exist )
839  if ( .not. exist ) then
840  !$omp parallel do collapse(2)
841  do j = 1, ja_org
842  do i = 1, ia_org
843  do k = 1, ka_org-2
844  qtrc_org(k+2,i,j,iq) = undef
845  end do
846  end do
847  end do
848  end if
849  end do
850 
851  else
852 
853  ! qv
854  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
855  qv_org(:,:,:), vars_atmos%get("QV"), &
856  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
857  exist = exist )
858  if ( exist ) then
859  rh2qv = .false.
860  else
861  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
862  rh_org(:,:,:), vars_atmos%get("RH"), &
863  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm )
864  rh2qv = .true.
865  end if
866 
867  ! qhyd
868  qnum_flag = .false.
869  do iq = 1, n_hyd
870  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
871  qhyd_org(:,:,:,iq), vars_atmos%get(hyd_name(iq)), &
872  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
873  exist = exist )
874  if ( .not. exist ) then
875  !$omp parallel do collapse(2)
876  do j = 1, ja_org
877  do i = 1, ia_org
878  do k = 1, ka_org-2
879  qhyd_org(k+2,i,j,iq) = 0.0_rp
880  end do
881  end do
882  end do
883  end if
884 
885  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
886  qnum_org(:,:,:,iq), vars_atmos%get(num_name(iq)), &
887  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
888  exist = exist )
889  if ( exist ) then
890  qnum_flag = .true.
891  else
892  !$omp parallel do collapse(2)
893  do j = 1, ja_org
894  do i = 1, ia_org
895  do k = 1, ka_org-2
896  qnum_org(k+2,i,j,iq) = undef
897  end do
898  end do
899  end do
900  end if
901  end do
902 
903  end if ! same_mp_type
904 
905  ! pressure
906  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
907  pres_org(:,:,:), vars_atmos%get("pressure"), &
908  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
909  exist = exist )
910  if ( .not. exist ) then
911  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
912  pres_org(:,:,:), vars_atmos%get("pbar"), &
913  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
914  exist = exist )
915  if ( exist ) then
916  call read3d( ka_org-2, ks_org, ke_org-2, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
917  work3d(:,:,:), vars_atmos%get("pdev"), &
918  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
919  exist = exist )
920  if ( .not. exist ) then
921  log_error("ParentAtmosInputNetCDF",*) '"pdev" is necessary if "pbar" exists'
922  call prc_abort
923  end if
924  !$omp parallel do collapse(2)
925  do j = 1, ja_org
926  do i = 1, ia_org
927  do k = 1, ka_org-2
928  pres_org(k+2,i,j) = pres_org(k+2,i,j) + work3d(k,i,j)
929  end do
930  end do
931  end do
932  end if
933  end if
934  if ( exist ) then
935  nopres = .false.
936  else
937  nopres = .true.
938  end if
939 
940  ! density
941  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
942  dens_org(:,:,:), vars_atmos%get("DENS"), &
943  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
944  exist = exist )
945  nodens = ( .not. exist )
946 
947  ! pt
948  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
949  pt_org(:,:,:), vars_atmos%get("PT"), &
950  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
951  exist = exist )
952  if ( exist ) then
953  temp2pt = .false.
954  else
955  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
956  pt_org(:,:,:), vars_atmos%get("RHOT"), &
957  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
958  exist = exist )
959  if ( exist ) then
960  if ( nodens ) then
961  if ( first_atm ) &
962  log_warn("ParentAtmosInputNetCDF",*) "RHOT is not used because DENS does not exist"
963  exist = .false.
964  else
965  !$omp parallel do collapse(2)
966  do j = 1, ja_org
967  do i = 1, ia_org
968  do k = 1, ka_org-2
969  pt_org(k+2,i,j) = pt_org(k+2,i,j) / dens_org(k+2,i,j)
970  end do
971  end do
972  end do
973  temp2pt = .false.
974  end if
975  end if
976  end if
977  if ( .not. exist ) then
978  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
979  temp_org(:,:,:), vars_atmos%get("T"), &
980  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
981  exist = exist )
982  temp2pt = .true.
983  end if
984  if ( .not. exist ) then
985  log_error("ParentAtmosInputNetCDF",*) 'Either "PT", "RHOT", or "T" is necessary'
986  call prc_abort
987  end if
988 
989  ! W
990  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
991  w_org(:,:,:), vars_atmos%get("W"), &
992  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
993  exist = exist )
994  if ( .not. exist ) then
995  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
996  w_org(:,:,:), vars_atmos%get("MOMZ"), &
997  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
998  exist = exist )
999  if ( exist ) then
1000  if ( nodens ) then
1001  if ( first_atm ) &
1002  log_warn("ParentAtmosInputNetCDF",*) "MOMZ is not used because DENS does not exist"
1003  exist = .false.
1004  else
1005  !$omp parallel do collapse(2)
1006  do j = 1, ja_org
1007  do i = 1, ia_org
1008  do k = 1, ka_org-2
1009  w_org(k+2,i,j) = w_org(k+2,i,j) / dens_org(k+2,i,j)
1010  end do
1011  end do
1012  end do
1013  end if
1014  end if
1015  end if
1016  if ( .not. exist ) then
1017  !$omp parallel do collapse(2)
1018  do j = 1, ja_org
1019  do i = 1, ia_org
1020  do k = 1, ka_org-2
1021  w_org(k+2,i,j) = 0.0_rp
1022  end do
1023  end do
1024  end do
1025  end if
1026 
1027  ! U
1028  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1029  u_org(:,:,:), vars_atmos%get("Umet"), &
1030  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1031  exist = exist )
1032  if ( exist ) then
1033  uvmet = .true.
1034  else
1035  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1036  u_org(:,:,:), vars_atmos%get("U"), &
1037  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1038  exist = exist )
1039  if ( exist ) then
1040  uvmet = .false.
1041  end if
1042  end if
1043  if ( .not. exist ) then
1044  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1045  u_org(:,:,:), vars_atmos%get("MOMX"), &
1046  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1047  exist = exist )
1048  if ( exist ) then
1049  if ( nodens ) then
1050  if ( first_atm ) &
1051  log_warn("ParentAtmosInputNetCDF",*) "MOMX is not used because DENS does not exist"
1052  exist = .false.
1053  else
1054  !$omp parallel do collapse(2)
1055  do j = 1, ja_org
1056  do i = 1, ia_org
1057  do k = 1, ka_org-2
1058  u_org(k+2,i,j) = u_org(k+2,i,j) / dens_org(k+2,i,j)
1059  end do
1060  end do
1061  end do
1062  end if
1063  uvmet = .false.
1064  end if
1065  end if
1066  if ( .not. exist ) then
1067  log_error("ParentAtmosInputNetCDF",*) 'Either "Ument", "U", or "MOMX" is necessary'
1068  call prc_abort
1069  end if
1070 
1071  ! V
1072  if ( uvmet ) then
1073  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1074  v_org(:,:,:), vars_atmos%get("Vmet"), &
1075  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1076  exist = exist )
1077  if ( .not. exist ) then
1078  log_error("ParentAtmosInputNetCDF",*) "Vmet is required when Umet exists"
1079  call prc_abort
1080  end if
1081  else
1082  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1083  v_org(:,:,:), vars_atmos%get("V"), &
1084  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1085  exist = exist )
1086  if ( .not. exist ) then
1087  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1088  v_org(:,:,:), vars_atmos%get("MOMY"), &
1089  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1090  exist = exist )
1091  if ( exist ) then
1092  if ( nodens ) then
1093  if ( first_atm ) &
1094  log_warn("ParentAtmosInputNetCDF",*) "MOMY is not used because DENS does not exist"
1095  exist = .false.
1096  else
1097  !$omp parallel do collapse(2)
1098  do j = 1, ja_org
1099  do i = 1, ia_org
1100  do k = 1, ka_org-2
1101  v_org(k+2,i,j) = v_org(k+2,i,j) / dens_org(k+2,i,j)
1102  end do
1103  end do
1104  end do
1105  end if
1106  end if
1107  end if
1108  if ( .not. exist ) then
1109  log_error("ParentAtmosInputNetCDF",*) 'Either "V" or "MOMY" is required when "U" or "MOMX" exists'
1110  call prc_abort
1111  end if
1112  end if
1113 
1114 
1115  if ( sfc_diagnoses ) then
1116 
1117  !$omp parallel do
1118  do j = 1, ja_org
1119  do i = 1, ia_org
1120  cz_org(1,i,j) = 0.0_rp
1121  end do
1122  end do
1123 
1124 
1125  if ( first_atm .or. update_coord ) then
1126  ! topo
1127  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1128  work2d(:,:), vars_atmos%get("topo"), &
1129  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1130  exist = exist )
1131  if ( exist ) then
1132  !$omp parallel do
1133  do j = 1, ja_org
1134  do i = 1, ia_org
1135  cz_org(2,i,j) = work2d(i,j)
1136  end do
1137  end do
1138  else
1139  !$omp parallel do
1140  do j = 1, ja_org
1141  do i = 1, ia_org
1142  cz_org(2,i,j) = undef
1143  end do
1144  end do
1145  end if
1146  end if
1147 
1148  ! MSLP
1149  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1150  work2d(:,:), vars_atmos%get("MSLP"), &
1151  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1152  exist = exist )
1153  if ( exist ) then
1154  !$omp parallel do
1155  do j = 1, ja_org
1156  do i = 1, ia_org
1157  pres_org(1,i,j) = work2d(i,j)
1158  end do
1159  end do
1160  else
1161  !$omp parallel do
1162  do j = 1, ja_org
1163  do i = 1, ia_org
1164  pres_org(1,i,j) = undef
1165  end do
1166  end do
1167  end if
1168 
1169  ! SFC_PRES
1170  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1171  work2d(:,:), vars_atmos%get("SFC_PRES"), &
1172  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1173  exist = exist )
1174  if ( exist ) then
1175  !$omp parallel do
1176  do j = 1, ja_org
1177  do i = 1, ia_org
1178  pres_org(2,i,j) = work2d(i,j)
1179  end do
1180  end do
1181  else
1182  !$omp parallel do
1183  do j = 1, ja_org
1184  do i = 1, ia_org
1185  pres_org(2,i,j) = undef
1186  end do
1187  end do
1188  end if
1189 
1190  ! U10, V10
1191  if ( uvmet ) then
1192  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1193  work2d(:,:), vars_atmos%get("U10met"), &
1194  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1195  exist = exist )
1196  if ( exist ) then
1197  !$omp parallel do
1198  do j = 1, ja_org
1199  do i = 1, ia_org
1200  u_org(2,i,j) = work2d(i,j)
1201  end do
1202  end do
1203  else
1204  !$omp parallel do
1205  do j = 1, ja_org
1206  do i = 1, ia_org
1207  u_org(2,i,j) = undef
1208  end do
1209  end do
1210  end if
1211  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1212  work2d(:,:), vars_atmos%get("V10met"), &
1213  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1214  exist = exist )
1215  if ( exist ) then
1216  !$omp parallel do
1217  do j = 1, ja_org
1218  do i = 1, ia_org
1219  v_org(2,i,j) = work2d(i,j)
1220  end do
1221  end do
1222  else
1223  !$omp parallel do
1224  do j = 1, ja_org
1225  do i = 1, ia_org
1226  v_org(2,i,j) = undef
1227  end do
1228  end do
1229  end if
1230  else
1231  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1232  work2d(:,:), vars_atmos%get("U10"), &
1233  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1234  exist = exist )
1235  if ( exist ) then
1236  !$omp parallel do
1237  do j = 1, ja_org
1238  do i = 1, ia_org
1239  u_org(2,i,j) = work2d(i,j)
1240  end do
1241  end do
1242  else
1243  !$omp parallel do
1244  do j = 1, ja_org
1245  do i = 1, ia_org
1246  u_org(2,i,j) = undef
1247  end do
1248  end do
1249  end if
1250  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1251  work2d(:,:), vars_atmos%get("V10"), &
1252  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1253  exist = exist )
1254  if ( exist ) then
1255  !$omp parallel do
1256  do j = 1, ja_org
1257  do i = 1, ia_org
1258  v_org(2,i,j) = work2d(i,j)
1259  end do
1260  end do
1261  else
1262  !$omp parallel do
1263  do j = 1, ja_org
1264  do i = 1, ia_org
1265  v_org(2,i,j) = undef
1266  end do
1267  end do
1268  end if
1269  end if
1270 
1271  ! T2
1272  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1273  work2d(:,:), vars_atmos%get("T2"), &
1274  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1275  exist = exist )
1276  if ( exist ) then
1277  !$omp parallel do
1278  do j = 1, ja_org
1279  do i = 1, ia_org
1280  temp_org(2,i,j) = work2d(i,j)
1281  end do
1282  end do
1283  end if
1284 
1285  ! Q2
1286  if ( rh2qv ) then
1287  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1288  work2d(:,:), vars_atmos%get("RH2"), &
1289  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1290  exist = exist )
1291  if ( exist ) then
1292  !$omp parallel do
1293  do j = 1, ja_org
1294  do i = 1, ia_org
1295  rh_org(2,i,j) = work2d(i,j)
1296  end do
1297  end do
1298  else
1299  !$omp parallel do
1300  do j = 1, ja_org
1301  do i = 1, ia_org
1302  rh_org(2,i,j) = undef
1303  end do
1304  end do
1305  end if
1306  else
1307  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1308  work2d(:,:), vars_atmos%get("Q2"), &
1309  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1310  exist = exist )
1311  if ( exist ) then
1312  if ( same_mp_type ) then
1313  !$omp parallel do
1314  do j = 1, ja_org
1315  do i = 1, ia_org
1316  qtrc_org(2,i,j,qs_mp) = work2d(i,j)
1317  end do
1318  end do
1319  else
1320  !$omp parallel do
1321  do j = 1, ja_org
1322  do i = 1, ia_org
1323  qv_org(2,i,j) = work2d(i,j)
1324  end do
1325  end do
1326  end if
1327  else
1328  !$omp parallel do
1329  do j = 1, ja_org
1330  do i = 1, ia_org
1331  qv_org(2,i,j) = undef
1332  end do
1333  end do
1334  end if
1335  end if
1336 
1337 
1338  else
1339 
1340  !$omp parallel do
1341  do j = 1, ja_org
1342  do i = 1, ia_org
1343  cz_org(1,i,j) = 0.0_rp
1344  cz_org(2,i,j) = 0.0_rp
1345  pres_org(1,i,j) = undef
1346  pres_org(2,i,j) = undef
1347  u_org(2,i,j) = undef
1348  v_org(2,i,j) = undef
1349  temp_org(2,i,j) = undef
1350  pt_org(2,i,j) = undef
1351  qv_org(2,i,j) = undef
1352  rh_org(2,i,j) = undef
1353  end do
1354  end do
1355 
1356  end if
1357 
1358  first_atm = .false.
1359 
1360  return
1361  end subroutine parentatmosinputnetcdf
1362 
1363  !-----------------------------------------------------------------------------
1365  subroutine parentlandsetupnetcdf( &
1366  ldims, &
1367  timelen, &
1368  lon_all, &
1369  lat_all, &
1370  basename_org, &
1371  basename_num, &
1372  use_file_landwater, &
1373  serial, &
1374  do_read )
1375  use scale_const, only: &
1376  d2r => const_d2r
1377  use scale_file, only: &
1379  use scale_comm_cartesc, only: &
1380  comm_bcast
1381  use scale_comm_cartesc_nest, only: &
1384 
1385  implicit none
1386 
1387  integer, intent(out) :: ldims(3)
1388  integer, intent(out) :: timelen
1389  real(rp), allocatable, intent(out) :: lon_all(:,:)
1390  real(rp), allocatable, intent(out) :: lat_all(:,:)
1391 
1392  character(len=*), intent(in) :: basename_org
1393  character(len=*), intent(in) :: basename_num
1394  logical, intent(in) :: use_file_landwater
1395 
1396  logical, intent(inout) :: serial
1397  logical, intent(inout) :: do_read
1398 
1399  character(len=8) :: file_type = "AUTO" ! "SCALE-RM", "WRFARW", "NAMELIST", "AUTO"
1400  character(len=FILE_HLONG) :: nm_file
1401  logical :: scale_multi_file = .true.
1402  integer :: scale_parent_prc_num_x
1403  integer :: scale_parent_prc_num_y
1404  character(len=FILE_HLONG) :: scale_latlon_catalogue
1405 
1406  namelist / param_mkinit_real_land_netcdf / &
1407  file_type, &
1408  nm_file, &
1409  scale_parent_prc_num_x, &
1410  scale_parent_prc_num_y, &
1411  scale_latlon_catalogue
1412 
1413  character(len=FILE_HLONG) :: basename
1414  character(len=FILE_HLONG) :: fname
1415  integer :: nmfid
1416 
1417  character(len=32) :: items(vars_max)
1418  integer :: nvars
1419  type(vinfo), pointer :: var_info
1420  class(*), pointer :: v
1421 
1422  logical :: error, exist
1423  integer :: n, i
1424  integer :: ierr
1425  !---------------------------------------------------------------------------
1426 
1427  log_info("ParentLandSetupNetCDF",*) 'Real Case/Land Setup'
1428 
1429  file_type = "AUTO"
1430  nm_file = ""
1431  scale_parent_prc_num_x = -1
1432  scale_parent_prc_num_y = -1
1433  scale_latlon_catalogue = ""
1434 
1435  !--- read namelist
1436  rewind(io_fid_conf)
1437  read(io_fid_conf,nml=param_mkinit_real_land_netcdf,iostat=ierr)
1438  if( ierr > 0 ) then
1439  log_error("ParentLandSetupNetCDF",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_LAND_NetCDF. Check!'
1440  call prc_abort
1441  endif
1442  log_nml(param_mkinit_real_land_netcdf)
1443 
1444 
1445  basename = trim(basename_org) // trim(basename_num)
1446 
1447  fid_lnd = -1
1448  if ( do_read ) then
1449  call check_filetype(fid_lnd, file_type, basename, scale_tile_lnd, "ParentLandSetupNetCDF")
1450  end if
1451 
1452  call comm_bcast( file_type )
1453 
1454  if ( file_type == "SCALE-RM" ) then
1455  call comm_bcast( scale_tile_lnd )
1456  if ( scale_tile_lnd ) then
1457  do_read = .true.
1458  serial = .false.
1459  end if
1460  end if
1461 
1462  if ( do_read ) then
1463  vars_land = hash_table()
1464 
1465  select case( file_type )
1466  case ( "SCALE-RM" )
1467  zname = "lz"
1468  xname = "x"
1469  yname = "y"
1470  tname = "time"
1471 
1472  call vars_land%put("lon", vinfo("lon"))
1473  call vars_land%put("lat", vinfo("lat"))
1474  call vars_land%put("lz", vinfo("lz"))
1475 
1476  call vars_land%put("topo", vinfo("topo"))
1477  call vars_land%put("lsmask", vinfo("lsmask"))
1478 
1479  call vars_land%put("LAND_TEMP", vinfo("LAND_TEMP"))
1480  if ( use_file_landwater ) then
1481  call vars_land%put("LAND_WATER", vinfo("LAND_WATER"))
1482  end if
1483 
1484  call vars_land%put("LAND_SFC_TEMP", vinfo("LAND_SFC_TEMP"))
1485 
1486  call vars_land%put("LAND_SFC_ALB_IR_dir", vinfo("LAND_SFC_ALB_IR_dir"))
1487  call vars_land%put("LAND_SFC_ALB_IR_dif", vinfo("LAND_SFC_ALB_IR_dif"))
1488  call vars_land%put("LAND_SFC_ALB_NIR_dir", vinfo("LAND_SFC_ALB_NIR_dir"))
1489  call vars_land%put("LAND_SFC_ALB_NIR_dif", vinfo("LAND_SFC_ALB_NIR_dif"))
1490  call vars_land%put("LAND_SFC_ALB_VIS_dir", vinfo("LAND_SFC_ALB_VIS_dir"))
1491  call vars_land%put("LAND_SFC_ALB_VIS_dif", vinfo("LAND_SFC_ALB_VIS_dif"))
1492 
1493  call vars_land%put("URBAN_SFC_TEMP", vinfo("URBAN_SFC_TEMP"))
1494 
1495  case ( "WRFARW" )
1496  zname = "soil_layers_stag"
1497  xname = "west_east"
1498  yname = "south_north"
1499  tname = "Time"
1500 
1501  call vars_land%put("lon", vinfo("XLONG"))
1502  call vars_land%put("lat", vinfo("XLAT"))
1503  call vars_land%put("lz", vinfo("ZS"))
1504 
1505  call vars_land%put("topo", vinfo("HGT"))
1506  call vars_land%put("lsmask", vinfo("LANDMASK"))
1507 
1508  call vars_land%put("LAND_TEMP", vinfo("TSLB"))
1509  if ( use_file_landwater ) then
1510  call vars_land%put("LAND_WATER", vinfo("SH2O"))
1511  end if
1512 
1513  call vars_land%put("LAND_SFC_TEMP", vinfo("TSK"))
1514 
1515  call vars_land%put("LAND_SFC_ALB_VIS_dir", vinfo("ALBEDO"))
1516  call vars_land%put("LAND_SFC_EMIS_IR_dif", vinfo("EMISS"))
1517 
1518  call vars_land%put("URBAN_SFC_TEMP", vinfo("URBAN_SFC_TEMP"))
1519 
1520 ! call vars_land%put("SNOW_WATER", vinfo("SNOW"))
1521 ! call vars_land%put("SNOW_TEMP", vinfo("TSNAV"))
1522 
1523  case ( "NAMELIST" )
1524  case default
1525  log_error("ParentLandSetupNetCDF",*) 'FILE_TYPE must be "SCALE-RM", "WRFARW", or "AUTO", ', trim(file_type)
1526  call prc_abort
1527  end select
1528 
1529 
1530  !--- read namelist
1531  nmfid = -1
1532  if ( nm_file /= "" ) then
1533  nmfid = io_get_available_fid()
1534  call io_get_fname(fname, nm_file)
1535  open(nmfid, file=fname, form="formatted", status="old", action="read", iostat=ierr)
1536  if ( ierr /= 0 ) then
1537  log_error("ParentLandSetupNetCDF",*) 'namelist file is not found! ', trim(fname)
1538  call prc_abort
1539  end if
1540 
1541  rewind(nmfid)
1542  read(nmfid, nml=netcdf_dims, iostat=ierr)
1543  if( ierr > 0 ) then
1544  log_error("ParentLandSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_DIMS in ', trim(fname), '. Check!'
1545  call prc_abort
1546  end if
1547 
1548  ! items
1549  rewind(nmfid)
1550  nvars = 0
1551  do n = 1, vars_max
1552  read(nmfid, nml=netcdf_item, iostat=ierr)
1553  if( ierr > 0 ) then
1554  log_error("ParentLandSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_ITEM in ', trim(fname), '. Check!'
1555  call prc_abort
1556  else if( ierr < 0 ) then
1557  exit
1558  end if
1559  nvars = nvars + 1
1560  items(nvars) = item
1561  end do
1562  if ( nvars > vars_max ) then
1563  log_error("ParentLandSetupNetCDF",*) "The number of item in the namelist file exceeds the limit! ", nvars
1564  call prc_abort
1565  end if
1566  rewind(nmfid)
1567  do n = 1, nvars
1568  ! set default
1569  if ( vars_land%has_key(items(n)) ) then
1570  item = items(n)
1571  v => vars_land%get(item)
1572  select type(v)
1573  type is (vinfo)
1574  var_info => v
1575  end select
1576  name = var_info%name
1577  fact = var_info%fact
1578  offset = var_info%offset
1579  else
1580  item = items(n)
1581  name = items(n)
1582  fact = 1.0_rp
1583  offset = 0.0_rp
1584  end if
1585  zstg = .false.
1586  xstg = .false.
1587  ystg = .false.
1588  read(nmfid, nml=netcdf_item, iostat=ierr)
1589  if ( ierr /= 0 ) exit
1590  ! set params
1591  call vars_land%put(item, vinfo(name=name, zstg=zstg, xstg=xstg, ystg=ystg, fact=fact, offset=offset))
1592  end do
1593 
1594  else if ( file_type == "NAMELIST" ) then
1595  log_error("ParentLANDSetupNetCDF",*) 'NM_FILE is necessary'
1596  call prc_abort
1597  end if
1598 
1599  end if
1600 
1601  if ( scale_tile_lnd ) then
1602 
1604  scale_domid_lnd, & ! [OUT]
1605  basename, & ! [IN]
1606  scale_parent_prc_num_x, & ! [IN]
1607  scale_parent_prc_num_y, & ! [IN]
1608  scale_latlon_catalogue ) ! [IN]
1609 
1611  scale_domid_lnd, & ! [IN]
1612  lkmax=ldims(1), & ! [OUT]
1613  imaxg=ldims(2), & ! [OUT]
1614  jmaxg=ldims(3), & ! [OUT]
1615  num_tile=nfiles_lnd ) ! [OUT]
1616 
1617  allocate( fids_lnd(nfiles_lnd) )
1618  allocate( tile_id_lnd(nfiles_lnd) )
1619 
1621  scale_domid_lnd, & ! [IN]
1622  tile_id = tile_id_lnd ) ! [OUT]
1623 
1624  call parentlandopennetcdf( basename_org, basename_num )
1625 
1626  else if ( do_read ) then
1627 
1628  call file_get_dimlength( fid_lnd, zname, ldims(1) )
1629  call file_get_dimlength( fid_lnd, xname, ldims(2) )
1630  call file_get_dimlength( fid_lnd, yname, ldims(3) )
1631 
1632  end if
1633 
1634  if ( do_read ) then
1635 
1636  call file_get_dimlength( fid_lnd, tname, timelen, error=error )
1637  if ( error ) timelen = 1
1638 
1639  allocate( lon_all(ldims(2), ldims(3)) )
1640  allocate( lat_all(ldims(2), ldims(3)) )
1641 
1642  call read2d( ldims(2), 1, ldims(2), ldims(3), 1, ldims(3), &
1643  lon_all(:,:), vars_land%get("lon"), &
1644  1, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1645  lon_all(:,:) = lon_all(:,:) * d2r
1646  call read2d( ldims(2), 1, ldims(2), ldims(3), 1, ldims(3), &
1647  lat_all(:,:), vars_land%get("lat"), &
1648  1, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1649  lat_all(:,:) = lat_all(:,:) * d2r
1650 
1651  end if
1652 
1653  first_lnd = .true.
1654 
1655  return
1656  end subroutine parentlandsetupnetcdf
1657 
1658  !-----------------------------------------------------------------------------
1660  subroutine parentlandopennetcdf( &
1661  basename_org, basename_num )
1662  use scale_file, only: &
1663  file_open
1664  implicit none
1665  character(len=*), intent(in) :: basename_org
1666  character(len=*), intent(in) :: basename_num
1667 
1668  character(len=FILE_HLONG) :: basename
1669 
1670  integer :: n
1671  !---------------------------------------------------------------------------
1672 
1673  log_newline
1674  log_info("ParentLandOpenNetCDF",*) 'Real Case/Land Open'
1675 
1676  basename = trim(basename_org) // trim(basename_num)
1677 
1678  if ( scale_tile_lnd ) then
1679 
1680  do n = 1, nfiles_lnd
1681  call file_open( &
1682  basename, & ! [IN]
1683  fids_lnd(n), & ! [OUT]
1684  aggregate=.false., & ! [IN]
1685  allnodes=.false., & ! [IN]
1686  rankid=tile_id_lnd(n) ) ! [IN]
1687  end do
1688 
1689  fid_lnd = fids_lnd(1)
1690 
1691  else
1692 
1693  call file_open(basename, fid_lnd, postfix="")
1694 
1695  end if
1696 
1697  return
1698  end subroutine parentlandopennetcdf
1699 
1700  !-----------------------------------------------------------------------------
1702  subroutine parentlandfinalizenetcdf
1703  implicit none
1704 
1705  !---------------------------------------------------------------------------
1706 
1707  log_info("ParentLandFinalizeNetCDF",*) 'Real Case/Land Finalize'
1708 
1709  if ( allocated(fids_lnd) ) deallocate( fids_lnd )
1710  if ( allocated(tile_id_lnd) ) deallocate( tile_id_lnd )
1711 
1712  call vars_land%destroy()
1713 
1714  scale_domid_lnd = -1
1715  nfiles_lnd = 0
1716  fid_lnd = -1
1717 
1718  return
1719  end subroutine parentlandfinalizenetcdf
1720 
1721  !-----------------------------------------------------------------------------
1722  subroutine parentlandinputnetcdf( &
1723  KA_org, KS_org, KE_org, &
1724  IA_org, IS_org, IE_org, &
1725  JA_org, JS_org, JE_org, &
1726  tg_org, &
1727  strg_org, &
1728  lst_org, &
1729  ust_org, &
1730  albg_org, &
1731  topo_org, &
1732  lmask_org, &
1733  lz_org, &
1734  use_file_landwater, &
1735  ldims, &
1736  it )
1737  use scale_const, only: &
1738  d2r => const_d2r, &
1739  undef => const_undef
1740  implicit none
1741  integer, intent(in) :: ka_org, ks_org, ke_org
1742  integer, intent(in) :: ia_org, is_org, ie_org
1743  integer, intent(in) :: ja_org, js_org, je_org
1744 
1745  real(rp), intent(out) :: tg_org(ka_org,ia_org,ja_org)
1746  real(rp), intent(out) :: strg_org(ka_org,ia_org,ja_org)
1747  real(rp), intent(out) :: lst_org(ia_org,ja_org)
1748  real(rp), intent(out) :: ust_org(ia_org,ja_org)
1749  real(rp), intent(out) :: albg_org(ia_org,ja_org,n_rad_dir,n_rad_rgn)
1750 
1751  real(rp), intent(inout) :: topo_org(ia_org,ja_org)
1752  real(rp), intent(inout) :: lmask_org(ia_org,ja_org)
1753  real(rp), intent(inout) :: lz_org(ka_org)
1754 
1755  logical, intent(in) :: use_file_landwater ! use land water data from files
1756  integer, intent(in) :: ldims(3)
1757  integer, intent(in) :: it
1758 
1759  logical :: exist
1760  integer :: k, i, j
1761  !---------------------------------------------------------------------------
1762 
1763  if ( first_lnd ) then
1764  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1765  topo_org(:,:), vars_land%get("topo"), &
1766  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1767 
1768  call read1d( ka_org, lz_org(:), vars_land%get("lz"), it, fid_lnd )
1769 
1770  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1771  lmask_org(:,:), vars_land%get("lsmask"), &
1772  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1773  end if
1774 
1775  call read3d( ka_org, ks_org, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1776  tg_org(:,:,:), vars_land%get("LAND_TEMP"), &
1777  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1778 
1779 
1780  ! soil liquid water [m3 m-3] (no wrfout-default)
1781  if( use_file_landwater ) then
1782  call read3d( ka_org, ks_org, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1783  strg_org(:,:,:), vars_land%get("LAND_WATER"), &
1784  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1785  endif
1786 
1787  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1788  lst_org(:,:), vars_land%get("LAND_SFC_TEMP"), &
1789  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1790 
1791  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1792  ust_org(:,:), vars_land%get("URBAN_SFC_TEMP"), &
1793  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1794  exist = exist )
1795  if ( .not. exist ) then
1796  !$omp parallel do
1797  do j = 1, ja_org
1798  do i = 1, ia_org
1799  ust_org(i,j) = lst_org(i,j)
1800  end do
1801  end do
1802  end if
1803 
1804  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1805  albg_org(:,:,i_r_direct,i_r_vis), vars_land%get("LAND_SFC_ALB_VIS_dir"), &
1806  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1807  exist = exist )
1808  if ( .not. exist ) then
1809  !$omp parallel do
1810  do j = 1, ja_org
1811  do i = 1, ia_org
1812  albg_org(i,j,i_r_direct,i_r_vis) = undef
1813  end do
1814  end do
1815  end if
1816  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1817  albg_org(:,:,i_r_diffuse,i_r_vis), vars_land%get("LAND_SFC_ALB_VIS_dif"), &
1818  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1819  exist = exist )
1820  if ( .not. exist ) then
1821  !$omp parallel do
1822  do j = 1, ja_org
1823  do i = 1, ia_org
1824  albg_org(i,j,i_r_diffuse,i_r_vis) = albg_org(i,j,i_r_direct,i_r_vis)
1825  end do
1826  end do
1827  end if
1828  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1829  albg_org(:,:,i_r_direct,i_r_nir), vars_land%get("LAND_SFC_ALB_NIR_dir"), &
1830  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1831  exist = exist )
1832  if ( .not. exist ) then
1833  !$omp parallel do
1834  do j = 1, ja_org
1835  do i = 1, ia_org
1836  albg_org(i,j,i_r_direct,i_r_nir) = undef
1837  end do
1838  end do
1839  end if
1840  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1841  albg_org(:,:,i_r_diffuse,i_r_nir), vars_land%get("LAND_SFC_ALB_NIR_dif"), &
1842  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1843  exist = exist )
1844  if ( .not. exist ) then
1845  !$omp parallel do
1846  do j = 1, ja_org
1847  do i = 1, ia_org
1848  albg_org(i,j,i_r_diffuse,i_r_nir) = albg_org(i,j,i_r_direct,i_r_nir)
1849  end do
1850  end do
1851  end if
1852  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1853  albg_org(:,:,i_r_diffuse,i_r_ir), vars_land%get("LAND_SFC_ALB_IR_dif"), &
1854  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1855  exist = exist )
1856  if ( .not. exist ) then
1857  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1858  albg_org(:,:,i_r_diffuse,i_r_ir), vars_land%get("LAND_SFC_EMIS_IR_dif"), &
1859  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1860  exist = exist )
1861  if ( exist ) then
1862  !$omp parallel do
1863  do j = 1, ja_org
1864  do i = 1, ia_org
1865  albg_org(i,j,i_r_diffuse,i_r_ir) = 1.0_rp - albg_org(i,j,i_r_diffuse,i_r_ir)
1866  end do
1867  end do
1868  else
1869  !$omp parallel do
1870  do j = 1, ja_org
1871  do i = 1, ia_org
1872  albg_org(i,j,i_r_diffuse,i_r_ir) = undef
1873  end do
1874  end do
1875  end if
1876  end if
1877  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1878  albg_org(:,:,i_r_direct,i_r_ir), vars_land%get("LAND_SFC_ALB_IR_dir"), &
1879  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1880  exist = exist )
1881  if ( .not. exist ) then
1882  !$omp parallel do
1883  do j = 1, ja_org
1884  do i = 1, ia_org
1885  albg_org(i,j,i_r_direct,i_r_ir) = albg_org(i,j,i_r_diffuse,i_r_ir)
1886  end do
1887  end do
1888  end if
1889 
1890  first_lnd = .false.
1891 
1892  return
1893  end subroutine parentlandinputnetcdf
1894 
1895  !-----------------------------------------------------------------------------
1897  subroutine parentoceansetupnetcdf( &
1898  odims, &
1899  timelen, &
1900  lon_all, &
1901  lat_all, &
1902  basename_org, &
1903  basename_num, &
1904  serial, &
1905  do_read )
1906  use scale_const, only: &
1907  d2r => const_d2r
1908  use scale_file, only: &
1910  use scale_comm_cartesc, only: &
1911  comm_bcast
1912  use scale_comm_cartesc_nest, only: &
1915 
1916  implicit none
1917 
1918  integer, intent(out) :: odims(2)
1919  integer, intent(out) :: timelen
1920  real(rp), allocatable, intent(out) :: lon_all(:,:)
1921  real(rp), allocatable, intent(out) :: lat_all(:,:)
1922 
1923  character(len=*), intent(in) :: basename_org
1924  character(len=*), intent(in) :: basename_num
1925 
1926  logical, intent(inout) :: serial
1927  logical, intent(inout) :: do_read
1928 
1929  character(len=8) :: file_type = "AUTO" ! "SCALE-RM", "WRFARW", "NAMELIST", "AUTO"
1930  character(len=FILE_HLONG) :: nm_file
1931  logical :: scale_multi_file = .true.
1932  integer :: scale_parent_prc_num_x
1933  integer :: scale_parent_prc_num_y
1934  character(len=FILE_HLONG) :: scale_latlon_catalogue
1935 
1936  namelist / param_mkinit_real_ocean_netcdf / &
1937  file_type, &
1938  nm_file, &
1939  scale_multi_file, &
1940  scale_parent_prc_num_x, &
1941  scale_parent_prc_num_y, &
1942  scale_latlon_catalogue
1943 
1944  character(len=FILE_HLONG) :: basename
1945  character(len=FILE_HLONG) :: fname
1946  integer :: nmfid
1947 
1948  character(len=32) :: items(vars_max)
1949  integer :: nvars
1950  type(vinfo), pointer :: var_info
1951  class(*), pointer :: v
1952 
1953  integer :: n, i
1954  integer :: ierr
1955  logical :: error
1956  !---------------------------------------------------------------------------
1957 
1958  log_info("ParentOceanSetupNetCDF",*) 'Real Case/Ocean Setup'
1959 
1960  file_type = "AUTO"
1961  nm_file = ""
1962  scale_parent_prc_num_x = -1
1963  scale_parent_prc_num_y = -1
1964  scale_latlon_catalogue = ""
1965 
1966  !--- read namelist
1967  rewind(io_fid_conf)
1968  read(io_fid_conf,nml=param_mkinit_real_ocean_netcdf,iostat=ierr)
1969  if( ierr > 0 ) then
1970  log_error("ParentOceanSetupNetCDF",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_OCEAN_NetCDF. Check!'
1971  call prc_abort
1972  endif
1973  log_nml(param_mkinit_real_ocean_netcdf)
1974 
1975 
1976  basename = trim(basename_org) // trim(basename_num)
1977 
1978  fid_ocn = -1
1979  if ( do_read ) then
1980  call check_filetype(fid_ocn, file_type, basename, scale_tile_ocn, "ParentOceanSetupNetCDF")
1981  end if
1982 
1983  call comm_bcast( file_type )
1984 
1985  if ( file_type == "SCALE-RM" ) then
1986  call comm_bcast( scale_tile_ocn )
1987  if ( scale_tile_ocn ) then
1988  do_read = .true.
1989  serial = .false.
1990  end if
1991  end if
1992 
1993 
1994  if ( do_read ) then
1995  vars_ocean = hash_table()
1996 
1997  select case( file_type )
1998  case ( "SCALE-RM" )
1999  xname = "x"
2000  yname = "y"
2001  tname = "time"
2002 
2003  call vars_ocean%put("lon", vinfo("lon"))
2004  call vars_ocean%put("lat", vinfo("lat"))
2005 
2006  call vars_ocean%put("lsmask", vinfo("lsmask"))
2007 
2008  call vars_ocean%put("OCEAN_TEMP", vinfo("OCEAN_TEMP"))
2009 
2010  call vars_ocean%put("OCEAN_SFC_TEMP", vinfo("OCEAN_SFC_TEMP"))
2011  call vars_ocean%put("OCEAN_SFC_Z0M", vinfo("OCEAN_SFC_Z0M"))
2012 
2013  call vars_ocean%put("OCEAN_SFC_ALB_IR_dir", vinfo("OCEAN_SFC_ALB_IR_dir"))
2014  call vars_ocean%put("OCEAN_SFC_ALB_IR_dif", vinfo("OCEAN_SFC_ALB_IR_dif"))
2015  call vars_ocean%put("OCEAN_SFC_ALB_NIR_dir", vinfo("OCEAN_SFC_ALB_NIR_dir"))
2016  call vars_ocean%put("OCEAN_SFC_ALB_NIR_dif", vinfo("OCEAN_SFC_ALB_NIR_dif"))
2017  call vars_ocean%put("OCEAN_SFC_ALB_VIS_dir", vinfo("OCEAN_SFC_ALB_VIS_dir"))
2018  call vars_ocean%put("OCEAN_SFC_ALB_VIS_dif", vinfo("OCEAN_SFC_ALB_VIS_dif"))
2019 
2020  case ( "WRFARW" )
2021  xname = "west_east"
2022  yname = "south_north"
2023  tname = "Time"
2024 
2025  call vars_ocean%put("lon", vinfo("XLONG"))
2026  call vars_ocean%put("lat", vinfo("XLAT"))
2027  call vars_ocean%put("lz", vinfo("ZS"))
2028 
2029  call vars_ocean%put("topo", vinfo("HGT"))
2030  call vars_ocean%put("lsmask", vinfo("LANDMASK"))
2031 
2032  call vars_ocean%put("OCEAN_TEMP", vinfo("OCEAN_TEMP"))
2033 
2034  call vars_ocean%put("OCEAN_SFC_TEMP", vinfo("SST"))
2035  call vars_ocean%put("OCEAN_SFC_Z0M", vinfo("ZNT"))
2036 
2037  call vars_ocean%put("OCEAN_SFC_ALB_VIS_dir", vinfo("ALBEDO"))
2038  call vars_ocean%put("OCEAN_SFC_EMIS_IR_dif", vinfo("EMISS"))
2039 
2040  case ( "NAMELIST" )
2041  case default
2042  log_error("ParentOCEANSetupNetCDF",*) 'FILE_TYPE must be "SCALE-RM", "WRFARW", "NAMELIST", or "AUTO", ', trim(file_type)
2043  call prc_abort
2044  end select
2045 
2046 
2047  !--- read namelist
2048  nmfid = -1
2049  if ( nm_file /= "" ) then
2050  nmfid = io_get_available_fid()
2051  call io_get_fname(fname, nm_file)
2052  open(nmfid, file=fname, form="formatted", status="old", action="read", iostat=ierr)
2053  if ( ierr /= 0 ) then
2054  log_error("ParentOceanSetupNetCDF",*) 'namelist file is not found! ', trim(fname)
2055  call prc_abort
2056  end if
2057 
2058  rewind(nmfid)
2059  read(nmfid, nml=netcdf_dims, iostat=ierr)
2060  if( ierr > 0 ) then
2061  log_error("ParentOceanSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_DIMS in ', trim(fname), '. Check!'
2062  call prc_abort
2063  end if
2064 
2065  ! items
2066  rewind(nmfid)
2067  nvars = 0
2068  do n = 1, vars_max
2069  read(nmfid, nml=netcdf_item, iostat=ierr)
2070  if( ierr > 0 ) then
2071  log_error("ParentLandSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_ITEM in ', trim(fname), '. Check!'
2072  call prc_abort
2073  else if( ierr < 0 ) then
2074  exit
2075  end if
2076  nvars = nvars + 1
2077  items(nvars) = item
2078  end do
2079  if ( nvars > vars_max ) then
2080  log_error("ParentLandSetupNetCDF",*) "The number of item in the namelist file exceeds the limit! ", nvars
2081  call prc_abort
2082  end if
2083  rewind(nmfid)
2084  do n = 1, nvars
2085  ! set default
2086  if ( vars_ocean%has_key(items(n)) ) then
2087  item = items(n)
2088  v => vars_ocean%get(item)
2089  select type(v)
2090  type is (vinfo)
2091  var_info => v
2092  end select
2093  name = var_info%name
2094  fact = var_info%fact
2095  offset = var_info%offset
2096  else
2097  item = items(n)
2098  name = items(n)
2099  fact = 1.0_rp
2100  offset = 0.0_rp
2101  end if
2102  zstg = .false.
2103  xstg = .false.
2104  ystg = .false.
2105  read(nmfid, nml=netcdf_item, iostat=ierr)
2106  if ( ierr /= 0 ) exit
2107  ! set params
2108  call vars_ocean%put(item, vinfo(name=name, zstg=zstg, xstg=xstg, ystg=ystg, fact=fact, offset=offset))
2109  end do
2110 
2111  else if ( file_type == "NAMELIST" ) then
2112  log_error("ParentLANDSetupNetCDF",*) 'NM_FILE is necessary'
2113  call prc_abort
2114  end if
2115 
2116  end if
2117 
2118  if ( scale_tile_ocn ) then
2119 
2121  scale_domid_ocn, & ! [OUT]
2122  basename, & ! [IN]
2123  scale_parent_prc_num_x, & ! [IN]
2124  scale_parent_prc_num_y, & ! [IN]
2125  scale_latlon_catalogue ) ! [IN]
2126 
2128  scale_domid_ocn, & ! [IN]
2129  imaxg=odims(1), & ! [OUT]
2130  jmaxg=odims(2), & ! [OUT]
2131  num_tile=nfiles_ocn ) ! [OUT]
2132 
2133  allocate( fids_ocn(nfiles_ocn) )
2134  allocate( tile_id_ocn(nfiles_ocn) )
2135 
2137  scale_domid_ocn, & ! [IN]
2138  tile_id = tile_id_ocn ) ! [OUT]
2139 
2140  call parentoceanopennetcdf( basename_org, basename_num )
2141 
2142  else if ( do_read ) then
2143 
2144  call file_get_dimlength( fid_ocn, xname, odims(1) )
2145  call file_get_dimlength( fid_ocn, yname, odims(2) )
2146 
2147  end if
2148 
2149  if ( do_read ) then
2150 
2151  call file_get_dimlength( fid_ocn, tname, timelen, error=error )
2152  if ( error ) timelen = 1
2153 
2154  allocate( lon_all(odims(1),odims(2)) )
2155  allocate( lat_all(odims(1),odims(2)) )
2156 
2157  call read2d( odims(1), 1, odims(1), odims(2), 1, odims(2), &
2158  lon_all(:,:), vars_ocean%get("lon"), &
2159  1, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2160  lon_all(:,:) = lon_all(:,:) * d2r
2161  call read2d( odims(1), 1, odims(1), odims(2), 1, odims(2), &
2162  lat_all(:,:), vars_ocean%get("lat"), &
2163  1, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2164  lat_all(:,:) = lat_all(:,:) * d2r
2165 
2166  end if
2167 
2168  first_ocn = .true.
2169 
2170  return
2171  end subroutine parentoceansetupnetcdf
2172 
2173  !-----------------------------------------------------------------------------
2175  subroutine parentoceanopennetcdf( &
2176  basename_org, basename_num )
2177  use scale_file, only: &
2178  file_open
2179  implicit none
2180  character(len=*), intent(in) :: basename_org
2181  character(len=*), intent(in) :: basename_num
2182 
2183  character(len=FILE_HLONG) :: basename
2184 
2185  integer :: n
2186  !---------------------------------------------------------------------------
2187 
2188  log_newline
2189  log_info("ParentOceanOpenNetCDF",*) 'Real Case/Ocean Open'
2190 
2191  basename = trim(basename_org) // trim(basename_num)
2192 
2193  if ( scale_tile_ocn ) then
2194 
2195  do n = 1, nfiles_ocn
2196  call file_open( &
2197  basename, & ! [IN]
2198  fids_ocn(n), & ! [OUT]
2199  aggregate=.false., & ! [IN]
2200  allnodes=.false., & ! [IN]
2201  rankid=tile_id_ocn(n) ) ! [IN]
2202  end do
2203 
2204  fid_ocn = fids_ocn(1)
2205 
2206  else
2207 
2208  call file_open(basename, fid_ocn, postfix="")
2209 
2210  end if
2211 
2212  return
2213  end subroutine parentoceanopennetcdf
2214 
2215  !-----------------------------------------------------------------------------
2217  subroutine parentoceanfinalizenetcdf
2218  implicit none
2219  !---------------------------------------------------------------------------
2220 
2221  log_info("ParentOceanFinalizeNetCDF",*) 'Real Case/Ocean Finalize'
2222 
2223  if ( allocated(fids_ocn) ) deallocate( fids_ocn )
2224  if ( allocated(tile_id_ocn) ) deallocate( tile_id_ocn )
2225 
2226  call vars_ocean%destroy()
2227 
2228  scale_domid_ocn = -1
2229  nfiles_ocn = 0
2230  fid_ocn = -1
2231 
2232  return
2233  end subroutine parentoceanfinalizenetcdf
2234 
2235  !-----------------------------------------------------------------------------
2236  subroutine parentoceaninputnetcdf( &
2237  IA_org, IS_org, IE_org, &
2238  JA_org, JS_org, JE_org, &
2239  tw_org, &
2240  sst_org, &
2241  albw_org, &
2242  z0w_org, &
2243  omask_org, &
2244  odims, &
2245  it )
2246  use scale_const, only: &
2247  d2r => const_d2r, &
2248  undef => const_undef
2249  implicit none
2250  integer, intent(in) :: ia_org, is_org, ie_org
2251  integer, intent(in) :: ja_org, js_org, je_org
2252 
2253  real(rp), intent(out) :: tw_org(ia_org,ja_org)
2254  real(rp), intent(out) :: sst_org(ia_org,ja_org)
2255  real(rp), intent(out) :: albw_org(ia_org,ja_org,n_rad_dir,n_rad_rgn)
2256  real(rp), intent(out) :: z0w_org(ia_org,ja_org)
2257  real(rp), intent(inout) :: omask_org(ia_org,ja_org)
2258 
2259  integer, intent(in) :: odims(2)
2260  integer, intent(in) :: it
2261 
2262  logical :: exist
2263  integer :: i, j
2264  !---------------------------------------------------------------------------
2265 
2266  if ( first_ocn ) then
2267  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2268  omask_org(:,:), vars_ocean%get("lsmask"), &
2269  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2270  end if
2271 
2272  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2273  tw_org(:,:), vars_ocean%get("OCEAN_SFC_TEMP"), &
2274  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2275  sst_org(:,:) = tw_org(:,:)
2276 
2277  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2278  z0w_org(:,:), vars_ocean%get("OCEAN_SFC_Z0M"), &
2279  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2280  exist = exist )
2281  if ( .not. exist ) then
2282  !$omp parallel do
2283  do j = 1, ja_org
2284  do i = 1, ia_org
2285  z0w_org(:,:) = undef
2286  end do
2287  end do
2288  end if
2289 
2290  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2291  albw_org(:,:,i_r_direct,i_r_vis), vars_ocean%get("OCEAN_SFC_ALB_VIS_dir"), &
2292  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2293  exist = exist )
2294  if ( .not. exist ) then
2295  !$omp parallel do
2296  do j = 1, ja_org
2297  do i = 1, ia_org
2298  albw_org(i,j,i_r_direct,i_r_vis) = undef
2299  end do
2300  end do
2301  end if
2302  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2303  albw_org(:,:,i_r_diffuse,i_r_vis), vars_ocean%get("OCEAN_SFC_ALB_VIS_dif"), &
2304  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2305  exist = exist )
2306  if ( .not. exist ) then
2307  !$omp parallel do
2308  do j = 1, ja_org
2309  do i = 1, ia_org
2310  albw_org(i,j,i_r_diffuse,i_r_vis) = albw_org(i,j,i_r_direct,i_r_vis)
2311  end do
2312  end do
2313  end if
2314  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2315  albw_org(:,:,i_r_direct,i_r_nir), vars_ocean%get("OCEAN_SFC_ALB_NIR_dir"), &
2316  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2317  exist = exist )
2318  if ( .not. exist ) then
2319  !$omp parallel do
2320  do j = 1, ja_org
2321  do i = 1, ia_org
2322  albw_org(i,j,i_r_direct,i_r_nir) = undef
2323  end do
2324  end do
2325  end if
2326  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2327  albw_org(:,:,i_r_diffuse,i_r_nir), vars_ocean%get("OCEAN_SFC_ALB_NIR_dif"), &
2328  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2329  exist = exist )
2330  if ( .not. exist ) then
2331  !$omp parallel do
2332  do j = 1, ja_org
2333  do i = 1, ia_org
2334  albw_org(i,j,i_r_diffuse,i_r_nir) = albw_org(i,j,i_r_direct,i_r_nir)
2335  end do
2336  end do
2337  end if
2338  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2339  albw_org(:,:,i_r_diffuse,i_r_ir), vars_ocean%get("OCEAN_SFC_ALB_IR_dif"), &
2340  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2341  exist = exist )
2342  if ( .not. exist ) then
2343  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2344  albw_org(:,:,i_r_diffuse,i_r_ir), vars_ocean%get("OCEAN_SFC_EMIS_IR_dif"), &
2345  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2346  exist = exist )
2347  if ( exist ) then
2348  !$omp parallel do
2349  do j = 1, ja_org
2350  do i = 1, ia_org
2351  albw_org(i,j,i_r_diffuse,i_r_ir) = 1.0_rp - albw_org(i,j,i_r_diffuse,i_r_ir)
2352  end do
2353  end do
2354  else
2355  !$omp parallel do
2356  do j = 1, ja_org
2357  do i = 1, ia_org
2358  albw_org(i,j,i_r_diffuse,i_r_ir) = undef
2359  end do
2360  end do
2361  end if
2362  end if
2363  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2364  albw_org(:,:,i_r_direct,i_r_ir), vars_ocean%get("OCEAN_SFC_ALB_IR_dir"), &
2365  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2366  exist = exist )
2367  if ( .not. exist ) then
2368  !$omp parallel do
2369  do j = 1, ja_org
2370  do i = 1, ia_org
2371  albw_org(i,j,i_r_direct,i_r_ir) = albw_org(i,j,i_r_diffuse,i_r_ir)
2372  end do
2373  end do
2374  end if
2375 
2376  first_ocn = .false.
2377 
2378  return
2379  end subroutine parentoceaninputnetcdf
2380 
2381  ! private
2382 
2383  subroutine check_filetype(fid, FILE_TYPE, basename_org, SCALE_tile, subname)
2384  use scale_file, only: &
2385  file_open, &
2386  file_get_attribute
2387  integer, intent(out) :: fid
2388  character(len=*), intent(inout) :: FILE_TYPE
2389  logical, intent(out) :: SCALE_tile
2390  character(len=*), intent(in) :: basename_org
2391  character(len=*), intent(in) :: subname
2392 
2393  character(len=FILE_HLONG) :: fname
2394  character(len=32) :: att
2395  logical :: exist
2396  integer :: i
2397 
2398  fname = basename_org
2399  inquire(file=fname, exist=exist)
2400  if ( .not. exist ) then
2401  fname = trim(basename_org)//".nc"
2402  inquire(file=fname, exist=exist)
2403  end if
2404  if ( .not. exist ) then
2405  fname = trim(basename_org)//".pe000000.nc"
2406  inquire(file=fname, exist=exist)
2407  end if
2408  if ( .not. exist ) then
2409  log_error(subname,*) "file is not found: ", trim(basename_org)
2410  call prc_abort
2411  end if
2412  call file_open(fname, fid, postfix="", allnodes=.false.)
2413  if ( file_type == "AUTO" ) then
2414  call file_get_attribute( &
2415  fid, "global", "source", &
2416  att, &
2417  existed = exist )
2418  if ( exist .and. att(:8)=="SCALE-RM" ) then
2419  file_type = "SCALE-RM"
2420  log_info(subname,*) 'FILE-TYPE SCALE-RM was detected'
2421  else
2422  call file_get_attribute( &
2423  fid, "global", "TITLE", &
2424  att, &
2425  existed = exist )
2426  if ( exist .and. index(att, "WRF") > 0 ) then
2427  file_type = "WRFARW"
2428  log_info(subname,*) 'FILE-TYPE WRF was detected'
2429  else
2430  file_type = "NAMELIST"
2431  end if
2432  end if
2433  end if
2434 
2435  scale_tile = .false.
2436  if ( file_type == "SCALE-RM" ) then
2437  call file_get_attribute( &
2438  fid, "global", "scale_cartesC_prc_num_x", &
2439  i, &
2440  existed = exist )
2441  if ( exist .and. i > 1 ) then
2442  scale_tile = .true.
2443  log_info(subname,*) 'Multi files was detected'
2444  else
2445  call file_get_attribute( &
2446  fid, "global", "scale_cartesC_prc_num_y", &
2447  i, &
2448  existed = exist )
2449  if ( exist .and. i > 1 ) then
2450  scale_tile = .true.
2451  log_info(subname,*) 'Multi files was detected'
2452  end if
2453  end if
2454  end if
2455 
2456  return
2457  end subroutine check_filetype
2458 
2459  subroutine read3d( &
2460  KA_org, KS_org, KE_org, &
2461  IA_org, IS_org, IE_org, &
2462  JA_org, JS_org, JE_org, &
2463  val, &
2464  var, &
2465  it, &
2466  nfiles, fid, fids, &
2467  scale_tile, scale_domid, &
2468  exist )
2469  use scale_file, only: &
2470  file_get_datainfo, &
2471  file_get_shape, &
2472  file_read
2473  use scale_comm_cartesc_nest, only: &
2475  integer, intent(in) :: KA_org, KS_org, KE_org
2476  integer, intent(in) :: IA_org, IS_org, IE_org
2477  integer, intent(in) :: JA_org, JS_org, JE_org
2478 
2479  real(RP), intent(out), target :: val(KA_org,IA_org,JA_org)
2480 
2481  class(*), pointer, intent(in) :: var
2482  integer, intent(in) :: it
2483  integer, intent(in) :: nfiles
2484  integer, intent(in) :: fid, fids(nfiles)
2485  logical, intent(in) :: scale_tile
2486  integer, intent(in) :: scale_domid
2487 
2488  logical, intent(out), optional :: exist
2489 
2490  real(RP), allocatable :: buf3d(:,:,:)
2491  real(RP), pointer :: work(:,:,:)
2492  real(RP), allocatable, target :: work_t(:,:,:)
2493 
2494  integer :: dims(3)
2495  integer :: tilei, tilej
2496  integer :: kmax
2497  integer :: cxs, cxe, cys, cye
2498  integer :: pxs, pxe, pys, pye
2499  logical :: has_tdim
2500  logical :: transpose
2501  logical :: exist_
2502  integer :: i0, i1, j0, j1
2503  integer :: kst, ist, jst
2504  integer :: it_
2505  integer :: k, i, j, n
2506 
2507  if ( .not. associated(var) ) then
2508  if ( present(exist) ) then
2509  exist = .false.
2510  else
2511  log_error("read3d",*) 'data is not found '
2512  call prc_abort
2513  end if
2514  return
2515  end if
2516 
2517  select type(var)
2518  type is (vinfo)
2519  if ( var%name == "" ) then
2520  if ( present(exist) ) then
2521  exist = .false.
2522  else
2523  log_error("read3d",*) 'data is not found '
2524  call prc_abort
2525  end if
2526  return
2527  end if
2528 
2529  call file_get_datainfo( fid, var%name, has_tdim=has_tdim, existed=exist_ )
2530  if ( .not. exist_ ) then
2531  if ( present(exist) ) then
2532  exist = .false.
2533  return
2534  else
2535  log_error("read3d",*) 'data is not found: ', trim(var%name)
2536  call prc_abort
2537  end if
2538  end if
2539 
2540  if ( has_tdim ) then
2541  it_ = it
2542  else
2543  it_ = 1
2544  end if
2545 
2546  kmax = ke_org - ks_org + 1
2547 
2548  if ( var%zstg ) then
2549  kst = 1
2550  else
2551  kst = 0
2552  end if
2553  if ( var%xstg ) then
2554  ist = 1
2555  else
2556  ist = 0
2557  end if
2558  if ( var%ystg ) then
2559  jst = 1
2560  else
2561  jst = 0
2562  end if
2563 
2564  call file_get_shape( fid, var%name, dims(:) )
2565  transpose = dims(1) .ne. kmax+kst
2566 
2567  if ( scale_tile ) then
2568  if ( var%xstg .or. var%ystg ) then
2569  allocate( work_t(ka_org,ia_org+ist,ja_org+jst) )
2570  work => work_t
2571  else
2572  work => val
2573  end if
2574  do n = 1, nfiles
2576  tilei, tilej, &
2577  cxs, cxe, cys, cye, &
2578  pxs, pxe, pys, pye, &
2579  scale_domid, n, &
2580  xstg = var%xstg, &
2581  ystg = var%ystg )
2582  i0 = max(is_org - cxs, 0)
2583  i1 = max(cxe - ie_org - ist, 0)
2584  j0 = max(js_org - cys, 0)
2585  j1 = max(cye - je_org - jst, 0)
2586  if ( pxs+i0 > pxe-i1 .or. pys+j0 > pye-j1 ) cycle
2587  if ( transpose ) then
2588  allocate( buf3d(pxs+i0:pxe-i1,pys+j0:pye-j1,ks_org:ke_org+kst) )
2589  call file_read( fids(n), var%name, buf3d(:,:,:), &
2590  step=it_, start=(/pxs+i0,pys+j0,1/), count=(/pxe-pxs+1-i1-i0,pye-pys+1-j1-j0,kmax+kst/) )
2591  if ( var%zstg ) then
2592  !$omp parallel do
2593  do j = j0, pye-pys-j1
2594  do i = i0, pxe-pxs-i1
2595  do k = ks_org, ke_org
2596  work(k,cxs+i-is_org+1,cys+j-js_org+1) = ( buf3d(pxs+i,pys+j,k) + buf3d(pxs+i,pys+j,k+1) ) * 0.5_rp * var%fact + var%offset
2597  end do
2598  end do
2599  end do
2600  else
2601  !$omp parallel do
2602  do j = j0, pye-pys-j1
2603  do i = i0, pxe-pxs-i1
2604  do k = ks_org, ke_org
2605  work(k,cxs+i-is_org+1,cys+j-js_org+1) = buf3d(pxs+i,pys+j,k) * var%fact + var%offset
2606  end do
2607  end do
2608  end do
2609  end if
2610  deallocate( buf3d )
2611  if ( var%xstg .and. cxs==2 .and. is_org==1 ) then ! tentative
2612  !$omp parallel do
2613  do j = j0, pye-pys-j1
2614  do k = ks_org, ke_org
2615  work(k,1,cys+j-js_org+1) = work(k,2,cys+j-js_org+1)
2616  end do
2617  end do
2618  end if
2619  if ( var%ystg .and. cys==2 .and. js_org==1 ) then ! tentative
2620  !$omp parallel do
2621  do i = i0, pxe-pxs-i1
2622  do k = ks_org, ke_org
2623  work(k,cxs+i-is_org+1,1) = work(k,cxs+i-is_org+1,2)
2624  end do
2625  end do
2626  end if
2627  else
2628  allocate( buf3d(ks_org:ke_org+kst,pxs+i0:pxe-i1,pys+j0:pye-j1) )
2629  call file_read( fids(n), var%name, buf3d(:,:,:), &
2630  step=it_, start=(/1,pxs+i0,pys+j0/), count=(/kmax+kst,pxe-pxs+1-i1-i0,pye-pys+1-j1-j0/) )
2631  if ( var%zstg ) then
2632  !$omp parallel do
2633  do j = j0, pye-pys-j1
2634  do i = i0, pxe-pxs-i1
2635  do k = ks_org, ke_org
2636  work(k,cxs+i-is_org+1,cys+j-js_org+1) = ( buf3d(k,pxs+i,pys+j) + buf3d(k+1,pxs+i,pys+j) ) * 0.5_rp * var%fact + var%offset
2637  end do
2638  end do
2639  end do
2640  else
2641  !$omp parallel do
2642  do j = j0, pye-pys-j1
2643  do i = i0, pxe-pxs-i1
2644  do k = ks_org, ke_org
2645  work(k,cxs+i-is_org+1,cys+j-js_org+1) = buf3d(k,pxs+i,pys+j) * var%fact + var%offset
2646  end do
2647  end do
2648  end do
2649  end if
2650  deallocate( buf3d )
2651  end if
2652  end do
2653  if ( var%xstg ) then
2654  !$omp parallel do collapse(2)
2655  do j = 1, ja_org
2656  do i = 1, ia_org
2657  do k = ks_org, ke_org
2658  val(k,i,j) = ( work(k,i,j) + work(k,i+1,j) ) * 0.5_rp
2659  end do
2660  end do
2661  end do
2662  else if ( var%ystg ) then
2663  !$omp parallel do collapse(2)
2664  do j = 1, ja_org
2665  do i = 1, ia_org
2666  do k = ks_org, ke_org
2667  val(k,i,j) = ( work(k,i,j) + work(k,i,j+1) ) * 0.5_rp
2668  end do
2669  end do
2670  end do
2671  end if
2672  if ( var%xstg .or. var%ystg ) then
2673  deallocate( work_t )
2674  end if
2675  nullify( work )
2676  else
2677  if ( transpose ) then
2678  allocate( buf3d(is_org:ie_org+ist,js_org:je_org+jst,ks_org:ke_org+kst) )
2679  call file_read( &
2680  fid, var%name, &
2681  buf3d(:,:,:), &
2682  step=it_, &
2683  start=(/is_org,js_org,1/), &
2684  count=(/ia_org+ist,ja_org+jst,kmax+kst/))
2685  if ( var%zstg ) then
2686  !$omp parallel do collapse(2)
2687  do j = 1, ja_org
2688  do i = 1, ia_org
2689  do k = ks_org, ke_org
2690  val(k,i,j) = ( buf3d(i+is_org-1,j+js_org-1,k) + buf3d(i+is_org-1,j+js_org-1,k+1) ) * 0.5_rp * var%fact + var%offset
2691  end do
2692  end do
2693  end do
2694  else if ( var%xstg ) then
2695  !$omp parallel do collapse(2)
2696  do j = 1, ja_org
2697  do i = 1, ia_org
2698  do k = ks_org, ke_org
2699  val(k,i,j) = ( buf3d(i+is_org-1,j+js_org-1,k) + buf3d(i+is_org,j+js_org-1,k) ) * 0.5_rp * var%fact + var%offset
2700  end do
2701  end do
2702  end do
2703  else if ( var%ystg ) then
2704  !$omp parallel do collapse(2)
2705  do j = 1, ja_org
2706  do i = 1, ia_org
2707  do k = ks_org, ke_org
2708  val(k,i,j) = ( buf3d(i+is_org-1,j+js_org-1,k) + buf3d(i+is_org-1,j+js_org,k) ) * 0.5_rp * var%fact + var%offset
2709  end do
2710  end do
2711  end do
2712  else
2713  !$omp parallel do collapse(2)
2714  do j = 1, ja_org
2715  do i = 1, ia_org
2716  do k = ks_org, ke_org
2717  val(k,i,j) = buf3d(i+is_org-1,j+js_org-1,k) * var%fact + var%offset
2718  end do
2719  end do
2720  end do
2721  end if
2722  deallocate( buf3d )
2723  else
2724  allocate( buf3d(ks_org:ke_org+kst,is_org:ie_org+ist,js_org:je_org+jst) )
2725  call file_read( &
2726  fid, var%name, &
2727  buf3d(:,:,:), &
2728  step=it_, &
2729  start=(/1,is_org,js_org/), &
2730  count=(/kmax+kst,ia_org+ist,ja_org+jst/) )
2731  if ( var%zstg ) then
2732  !$omp parallel do collapse(2)
2733  do j = 1, ja_org
2734  do i = 1, ia_org
2735  do k = ks_org, ke_org
2736  val(k,i,j) = ( buf3d(k,i+is_org-1,j+js_org-1) + buf3d(k+1,i+is_org-1,j+js_org-1) ) * 0.5_rp * var%fact + var%offset
2737  end do
2738  end do
2739  end do
2740  else if ( var%xstg ) then
2741  !$omp parallel do collapse(2)
2742  do j = 1, ja_org
2743  do i = 1, ia_org
2744  do k = ks_org, ke_org
2745  val(k,i,j) = ( buf3d(k,i+is_org-1,j+js_org-1) + buf3d(k,i+is_org,j+js_org-1) ) * 0.5_rp * var%fact + var%offset
2746  end do
2747  end do
2748  end do
2749  else if ( var%ystg ) then
2750  !$omp parallel do collapse(2)
2751  do j = 1, ja_org
2752  do i = 1, ia_org
2753  do k = ks_org, ke_org
2754  val(k,i,j) = ( buf3d(k,i+is_org-1,j+js_org-1) + buf3d(k,i+is_org-1,j+js_org) ) * 0.5_rp * var%fact + var%offset
2755  end do
2756  end do
2757  end do
2758  else
2759  !$omp parallel do collapse(2)
2760  do j = 1, ja_org
2761  do i = 1, ia_org
2762  do k = ks_org, ke_org
2763  val(k,i,j) = buf3d(k,i+is_org-1,j+js_org-1) * var%fact + var%offset
2764  end do
2765  end do
2766  end do
2767  end if
2768  deallocate( buf3d )
2769  end if
2770  end if
2771 
2772  if ( present(exist) ) exist = .true.
2773  end select
2774 
2775  return
2776  end subroutine read3d
2777 
2778  subroutine read2d( &
2779  IA_org, IS_org, IE_org, &
2780  JA_org, JS_org, JE_org, &
2781  val, &
2782  var, &
2783  it, &
2784  nfiles, fid, fids, &
2785  scale_tile, scale_domid, &
2786  exist )
2787  use scale_file, only: &
2788  file_get_datainfo, &
2789  file_get_shape, &
2790  file_read
2791  use scale_comm_cartesc_nest, only: &
2793  integer, intent(in) :: IA_org, IS_org, IE_org
2794  integer, intent(in) :: JA_org, JS_org, JE_org
2795 
2796  real(RP), intent(out), target :: val(IA_org,JA_org)
2797 
2798  class(*), pointer, intent(in) :: var
2799  integer, intent(in) :: it
2800  integer, intent(in) :: nfiles
2801  integer, intent(in) :: fid, fids(nfiles)
2802  logical, intent(in) :: scale_tile
2803  integer, intent(in) :: scale_domid
2804 
2805  logical, intent(out), optional :: exist
2806 
2807  real(RP), allocatable :: buf2d(:,:)
2808  real(RP), pointer :: work(:,:)
2809  real(RP), allocatable, target :: work_t(:,:)
2810 
2811  integer :: tilei, tilej
2812  integer :: cxs, cxe, cys, cye
2813  integer :: pxs, pxe, pys, pye
2814 
2815  integer :: dims(2)
2816 
2817  integer :: i0, i1, j0, j1
2818  integer :: ist, jst
2819 
2820  logical :: has_tdim
2821  logical :: exist_
2822  integer :: it_
2823  integer :: n, i, j
2824 
2825  if ( .not. associated(var) ) then
2826  if ( present(exist) ) then
2827  exist = .false.
2828  else
2829  log_error("read2d",*) 'data is not found '
2830  call prc_abort
2831  end if
2832  return
2833  end if
2834 
2835  select type(var)
2836  type is (vinfo)
2837  if ( var%name == "" ) then
2838  if ( present(exist) ) then
2839  exist = .false.
2840  else
2841  log_error("read2d",*) 'data is not found '
2842  call prc_abort
2843  end if
2844  return
2845  end if
2846 
2847  call file_get_datainfo( fid, var%name, has_tdim=has_tdim, existed=exist_ )
2848  if ( .not. exist_ ) then
2849  if ( present(exist) ) then
2850  exist = .false.
2851  return
2852  else
2853  log_error("read2d",*) 'data is not found: ', trim(var%name)
2854  call prc_abort
2855  end if
2856  end if
2857 
2858  if ( has_tdim ) then
2859  it_ = it
2860  else
2861  it_ = 1
2862  end if
2863 
2864  if ( var%xstg ) then
2865  ist = 1
2866  else
2867  ist = 0
2868  end if
2869  if ( var%ystg ) then
2870  jst = 1
2871  else
2872  jst = 0
2873  end if
2874 
2875  if ( scale_domid > 0 ) then
2876  if ( var%xstg .or. var%ystg ) then
2877  allocate( work_t(ia_org+ist,ja_org+jst) )
2878  work => work_t
2879  else
2880  work => val
2881  end if
2882  do n = 1, nfiles
2884  tilei, tilej, &
2885  cxs, cxe, cys, cye, &
2886  pxs, pxe, pys, pye, &
2887  scale_domid, n, &
2888  xstg = var%xstg, &
2889  ystg = var%ystg )
2890  i0 = max(is_org - cxs, 0)
2891  i1 = max(cxe - ie_org - ist, 0)
2892  j0 = max(js_org - cys, 0)
2893  j1 = max(cye - je_org - jst, 0)
2894  if ( pxs+i0 > pxe-i1 .or. pys+j0 > pye-j1 ) cycle
2895  allocate( buf2d(pxs+i0:pxe-i1,pys+j0:pye-j1) )
2896  call file_read( fids(n), var%name, buf2d(:,:), &
2897  step=it_, start=(/pxs+i0,pys+j0/), count=(/pxe-pxs+1-i1-i0,pye-pys+1-j1-j0/) )
2898  !$omp parallel do
2899  do j = j0, pye-pys-j1
2900  do i = i0, pxe-pxs-i1
2901  work(cxs+i-is_org+1,cys+j-js_org+1) = buf2d(pxs+i,pys+j) * var%fact + var%offset
2902  end do
2903  end do
2904  deallocate( buf2d )
2905  if ( var%xstg .and. cxs==2 .and. is_org==1 ) then ! tentative
2906  !$omp parallel do
2907  do j = j0, pye-pys-j1
2908  work(1,cys+j-js_org+1) = work(2,cys+j-js_org+1)
2909  end do
2910  end if
2911  if ( var%ystg .and. cys==2 .and. js_org==1 ) then ! tentative
2912  !$omp parallel do
2913  do i = i0, pxe-pxs-i1
2914  work(cxs+i-is_org+1,1) = work(cxs+i-is_org+1,2)
2915  end do
2916  end if
2917  end do
2918  if ( var%xstg ) then
2919  !$omp parallel do
2920  do j = 1, ja_org
2921  do i = 1, ia_org
2922  val(i,j) = ( work(i,j) + work(i+1,j) ) * 0.5_rp
2923  end do
2924  end do
2925  else if ( var%ystg ) then
2926  !$omp parallel do
2927  do j = 1, ja_org
2928  do i = 1, ia_org
2929  val(i,j) = ( work(i,j) + work(i,j+1) ) * 0.5_rp
2930  end do
2931  end do
2932  end if
2933  if ( var%xstg .or. var%ystg ) then
2934  deallocate( work_t )
2935  end if
2936  nullify( work )
2937  else
2938  if ( var%xstg .or. var%ystg ) then
2939  allocate( work_t(is_org:ie_org+ist,js_org:je_org+jst) )
2940  work => work_t
2941  else
2942  work => val
2943  end if
2944  call file_read( &
2945  fid, var%name, &
2946  work(:,:), &
2947  step=it_, &
2948  start=(/is_org,js_org/), &
2949  count=(/ia_org+ist,ja_org+jst/) )
2950  if ( var%xstg ) then
2951  !$omp parallel do
2952  do j = 1, ja_org
2953  do i = 1, ia_org
2954  val(i,j) = ( work(i,j) + work(i+1,j) ) * 0.5_rp * var%fact + var%offset
2955  end do
2956  end do
2957  else if ( var%ystg ) then
2958  !$omp parallel do
2959  do j = 1, ja_org
2960  do i = 1, ia_org
2961  val(i,j) = ( work(i,j) + work(i,j+1) ) * 0.5_rp * var%fact + var%offset
2962  end do
2963  end do
2964  else if ( var%fact .ne. 1.0_rp .or. var%offset .ne. 0.0_rp ) then
2965  !$omp parallel do
2966  do j = 1, ja_org
2967  do i = 1, ia_org
2968  val(i,j) = val(i,j) * var%fact + var%offset
2969  end do
2970  end do
2971  end if
2972  if ( var%xstg .or. var%ystg ) then
2973  deallocate( work_t )
2974  end if
2975  nullify( work )
2976  end if
2977 
2978  if ( present(exist) ) exist = .true.
2979  end select
2980 
2981  return
2982  end subroutine read2d
2983 
2984  subroutine read1d( &
2985  KA_org, &
2986  val, &
2987  var, &
2988  it, &
2989  fid, &
2990  exist )
2991  use scale_file, only: &
2992  file_get_datainfo, &
2993  file_read
2994  integer, intent(in) :: KA_org
2995 
2996  real(RP), intent(out) :: val(KA_org)
2997 
2998  class(*), pointer, intent(in) :: var
2999  integer, intent(in) :: it
3000  integer, intent(in) :: fid
3001 
3002  logical, intent(out), optional :: exist
3003 
3004  logical :: has_tdim
3005  logical :: exist_
3006 
3007  integer :: it_
3008 
3009  if ( .not. associated(var) ) then
3010  if ( present(exist) ) then
3011  exist = .false.
3012  else
3013  log_error("read1d",*) 'data is not found '
3014  call prc_abort
3015  end if
3016  return
3017  end if
3018 
3019  select type(var)
3020  type is ( vinfo )
3021  if ( var%name == "" ) then
3022  if ( present(exist) ) then
3023  exist = .false.
3024  else
3025  log_error("read1d",*) 'data is not found '
3026  call prc_abort
3027  end if
3028  return
3029  end if
3030 
3031  call file_get_datainfo( fid, var%name, has_tdim=has_tdim, existed=exist_ )
3032  if ( .not. exist_ ) then
3033  if ( present(exist) ) then
3034  exist = .false.
3035  return
3036  else
3037  log_error("read1d",*) 'data is not found: ', trim(var%name)
3038  call prc_abort
3039  end if
3040  end if
3041 
3042  if ( has_tdim ) then
3043  it_ = it
3044  else
3045  it_ = 1
3046  end if
3047 
3048  call file_read( fid, var%name, val(:), step=it_ )
3049  if ( var%fact .ne. 1.0_rp .or. var%offset .ne. 0.0_rp ) then
3050  val(:) = val(:) * var%fact + var%offset
3051  end if
3052 
3053  if ( present(exist) ) exist = .true.
3054  end select
3055 
3056  return
3057  end subroutine read1d
3058 
3059 end module mod_realinput_netcdf
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
mod_realinput_netcdf::parentoceanopennetcdf
subroutine, public parentoceanopennetcdf(basename_org, basename_num)
Ocean Open.
Definition: mod_realinput_netcdf.F90:2177
scale_cpl_sfc_index::n_rad_dir
integer, parameter, public n_rad_dir
Definition: scale_cpl_sfc_index.F90:36
scale_atmos_hydrometeor::num_name
character(len=h_short), dimension(n_hyd), parameter, public num_name
Definition: scale_atmos_hydrometeor.F90:108
mod_atmos_phy_mp_vars
module Atmosphere / Physics Cloud Microphysics
Definition: mod_atmos_phy_mp_vars.F90:12
mod_realinput_netcdf::parentlandopennetcdf
subroutine, public parentlandopennetcdf(basename_org, basename_num)
Land Open.
Definition: mod_realinput_netcdf.F90:1662
scale_comm_cartesc_nest::comm_cartesc_nest_domain_shape
subroutine, public comm_cartesc_nest_domain_shape(tilei, tilej, cxs, cxe, cys, cye, pxs, pxe, pys, pye, dom_id, iloc, xstg, ystg)
Return shape of ParentDomain at the specified rank (for offline)
Definition: scale_comm_cartesC_nest.F90:1187
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_cpl_sfc_index::i_r_direct
integer, parameter, public i_r_direct
Definition: scale_cpl_sfc_index.F90:37
mod_realinput_netcdf::parentoceanfinalizenetcdf
subroutine, public parentoceanfinalizenetcdf
Ocean Finalize.
Definition: mod_realinput_netcdf.F90:2218
mod_realinput_netcdf::parentlandfinalizenetcdf
subroutine, public parentlandfinalizenetcdf
Land Finalize.
Definition: mod_realinput_netcdf.F90:1703
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
scale_hash
module CONSTANT
Definition: scale_hash.F90:11
scale_atmos_hydrometeor::hyd_name
character(len=h_short), dimension(n_hyd), parameter, public hyd_name
Definition: scale_atmos_hydrometeor.F90:104
scale_file::file_open
subroutine, public file_open(basename, fid, mode, single, allnodes, aggregate, rankid, postfix)
Definition: scale_file.F90:536
scale_cpl_sfc_index::i_r_diffuse
integer, parameter, public i_r_diffuse
Definition: scale_cpl_sfc_index.F90:38
mod_realinput_netcdf::parentlandinputnetcdf
subroutine, public parentlandinputnetcdf(KA_org, KS_org, KE_org, IA_org, IS_org, IE_org, JA_org, JS_org, JE_org, tg_org, strg_org, lst_org, ust_org, albg_org, topo_org, lmask_org, lz_org, use_file_landwater, ldims, it)
Definition: mod_realinput_netcdf.F90:1737
scale_precision
module PRECISION
Definition: scale_precision.F90:14
mod_atmos_phy_mp_vars::qs_mp
integer, public qs_mp
Definition: mod_atmos_phy_mp_vars.F90:79
scale_comm_cartesc_nest
module Communication CartesianC nesting
Definition: scale_comm_cartesC_nest.F90:12
scale_cpl_sfc_index::i_r_ir
integer, parameter, public i_r_ir
Definition: scale_cpl_sfc_index.F90:29
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:91
mod_realinput_netcdf::parentatmosinputnetcdf
subroutine, public parentatmosinputnetcdf(KA_org, KS_org, KE_org, IA_org, IS_org, IE_org, JA_org, JS_org, JE_org, QA, cz_org, w_org, u_org, v_org, pres_org, dens_org, temp_org, pt_org, qtrc_org, qv_org, rh_org, qhyd_org, qnum_org, nopres, nodens, uvmet, temp2pt, rh2qv, qnum_flag, same_mp_type, sfc_diagnoses, update_coord, dims, it)
Definition: mod_realinput_netcdf.F90:727
scale_file::file_get_dimlength
subroutine, public file_get_dimlength(fid, dimname, len, error)
get length of dimension
Definition: scale_file.F90:633
scale_io::io_get_available_fid
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:373
scale_comm_cartesc_nest::comm_cartesc_nest_parent_info
subroutine, public comm_cartesc_nest_parent_info(dom_id, KMAX, LKMAX, IMAXG, JMAXG, num_tile, tile_id)
Return infomation of parent domain (for offline)
Definition: scale_comm_cartesC_nest.F90:1130
mod_realinput_netcdf::read1d
subroutine read1d(KA_org, val, var, it, fid, exist)
Definition: mod_realinput_netcdf.F90:2991
scale_file
module file
Definition: scale_file.F90:15
scale_hash::hash_table
Definition: scale_hash.F90:45
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_realinput_netcdf::check_filetype
subroutine check_filetype(fid, FILE_TYPE, basename_org, SCALE_tile, subname)
Definition: mod_realinput_netcdf.F90:2384
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
scale_io::io_get_fname
subroutine, public io_get_fname(outstr, instr, rank, ext, len)
generate process specific filename
Definition: scale_io.F90:421
scale_io
module STDIO
Definition: scale_io.F90:10
scale_cpl_sfc_index::i_r_nir
integer, parameter, public i_r_nir
Definition: scale_cpl_sfc_index.F90:30
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_comm_cartesc_nest::comm_cartesc_nest_domain_regist_file
subroutine, public comm_cartesc_nest_domain_regist_file(dom_id, PARENT_BASENAME, PARENT_PRC_NUM_X, PARENT_PRC_NUM_Y, LATLON_CATALOGUE_FNAME)
offline setup
Definition: scale_comm_cartesC_nest.F90:708
mod_realinput_netcdf::read3d
subroutine read3d(KA_org, KS_org, KE_org, IA_org, IS_org, IE_org, JA_org, JS_org, JE_org, val, var, it, nfiles, fid, fids, scale_tile, scale_domid, exist)
Definition: mod_realinput_netcdf.F90:2469
mod_realinput_netcdf::parentatmossetupnetcdf
subroutine, public parentatmossetupnetcdf(dims, timelen, mixing_ratio, update_coord, mapping_info, qtrc_flag, lon_all, lat_all, basename_org, basename_num, same_mp_type, pt_dry, serial, do_read)
Atmos Setup.
Definition: mod_realinput_netcdf.F90:150
scale_tracer::tracer_name
character(len=h_short), dimension(qa_max), public tracer_name
Definition: scale_tracer.F90:39
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_file_h
module file_h
Definition: scale_file_h.F90:11
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
mod_realinput_netcdf::parentoceansetupnetcdf
subroutine, public parentoceansetupnetcdf(odims, timelen, lon_all, lat_all, basename_org, basename_num, serial, do_read)
Ocean Setup.
Definition: mod_realinput_netcdf.F90:1906
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_cpl_sfc_index
module coupler / surface-atmospehre
Definition: scale_cpl_sfc_index.F90:11
scale_cpl_sfc_index::i_r_vis
integer, parameter, public i_r_vis
Definition: scale_cpl_sfc_index.F90:31
mod_realinput_netcdf
module REAL input netCDF
Definition: mod_realinput_netcdf.F90:11
mod_realinput_netcdf::parentatmosfinalizenetcdf
subroutine, public parentatmosfinalizenetcdf
Atmos Finalize.
Definition: mod_realinput_netcdf.F90:684
mod_realinput_netcdf::parentlandsetupnetcdf
subroutine, public parentlandsetupnetcdf(ldims, timelen, lon_all, lat_all, basename_org, basename_num, use_file_landwater, serial, do_read)
Land Setup.
Definition: mod_realinput_netcdf.F90:1375
mod_atmos_phy_mp_vars::qe_mp
integer, public qe_mp
Definition: mod_atmos_phy_mp_vars.F90:80
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
scale_cpl_sfc_index::n_rad_rgn
integer, parameter, public n_rad_rgn
Definition: scale_cpl_sfc_index.F90:28
mod_realinput_netcdf::read2d
subroutine read2d(IA_org, IS_org, IE_org, JA_org, JS_org, JE_org, val, var, it, nfiles, fid, fids, scale_tile, scale_domid, exist)
Definition: mod_realinput_netcdf.F90:2787
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:95
mod_realinput_netcdf::parentoceaninputnetcdf
subroutine, public parentoceaninputnetcdf(IA_org, IS_org, IE_org, JA_org, JS_org, JE_org, tw_org, sst_org, albw_org, z0w_org, omask_org, odims, it)
Definition: mod_realinput_netcdf.F90:2246
scale_mapprojection::mappinginfo
Definition: scale_mapprojection.F90:93
scale_prc::prc_ismaster
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:92
mod_realinput_netcdf::parentatmosopennetcdf
subroutine, public parentatmosopennetcdf(basename_org, basename_num)
Atmos Open.
Definition: mod_realinput_netcdf.F90:642