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_mptype, &
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_mptype
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=16) :: 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_easting = 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_mptype ) 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(:), existed=exist )
367  if ( .not. exist ) &
368  call file_get_attribute( fid_atm, map, "standard_parallel", standard_parallel(1:1), 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_mptype ) then
418  log_error("ParentAtmosSetupNetCDF",*) 'same_mptype 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_mptype , &
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_mptype
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_mptype ) 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
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  if ( exist ) then
946  nodens = .false.
947  else
948  nodens = .true.
949  end if
950 
951  ! pt
952  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
953  pt_org(:,:,:), vars_atmos%get("PT"), &
954  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
955  exist = exist )
956  if ( exist ) then
957  temp2pt = .false.
958  else
959  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
960  pt_org(:,:,:), vars_atmos%get("RHOT"), &
961  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
962  exist = exist )
963  if ( exist ) then
964  if ( nodens ) then
965  log_error("ParentAtmosInputNetCDF",*) "DENS is necessary to calculate PT from RHOT"
966  call prc_abort
967  end if
968  !$omp parallel do collapse(2)
969  do j = 1, ja_org
970  do i = 1, ia_org
971  do k = 1, ka_org-2
972  pt_org(k+2,i,j) = pt_org(k+2,i,j) / dens_org(k+2,i,j)
973  end do
974  end do
975  end do
976  temp2pt = .false.
977  else
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  if ( .not. exist ) then
983  log_error("ParentAtmosInputNetCDF",*) '"PT", "RHOT", or "T" is necessary'
984  call prc_abort
985  end if
986  temp2pt = .true.
987  end if
988  end if
989 
990  ! W
991  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
992  w_org(:,:,:), vars_atmos%get("W"), &
993  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
994  exist = exist )
995  if ( .not. exist ) then
996  if ( nodens ) then
997  log_error("ParentAtmosInputNetCDF",*) "DENS is necessary to use MOMZ"
998  call prc_abort
999  end if
1000  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1001  w_org(:,:,:), vars_atmos%get("MOMZ"), &
1002  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1003  exist = exist )
1004  if ( exist ) then
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  else
1014  !$omp parallel do collapse(2)
1015  do j = 1, ja_org
1016  do i = 1, ia_org
1017  do k = 1, ka_org-2
1018  w_org(k+2,i,j) = 0.0_rp
1019  end do
1020  end do
1021  end do
1022  endif
1023  end if
1024 
1025  ! U
1026  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1027  u_org(:,:,:), vars_atmos%get("Umet"), &
1028  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1029  exist = exist )
1030  if ( exist ) then
1031  uvmet = .true.
1032  else
1033  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1034  u_org(:,:,:), vars_atmos%get("U"), &
1035  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1036  exist = exist )
1037  if ( .not. exist ) then
1038  if ( nodens ) then
1039  log_error("ParentAtmosInputNetCDF",*) "DENS is necessary to use MOMX"
1040  call prc_abort
1041  end if
1042  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1043  u_org(:,:,:), vars_atmos%get("MOMX"), &
1044  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1045  exist = exist )
1046  if ( .not. exist ) then
1047  log_error("ParentAtmosInputNetCDF",*) '"Ument", "U", or "MOMX" is necessary'
1048  call prc_abort
1049  end if
1050  !$omp parallel do collapse(2)
1051  do j = 1, ja_org
1052  do i = 1, ia_org
1053  do k = 1, ka_org-2
1054  u_org(k+2,i,j) = u_org(k+2,i,j) / dens_org(k+2,i,j)
1055  end do
1056  end do
1057  end do
1058  end if
1059  uvmet = .false.
1060  end if
1061 
1062  ! V
1063  if ( uvmet ) then
1064  call read3d( ka_org, ks_org+2, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1065  v_org(:,:,:), vars_atmos%get("Vmet"), &
1066  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1067  exist = exist )
1068  if ( .not. exist ) then
1069  log_error("ParentAtmosInputNetCDF",*) "Vmet is required when Umet exists"
1070  call prc_abort
1071  end if
1072  else
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("V"), &
1075  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1076  exist = exist )
1077  if ( .not. exist ) then
1078  if ( nodens ) then
1079  log_error("ParentAtmosInputNetCDF",*) "DENS is necessary to use MOMY"
1080  call prc_abort
1081  end if
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("MOMY"), &
1084  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1085  exist = exist )
1086  if ( .not. exist ) then
1087  log_error("ParentAtmosInputNetCDF",*) '"V" or "MOMY" is required when "U" or "MOMX" exists'
1088  call prc_abort
1089  end if
1090  !$omp parallel do collapse(2)
1091  do j = 1, ja_org
1092  do i = 1, ia_org
1093  do k = 1, ka_org-2
1094  v_org(k+2,i,j) = v_org(k+2,i,j) / dens_org(k+2,i,j)
1095  end do
1096  end do
1097  end do
1098  end if
1099  end if
1100 
1101 
1102  if ( sfc_diagnoses ) then
1103 
1104  !$omp parallel do
1105  do j = 1, ja_org
1106  do i = 1, ia_org
1107  cz_org(1,i,j) = 0.0_rp
1108  end do
1109  end do
1110 
1111 
1112  if ( first_atm .or. update_coord ) then
1113  ! topo
1114  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1115  work2d(:,:), vars_atmos%get("topo"), &
1116  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1117  exist = exist )
1118  if ( exist ) then
1119  !$omp parallel do
1120  do j = 1, ja_org
1121  do i = 1, ia_org
1122  cz_org(2,i,j) = work2d(i,j)
1123  end do
1124  end do
1125  else
1126  !$omp parallel do
1127  do j = 1, ja_org
1128  do i = 1, ia_org
1129  cz_org(2,i,j) = undef
1130  end do
1131  end do
1132  end if
1133  end if
1134 
1135  ! MSLP
1136  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1137  work2d(:,:), vars_atmos%get("MSLP"), &
1138  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1139  exist = exist )
1140  if ( exist ) then
1141  !$omp parallel do
1142  do j = 1, ja_org
1143  do i = 1, ia_org
1144  pres_org(1,i,j) = work2d(i,j)
1145  end do
1146  end do
1147  else
1148  !$omp parallel do
1149  do j = 1, ja_org
1150  do i = 1, ia_org
1151  pres_org(1,i,j) = undef
1152  end do
1153  end do
1154  end if
1155 
1156  ! SFC_PRES
1157  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1158  work2d(:,:), vars_atmos%get("SFC_PRES"), &
1159  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1160  exist = exist )
1161  if ( exist ) then
1162  !$omp parallel do
1163  do j = 1, ja_org
1164  do i = 1, ia_org
1165  pres_org(2,i,j) = work2d(i,j)
1166  end do
1167  end do
1168  else
1169  !$omp parallel do
1170  do j = 1, ja_org
1171  do i = 1, ia_org
1172  pres_org(2,i,j) = undef
1173  end do
1174  end do
1175  end if
1176 
1177  ! U10, V10
1178  if ( uvmet ) then
1179  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1180  work2d(:,:), vars_atmos%get("U10met"), &
1181  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1182  exist = exist )
1183  if ( exist ) then
1184  !$omp parallel do
1185  do j = 1, ja_org
1186  do i = 1, ia_org
1187  u_org(2,i,j) = work2d(i,j)
1188  end do
1189  end do
1190  else
1191  !$omp parallel do
1192  do j = 1, ja_org
1193  do i = 1, ia_org
1194  u_org(2,i,j) = undef
1195  end do
1196  end do
1197  end if
1198  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1199  work2d(:,:), vars_atmos%get("V10met"), &
1200  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1201  exist = exist )
1202  if ( exist ) then
1203  !$omp parallel do
1204  do j = 1, ja_org
1205  do i = 1, ia_org
1206  v_org(2,i,j) = work2d(i,j)
1207  end do
1208  end do
1209  else
1210  !$omp parallel do
1211  do j = 1, ja_org
1212  do i = 1, ia_org
1213  v_org(2,i,j) = undef
1214  end do
1215  end do
1216  end if
1217  else
1218  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1219  work2d(:,:), vars_atmos%get("U10"), &
1220  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1221  exist = exist )
1222  if ( exist ) then
1223  !$omp parallel do
1224  do j = 1, ja_org
1225  do i = 1, ia_org
1226  u_org(2,i,j) = work2d(i,j)
1227  end do
1228  end do
1229  else
1230  !$omp parallel do
1231  do j = 1, ja_org
1232  do i = 1, ia_org
1233  u_org(2,i,j) = undef
1234  end do
1235  end do
1236  end if
1237  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1238  work2d(:,:), vars_atmos%get("V10"), &
1239  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1240  exist = exist )
1241  if ( exist ) then
1242  !$omp parallel do
1243  do j = 1, ja_org
1244  do i = 1, ia_org
1245  v_org(2,i,j) = work2d(i,j)
1246  end do
1247  end do
1248  else
1249  !$omp parallel do
1250  do j = 1, ja_org
1251  do i = 1, ia_org
1252  v_org(2,i,j) = undef
1253  end do
1254  end do
1255  end if
1256  end if
1257 
1258  ! T2
1259  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1260  work2d(:,:), vars_atmos%get("T2"), &
1261  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1262  exist = exist )
1263  if ( exist ) then
1264  !$omp parallel do
1265  do j = 1, ja_org
1266  do i = 1, ia_org
1267  temp_org(2,i,j) = work2d(i,j)
1268  end do
1269  end do
1270  end if
1271 
1272  ! Q2
1273  if ( rh2qv ) then
1274  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1275  work2d(:,:), vars_atmos%get("RH2"), &
1276  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1277  exist = exist )
1278  if ( exist ) then
1279  !$omp parallel do
1280  do j = 1, ja_org
1281  do i = 1, ia_org
1282  rh_org(2,i,j) = work2d(i,j)
1283  end do
1284  end do
1285  else
1286  !$omp parallel do
1287  do j = 1, ja_org
1288  do i = 1, ia_org
1289  rh_org(2,i,j) = undef
1290  end do
1291  end do
1292  end if
1293  else
1294  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1295  work2d(:,:), vars_atmos%get("Q2"), &
1296  it, nfiles_atm, fid_atm, fids_atm, scale_tile_atm, scale_domid_atm, &
1297  exist = exist )
1298  if ( exist ) then
1299  if ( same_mptype ) then
1300  !$omp parallel do
1301  do j = 1, ja_org
1302  do i = 1, ia_org
1303  qtrc_org(2,i,j,qs_mp) = work2d(i,j)
1304  end do
1305  end do
1306  else
1307  !$omp parallel do
1308  do j = 1, ja_org
1309  do i = 1, ia_org
1310  qv_org(2,i,j) = work2d(i,j)
1311  end do
1312  end do
1313  end if
1314  else
1315  !$omp parallel do
1316  do j = 1, ja_org
1317  do i = 1, ia_org
1318  qv_org(2,i,j) = undef
1319  end do
1320  end do
1321  end if
1322  end if
1323 
1324 
1325  else
1326 
1327  !$omp parallel do
1328  do j = 1, ja_org
1329  do i = 1, ia_org
1330  cz_org(1,i,j) = 0.0_rp
1331  cz_org(2,i,j) = 0.0_rp
1332  pres_org(1,i,j) = undef
1333  pres_org(2,i,j) = undef
1334  u_org(2,i,j) = undef
1335  v_org(2,i,j) = undef
1336  temp_org(2,i,j) = undef
1337  pt_org(2,i,j) = undef
1338  qv_org(2,i,j) = undef
1339  rh_org(2,i,j) = undef
1340  end do
1341  end do
1342 
1343  end if
1344 
1345  first_atm = .false.
1346 
1347  return
1348  end subroutine parentatmosinputnetcdf
1349 
1350  !-----------------------------------------------------------------------------
1352  subroutine parentlandsetupnetcdf( &
1353  ldims, &
1354  timelen, &
1355  lon_all, &
1356  lat_all, &
1357  basename_org, &
1358  basename_num, &
1359  use_file_landwater, &
1360  serial, &
1361  do_read )
1362  use scale_const, only: &
1363  d2r => const_d2r
1364  use scale_file, only: &
1366  use scale_comm_cartesc, only: &
1367  comm_bcast
1368  use scale_comm_cartesc_nest, only: &
1371 
1372  implicit none
1373 
1374  integer, intent(out) :: ldims(3)
1375  integer, intent(out) :: timelen
1376  real(rp), allocatable, intent(out) :: lon_all(:,:)
1377  real(rp), allocatable, intent(out) :: lat_all(:,:)
1378 
1379  character(len=*), intent(in) :: basename_org
1380  character(len=*), intent(in) :: basename_num
1381  logical, intent(in) :: use_file_landwater
1382 
1383  logical, intent(inout) :: serial
1384  logical, intent(inout) :: do_read
1385 
1386  character(len=8) :: file_type = "AUTO" ! "SCALE-RM", "WRFARW", "NAMELIST", "AUTO"
1387  character(len=FILE_HLONG) :: nm_file
1388  logical :: scale_multi_file = .true.
1389  integer :: scale_parent_prc_num_x
1390  integer :: scale_parent_prc_num_y
1391  character(len=FILE_HLONG) :: scale_latlon_catalogue
1392 
1393  namelist / param_mkinit_real_land_netcdf / &
1394  file_type, &
1395  nm_file, &
1396  scale_parent_prc_num_x, &
1397  scale_parent_prc_num_y, &
1398  scale_latlon_catalogue
1399 
1400  character(len=FILE_HLONG) :: basename
1401  character(len=FILE_HLONG) :: fname
1402  integer :: nmfid
1403 
1404  character(len=32) :: items(vars_max)
1405  integer :: nvars
1406  type(vinfo), pointer :: var_info
1407  class(*), pointer :: v
1408 
1409  logical :: error, exist
1410  integer :: n, i
1411  integer :: ierr
1412  !---------------------------------------------------------------------------
1413 
1414  log_info("ParentLandSetupNetCDF",*) 'Real Case/Land Setup'
1415 
1416  file_type = "AUTO"
1417  nm_file = ""
1418  scale_parent_prc_num_x = -1
1419  scale_parent_prc_num_y = -1
1420  scale_latlon_catalogue = ""
1421 
1422  !--- read namelist
1423  rewind(io_fid_conf)
1424  read(io_fid_conf,nml=param_mkinit_real_land_netcdf,iostat=ierr)
1425  if( ierr > 0 ) then
1426  log_error("ParentLandSetupNetCDF",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_LAND_NetCDF. Check!'
1427  call prc_abort
1428  endif
1429  log_nml(param_mkinit_real_land_netcdf)
1430 
1431 
1432  basename = trim(basename_org) // trim(basename_num)
1433 
1434  fid_lnd = -1
1435  if ( do_read ) then
1436  call check_filetype(fid_lnd, file_type, basename, scale_tile_lnd, "ParentLandSetupNetCDF")
1437  end if
1438 
1439  call comm_bcast( file_type )
1440 
1441  if ( file_type == "SCALE-RM" ) then
1442  call comm_bcast( scale_tile_lnd )
1443  if ( scale_tile_lnd ) then
1444  do_read = .true.
1445  serial = .false.
1446  end if
1447  end if
1448 
1449  if ( do_read ) then
1450  vars_land = hash_table()
1451 
1452  select case( file_type )
1453  case ( "SCALE-RM" )
1454  zname = "lz"
1455  xname = "x"
1456  yname = "y"
1457  tname = "time"
1458 
1459  call vars_land%put("lon", vinfo("lon"))
1460  call vars_land%put("lat", vinfo("lat"))
1461  call vars_land%put("lz", vinfo("lz"))
1462 
1463  call vars_land%put("topo", vinfo("topo"))
1464  call vars_land%put("lsmask", vinfo("lsmask"))
1465 
1466  call vars_land%put("LAND_TEMP", vinfo("LAND_TEMP"))
1467  if ( use_file_landwater ) then
1468  call vars_land%put("LAND_WATER", vinfo("LAND_WATER"))
1469  end if
1470 
1471  call vars_land%put("LAND_SFC_TEMP", vinfo("LAND_SFC_TEMP"))
1472 
1473  call vars_land%put("LAND_SFC_ALB_IR_dir", vinfo("LAND_SFC_ALB_IR_dir"))
1474  call vars_land%put("LAND_SFC_ALB_IR_dif", vinfo("LAND_SFC_ALB_IR_dif"))
1475  call vars_land%put("LAND_SFC_ALB_NIR_dir", vinfo("LAND_SFC_ALB_NIR_dir"))
1476  call vars_land%put("LAND_SFC_ALB_NIR_dif", vinfo("LAND_SFC_ALB_NIR_dif"))
1477  call vars_land%put("LAND_SFC_ALB_VIS_dir", vinfo("LAND_SFC_ALB_VIS_dir"))
1478  call vars_land%put("LAND_SFC_ALB_VIS_dif", vinfo("LAND_SFC_ALB_VIS_dif"))
1479 
1480  call vars_land%put("URBAN_SFC_TEMP", vinfo("URBAN_SFC_TEMP"))
1481 
1482  case ( "WRFARW" )
1483  zname = "soil_layers_stag"
1484  xname = "west_east"
1485  yname = "south_north"
1486  tname = "Time"
1487 
1488  call vars_land%put("lon", vinfo("XLONG"))
1489  call vars_land%put("lat", vinfo("XLAT"))
1490  call vars_land%put("lz", vinfo("ZS"))
1491 
1492  call vars_land%put("topo", vinfo("HGT"))
1493  call vars_land%put("lsmask", vinfo("LANDMASK"))
1494 
1495  call vars_land%put("LAND_TEMP", vinfo("TSLB"))
1496  if ( use_file_landwater ) then
1497  call vars_land%put("LAND_WATER", vinfo("SH2O"))
1498  end if
1499 
1500  call vars_land%put("LAND_SFC_TEMP", vinfo("TSK"))
1501 
1502  call vars_land%put("LAND_SFC_ALB_VIS_dir", vinfo("ALBEDO"))
1503  call vars_land%put("LAND_SFC_EMIS_IR_dif", vinfo("EMISS"))
1504 
1505  call vars_land%put("URBAN_SFC_TEMP", vinfo("URBAN_SFC_TEMP"))
1506 
1507 ! call vars_land%put("SNOW_WATER", vinfo("SNOW"))
1508 ! call vars_land%put("SNOW_TEMP", vinfo("TSNAV"))
1509 
1510  case ( "NAMELIST" )
1511  case default
1512  log_error("ParentLandSetupNetCDF",*) 'FILE_TYPE must be "SCALE-RM", "WRFARW", or "AUTO", ', trim(file_type)
1513  call prc_abort
1514  end select
1515 
1516 
1517  !--- read namelist
1518  nmfid = -1
1519  if ( nm_file /= "" ) then
1520  nmfid = io_get_available_fid()
1521  call io_get_fname(fname, nm_file)
1522  open(nmfid, file=fname, form="formatted", status="old", action="read", iostat=ierr)
1523  if ( ierr /= 0 ) then
1524  log_error("ParentLandSetupNetCDF",*) 'namelist file is not found! ', trim(fname)
1525  call prc_abort
1526  end if
1527 
1528  rewind(nmfid)
1529  read(nmfid, nml=netcdf_dims, iostat=ierr)
1530  if( ierr > 0 ) then
1531  log_error("ParentLandSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_DIMS in ', trim(fname), '. Check!'
1532  call prc_abort
1533  end if
1534 
1535  ! items
1536  rewind(nmfid)
1537  nvars = 0
1538  do n = 1, vars_max
1539  read(nmfid, nml=netcdf_item, iostat=ierr)
1540  if( ierr > 0 ) then
1541  log_error("ParentLandSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_ITEM in ', trim(fname), '. Check!'
1542  call prc_abort
1543  else if( ierr < 0 ) then
1544  exit
1545  end if
1546  nvars = nvars + 1
1547  items(nvars) = item
1548  end do
1549  if ( nvars > vars_max ) then
1550  log_error("ParentLandSetupNetCDF",*) "The number of item in the namelist file exceeds the limit! ", nvars
1551  call prc_abort
1552  end if
1553  rewind(nmfid)
1554  do n = 1, nvars
1555  ! set default
1556  if ( vars_land%has_key(items(n)) ) then
1557  item = items(n)
1558  v => vars_land%get(item)
1559  select type(v)
1560  type is (vinfo)
1561  var_info => v
1562  end select
1563  name = var_info%name
1564  fact = var_info%fact
1565  offset = var_info%offset
1566  else
1567  item = items(n)
1568  name = items(n)
1569  fact = 1.0_rp
1570  offset = 0.0_rp
1571  end if
1572  zstg = .false.
1573  xstg = .false.
1574  ystg = .false.
1575  read(nmfid, nml=netcdf_item, iostat=ierr)
1576  if ( ierr /= 0 ) exit
1577  ! set params
1578  call vars_land%put(item, vinfo(name=name, zstg=zstg, xstg=xstg, ystg=ystg, fact=fact, offset=offset))
1579  end do
1580 
1581  else if ( file_type == "NAMELIST" ) then
1582  log_error("ParentLANDSetupNetCDF",*) 'NM_FILE is necessary'
1583  call prc_abort
1584  end if
1585 
1586  end if
1587 
1588  if ( scale_tile_lnd ) then
1589 
1591  scale_domid_lnd, & ! [OUT]
1592  basename, & ! [IN]
1593  scale_parent_prc_num_x, & ! [IN]
1594  scale_parent_prc_num_y, & ! [IN]
1595  scale_latlon_catalogue ) ! [IN]
1596 
1598  scale_domid_lnd, & ! [IN]
1599  lkmax=ldims(1), & ! [OUT]
1600  imaxg=ldims(2), & ! [OUT]
1601  jmaxg=ldims(3), & ! [OUT]
1602  num_tile=nfiles_lnd ) ! [OUT]
1603 
1604  allocate( fids_lnd(nfiles_lnd) )
1605  allocate( tile_id_lnd(nfiles_lnd) )
1606 
1608  scale_domid_lnd, & ! [IN]
1609  tile_id = tile_id_lnd ) ! [OUT]
1610 
1611  call parentlandopennetcdf( basename_org, basename_num )
1612 
1613  else if ( do_read ) then
1614 
1615  call file_get_dimlength( fid_lnd, zname, ldims(1) )
1616  call file_get_dimlength( fid_lnd, xname, ldims(2) )
1617  call file_get_dimlength( fid_lnd, yname, ldims(3) )
1618 
1619  end if
1620 
1621  if ( do_read ) then
1622 
1623  call file_get_dimlength( fid_lnd, tname, timelen, error=error )
1624  if ( error ) timelen = 1
1625 
1626  allocate( lon_all(ldims(2), ldims(3)) )
1627  allocate( lat_all(ldims(2), ldims(3)) )
1628 
1629  call read2d( ldims(2), 1, ldims(2), ldims(3), 1, ldims(3), &
1630  lon_all(:,:), vars_land%get("lon"), &
1631  1, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1632  lon_all(:,:) = lon_all(:,:) * d2r
1633  call read2d( ldims(2), 1, ldims(2), ldims(3), 1, ldims(3), &
1634  lat_all(:,:), vars_land%get("lat"), &
1635  1, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1636  lat_all(:,:) = lat_all(:,:) * d2r
1637 
1638  end if
1639 
1640  first_lnd = .true.
1641 
1642  return
1643  end subroutine parentlandsetupnetcdf
1644 
1645  !-----------------------------------------------------------------------------
1647  subroutine parentlandopennetcdf( &
1648  basename_org, basename_num )
1649  use scale_file, only: &
1650  file_open
1651  implicit none
1652  character(len=*), intent(in) :: basename_org
1653  character(len=*), intent(in) :: basename_num
1654 
1655  character(len=FILE_HLONG) :: basename
1656 
1657  integer :: n
1658  !---------------------------------------------------------------------------
1659 
1660  log_newline
1661  log_info("ParentLandOpenNetCDF",*) 'Real Case/Land Open'
1662 
1663  basename = trim(basename_org) // trim(basename_num)
1664 
1665  if ( scale_tile_lnd ) then
1666 
1667  do n = 1, nfiles_lnd
1668  call file_open( &
1669  basename, & ! [IN]
1670  fids_lnd(n), & ! [OUT]
1671  aggregate=.false., & ! [IN]
1672  allnodes=.false., & ! [IN]
1673  rankid=tile_id_lnd(n) ) ! [IN]
1674  end do
1675 
1676  fid_lnd = fids_lnd(1)
1677 
1678  else
1679 
1680  call file_open(basename, fid_lnd, postfix="")
1681 
1682  end if
1683 
1684  return
1685  end subroutine parentlandopennetcdf
1686 
1687  !-----------------------------------------------------------------------------
1689  subroutine parentlandfinalizenetcdf
1690  implicit none
1691 
1692  !---------------------------------------------------------------------------
1693 
1694  log_info("ParentLandFinalizeNetCDF",*) 'Real Case/Land Finalize'
1695 
1696  if ( allocated(fids_lnd) ) deallocate( fids_lnd )
1697  if ( allocated(tile_id_lnd) ) deallocate( tile_id_lnd )
1698 
1699  call vars_land%destroy()
1700 
1701  scale_domid_lnd = -1
1702  nfiles_lnd = 0
1703  fid_lnd = -1
1704 
1705  return
1706  end subroutine parentlandfinalizenetcdf
1707 
1708  !-----------------------------------------------------------------------------
1709  subroutine parentlandinputnetcdf( &
1710  KA_org, KS_org, KE_org, &
1711  IA_org, IS_org, IE_org, &
1712  JA_org, JS_org, JE_org, &
1713  tg_org, &
1714  strg_org, &
1715  lst_org, &
1716  ust_org, &
1717  albg_org, &
1718  topo_org, &
1719  lmask_org, &
1720  lz_org, &
1721  use_file_landwater, &
1722  ldims, &
1723  it )
1724  use scale_const, only: &
1725  d2r => const_d2r, &
1726  undef => const_undef
1727  implicit none
1728  integer, intent(in) :: ka_org, ks_org, ke_org
1729  integer, intent(in) :: ia_org, is_org, ie_org
1730  integer, intent(in) :: ja_org, js_org, je_org
1731 
1732  real(rp), intent(out) :: tg_org(ka_org,ia_org,ja_org)
1733  real(rp), intent(out) :: strg_org(ka_org,ia_org,ja_org)
1734  real(rp), intent(out) :: lst_org(ia_org,ja_org)
1735  real(rp), intent(out) :: ust_org(ia_org,ja_org)
1736  real(rp), intent(out) :: albg_org(ia_org,ja_org,n_rad_dir,n_rad_rgn)
1737 
1738  real(rp), intent(inout) :: topo_org(ia_org,ja_org)
1739  real(rp), intent(inout) :: lmask_org(ia_org,ja_org)
1740  real(rp), intent(inout) :: lz_org(ka_org)
1741 
1742  logical, intent(in) :: use_file_landwater ! use land water data from files
1743  integer, intent(in) :: ldims(3)
1744  integer, intent(in) :: it
1745 
1746  logical :: exist
1747  integer :: k, i, j
1748  !---------------------------------------------------------------------------
1749 
1750  if ( first_lnd ) then
1751  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1752  topo_org(:,:), vars_land%get("topo"), &
1753  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1754 
1755  call read1d( ka_org, lz_org(:), vars_land%get("lz"), it, fid_lnd )
1756 
1757  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1758  lmask_org(:,:), vars_land%get("lsmask"), &
1759  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1760  end if
1761 
1762  call read3d( ka_org, ks_org, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1763  tg_org(:,:,:), vars_land%get("LAND_TEMP"), &
1764  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1765 
1766 
1767  ! soil liquid water [m3 m-3] (no wrfout-default)
1768  if( use_file_landwater ) then
1769  call read3d( ka_org, ks_org, ke_org, ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1770  strg_org(:,:,:), vars_land%get("LAND_WATER"), &
1771  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1772  endif
1773 
1774  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1775  lst_org(:,:), vars_land%get("LAND_SFC_TEMP"), &
1776  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1777 
1778  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1779  ust_org(:,:), vars_land%get("URBAN_SFC_TEMP"), &
1780  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1781  exist = exist )
1782  if ( .not. exist ) then
1783  !$omp parallel do
1784  do j = 1, ja_org
1785  do i = 1, ia_org
1786  ust_org(i,j) = lst_org(i,j)
1787  end do
1788  end do
1789  end if
1790 
1791  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1792  albg_org(:,:,i_r_direct,i_r_vis), vars_land%get("LAND_SFC_ALB_VIS_dir"), &
1793  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd )
1794  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1795  albg_org(:,:,i_r_diffuse,i_r_vis), vars_land%get("LAND_SFC_ALB_VIS_dif"), &
1796  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1797  exist = exist )
1798  if ( .not. exist ) then
1799  !$omp parallel do
1800  do j = 1, ja_org
1801  do i = 1, ia_org
1802  albg_org(i,j,i_r_diffuse,i_r_vis) = albg_org(i,j,i_r_direct,i_r_vis)
1803  end do
1804  end do
1805  end if
1806  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1807  albg_org(:,:,i_r_direct,i_r_nir), vars_land%get("LAND_SFC_ALB_NIR_dir"), &
1808  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1809  exist = exist )
1810  if ( .not. exist ) then
1811  !$omp parallel do
1812  do j = 1, ja_org
1813  do i = 1, ia_org
1814  albg_org(i,j,i_r_direct,i_r_nir) = albg_org(i,j,i_r_direct,i_r_vis)
1815  end do
1816  end do
1817  end if
1818  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1819  albg_org(:,:,i_r_diffuse,i_r_nir), vars_land%get("LAND_SFC_ALB_NIR_dif"), &
1820  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1821  exist = exist )
1822  if ( .not. exist ) then
1823  !$omp parallel do
1824  do j = 1, ja_org
1825  do i = 1, ia_org
1826  albg_org(i,j,i_r_diffuse,i_r_nir) = albg_org(i,j,i_r_diffuse,i_r_vis)
1827  end do
1828  end do
1829  end if
1830  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1831  albg_org(:,:,i_r_diffuse,i_r_ir), vars_land%get("LAND_SFC_ALB_IR_dif"), &
1832  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1833  exist = exist )
1834  if ( .not. exist ) then
1835  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1836  albg_org(:,:,i_r_diffuse,i_r_ir), vars_land%get("LAND_SFC_EMIS_IR_dif"), &
1837  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1838  exist = exist )
1839  if ( .not. exist ) then
1840  log_error("ParentLandInputNetCDF",*) '"LAND_SFC_ALB_IR_dif" or "LAND_SFC_EMIS_IR_dif" is necessary'
1841  call prc_abort
1842  end if
1843  !$omp parallel do
1844  do j = 1, ja_org
1845  do i = 1, ia_org
1846  albg_org(i,j,i_r_diffuse,i_r_ir) = 1.0_rp - albg_org(i,j,i_r_diffuse,i_r_ir)
1847  end do
1848  end do
1849  end if
1850  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
1851  albg_org(:,:,i_r_direct,i_r_ir), vars_land%get("LAND_SFC_ALB_IR_dir"), &
1852  it, nfiles_lnd, fid_lnd, fids_lnd, scale_tile_lnd, scale_domid_lnd, &
1853  exist = exist )
1854  if ( .not. exist ) then
1855  !$omp parallel do
1856  do j = 1, ja_org
1857  do i = 1, ia_org
1858  albg_org(i,j,i_r_direct,i_r_ir) = albg_org(i,j,i_r_diffuse,i_r_ir)
1859  end do
1860  end do
1861  end if
1862 
1863  first_lnd = .false.
1864 
1865  return
1866  end subroutine parentlandinputnetcdf
1867 
1868  !-----------------------------------------------------------------------------
1870  subroutine parentoceansetupnetcdf( &
1871  odims, &
1872  timelen, &
1873  lon_all, &
1874  lat_all, &
1875  basename_org, &
1876  basename_num, &
1877  serial, &
1878  do_read )
1879  use scale_const, only: &
1880  d2r => const_d2r
1881  use scale_file, only: &
1883  use scale_comm_cartesc, only: &
1884  comm_bcast
1885  use scale_comm_cartesc_nest, only: &
1888 
1889  implicit none
1890 
1891  integer, intent(out) :: odims(2)
1892  integer, intent(out) :: timelen
1893  real(rp), allocatable, intent(out) :: lon_all(:,:)
1894  real(rp), allocatable, intent(out) :: lat_all(:,:)
1895 
1896  character(len=*), intent(in) :: basename_org
1897  character(len=*), intent(in) :: basename_num
1898 
1899  logical, intent(inout) :: serial
1900  logical, intent(inout) :: do_read
1901 
1902  character(len=8) :: file_type = "AUTO" ! "SCALE-RM", "WRFARW", "NAMELIST", "AUTO"
1903  character(len=FILE_HLONG) :: nm_file
1904  logical :: scale_multi_file = .true.
1905  integer :: scale_parent_prc_num_x
1906  integer :: scale_parent_prc_num_y
1907  character(len=FILE_HLONG) :: scale_latlon_catalogue
1908 
1909  namelist / param_mkinit_real_ocean_netcdf / &
1910  file_type, &
1911  nm_file, &
1912  scale_multi_file, &
1913  scale_parent_prc_num_x, &
1914  scale_parent_prc_num_y, &
1915  scale_latlon_catalogue
1916 
1917  character(len=FILE_HLONG) :: basename
1918  character(len=FILE_HLONG) :: fname
1919  integer :: nmfid
1920 
1921  character(len=32) :: items(vars_max)
1922  integer :: nvars
1923  type(vinfo), pointer :: var_info
1924  class(*), pointer :: v
1925 
1926  integer :: n, i
1927  integer :: ierr
1928  logical :: error
1929  !---------------------------------------------------------------------------
1930 
1931  log_info("ParentOceanSetupNetCDF",*) 'Real Case/Ocean Setup'
1932 
1933  file_type = "AUTO"
1934  nm_file = ""
1935  scale_parent_prc_num_x = -1
1936  scale_parent_prc_num_y = -1
1937  scale_latlon_catalogue = ""
1938 
1939  !--- read namelist
1940  rewind(io_fid_conf)
1941  read(io_fid_conf,nml=param_mkinit_real_ocean_netcdf,iostat=ierr)
1942  if( ierr > 0 ) then
1943  log_error("ParentOceanSetupNetCDF",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_OCEAN_NetCDF. Check!'
1944  call prc_abort
1945  endif
1946  log_nml(param_mkinit_real_ocean_netcdf)
1947 
1948 
1949  basename = trim(basename_org) // trim(basename_num)
1950 
1951  fid_ocn = -1
1952  if ( do_read ) then
1953  call check_filetype(fid_ocn, file_type, basename, scale_tile_ocn, "ParentOceanSetupNetCDF")
1954  end if
1955 
1956  call comm_bcast( file_type )
1957 
1958  if ( file_type == "SCALE-RM" ) then
1959  call comm_bcast( scale_tile_ocn )
1960  if ( scale_tile_ocn ) then
1961  do_read = .true.
1962  serial = .false.
1963  end if
1964  end if
1965 
1966 
1967  if ( do_read ) then
1968  vars_ocean = hash_table()
1969 
1970  select case( file_type )
1971  case ( "SCALE-RM" )
1972  xname = "x"
1973  yname = "y"
1974  tname = "time"
1975 
1976  call vars_ocean%put("lon", vinfo("lon"))
1977  call vars_ocean%put("lat", vinfo("lat"))
1978 
1979  call vars_ocean%put("lsmask", vinfo("lsmask"))
1980 
1981  call vars_ocean%put("OCEAN_TEMP", vinfo("OCEAN_TEMP"))
1982 
1983  call vars_ocean%put("OCEAN_SFC_TEMP", vinfo("OCEAN_SFC_TEMP"))
1984  call vars_ocean%put("OCEAN_SFC_Z0M", vinfo("OCEAN_SFC_Z0M"))
1985 
1986  call vars_ocean%put("OCEAN_SFC_ALB_IR_dir", vinfo("OCEAN_SFC_ALB_IR_dir"))
1987  call vars_ocean%put("OCEAN_SFC_ALB_IR_dif", vinfo("OCEAN_SFC_ALB_IR_dif"))
1988  call vars_ocean%put("OCEAN_SFC_ALB_NIR_dir", vinfo("OCEAN_SFC_ALB_NIR_dir"))
1989  call vars_ocean%put("OCEAN_SFC_ALB_NIR_dif", vinfo("OCEAN_SFC_ALB_NIR_dif"))
1990  call vars_ocean%put("OCEAN_SFC_ALB_VIS_dir", vinfo("OCEAN_SFC_ALB_VIS_dir"))
1991  call vars_ocean%put("OCEAN_SFC_ALB_VIS_dif", vinfo("OCEAN_SFC_ALB_VIS_dif"))
1992 
1993  case ( "WRFARW" )
1994  xname = "west_east"
1995  yname = "south_north"
1996  tname = "Time"
1997 
1998  call vars_ocean%put("lon", vinfo("XLONG"))
1999  call vars_ocean%put("lat", vinfo("XLAT"))
2000  call vars_ocean%put("lz", vinfo("ZS"))
2001 
2002  call vars_ocean%put("topo", vinfo("HGT"))
2003  call vars_ocean%put("lsmask", vinfo("LANDMASK"))
2004 
2005  call vars_ocean%put("OCEAN_TEMP", vinfo("OCEAN_TEMP"))
2006 
2007  call vars_ocean%put("OCEAN_SFC_TEMP", vinfo("SST"))
2008  call vars_ocean%put("OCEAN_SFC_Z0M", vinfo("ZNT"))
2009 
2010  call vars_ocean%put("OCEAN_SFC_ALB_VIS_dir", vinfo("ALBEDO"))
2011  call vars_ocean%put("OCEAN_SFC_EMIS_IR_dif", vinfo("EMISS"))
2012 
2013  case ( "NAMELIST" )
2014  case default
2015  log_error("ParentOCEANSetupNetCDF",*) 'FILE_TYPE must be "SCALE-RM", "WRFARW", "NAMELIST", or "AUTO", ', trim(file_type)
2016  call prc_abort
2017  end select
2018 
2019 
2020  !--- read namelist
2021  nmfid = -1
2022  if ( nm_file /= "" ) then
2023  nmfid = io_get_available_fid()
2024  call io_get_fname(fname, nm_file)
2025  open(nmfid, file=fname, form="formatted", status="old", action="read", iostat=ierr)
2026  if ( ierr /= 0 ) then
2027  log_error("ParentOceanSetupNetCDF",*) 'namelist file is not found! ', trim(fname)
2028  call prc_abort
2029  end if
2030 
2031  rewind(nmfid)
2032  read(nmfid, nml=netcdf_dims, iostat=ierr)
2033  if( ierr > 0 ) then
2034  log_error("ParentOceanSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_DIMS in ', trim(fname), '. Check!'
2035  call prc_abort
2036  end if
2037 
2038  ! items
2039  rewind(nmfid)
2040  nvars = 0
2041  do n = 1, vars_max
2042  read(nmfid, nml=netcdf_item, iostat=ierr)
2043  if( ierr > 0 ) then
2044  log_error("ParentLandSetupNetCDF",*) 'Not appropriate names in namelist NetCDF_ITEM in ', trim(fname), '. Check!'
2045  call prc_abort
2046  else if( ierr < 0 ) then
2047  exit
2048  end if
2049  nvars = nvars + 1
2050  items(nvars) = item
2051  end do
2052  if ( nvars > vars_max ) then
2053  log_error("ParentLandSetupNetCDF",*) "The number of item in the namelist file exceeds the limit! ", nvars
2054  call prc_abort
2055  end if
2056  rewind(nmfid)
2057  do n = 1, nvars
2058  ! set default
2059  if ( vars_ocean%has_key(items(n)) ) then
2060  item = items(n)
2061  v => vars_ocean%get(item)
2062  select type(v)
2063  type is (vinfo)
2064  var_info => v
2065  end select
2066  name = var_info%name
2067  fact = var_info%fact
2068  offset = var_info%offset
2069  else
2070  item = items(n)
2071  name = items(n)
2072  fact = 1.0_rp
2073  offset = 0.0_rp
2074  end if
2075  zstg = .false.
2076  xstg = .false.
2077  ystg = .false.
2078  read(nmfid, nml=netcdf_item, iostat=ierr)
2079  if ( ierr /= 0 ) exit
2080  ! set params
2081  call vars_ocean%put(item, vinfo(name=name, zstg=zstg, xstg=xstg, ystg=ystg, fact=fact, offset=offset))
2082  end do
2083 
2084  else if ( file_type == "NAMELIST" ) then
2085  log_error("ParentLANDSetupNetCDF",*) 'NM_FILE is necessary'
2086  call prc_abort
2087  end if
2088 
2089  end if
2090 
2091  if ( scale_tile_ocn ) then
2092 
2094  scale_domid_ocn, & ! [OUT]
2095  basename, & ! [IN]
2096  scale_parent_prc_num_x, & ! [IN]
2097  scale_parent_prc_num_y, & ! [IN]
2098  scale_latlon_catalogue ) ! [IN]
2099 
2101  scale_domid_ocn, & ! [IN]
2102  imaxg=odims(1), & ! [OUT]
2103  jmaxg=odims(2), & ! [OUT]
2104  num_tile=nfiles_ocn ) ! [OUT]
2105 
2106  allocate( fids_ocn(nfiles_ocn) )
2107  allocate( tile_id_ocn(nfiles_ocn) )
2108 
2110  scale_domid_ocn, & ! [IN]
2111  tile_id = tile_id_ocn ) ! [OUT]
2112 
2113  call parentoceanopennetcdf( basename_org, basename_num )
2114 
2115  else if ( do_read ) then
2116 
2117  call file_get_dimlength( fid_ocn, xname, odims(1) )
2118  call file_get_dimlength( fid_ocn, yname, odims(2) )
2119 
2120  end if
2121 
2122  if ( do_read ) then
2123 
2124  call file_get_dimlength( fid_ocn, tname, timelen, error=error )
2125  if ( error ) timelen = 1
2126 
2127  allocate( lon_all(odims(1),odims(2)) )
2128  allocate( lat_all(odims(1),odims(2)) )
2129 
2130  call read2d( odims(1), 1, odims(1), odims(2), 1, odims(2), &
2131  lon_all(:,:), vars_ocean%get("lon"), &
2132  1, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2133  lon_all(:,:) = lon_all(:,:) * d2r
2134  call read2d( odims(1), 1, odims(1), odims(2), 1, odims(2), &
2135  lat_all(:,:), vars_ocean%get("lat"), &
2136  1, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2137  lat_all(:,:) = lat_all(:,:) * d2r
2138 
2139  end if
2140 
2141  first_ocn = .true.
2142 
2143  return
2144  end subroutine parentoceansetupnetcdf
2145 
2146  !-----------------------------------------------------------------------------
2148  subroutine parentoceanopennetcdf( &
2149  basename_org, basename_num )
2150  use scale_file, only: &
2151  file_open
2152  implicit none
2153  character(len=*), intent(in) :: basename_org
2154  character(len=*), intent(in) :: basename_num
2155 
2156  character(len=FILE_HLONG) :: basename
2157 
2158  integer :: n
2159  !---------------------------------------------------------------------------
2160 
2161  log_newline
2162  log_info("ParentOceanOpenNetCDF",*) 'Real Case/Ocean Open'
2163 
2164  basename = trim(basename_org) // trim(basename_num)
2165 
2166  if ( scale_tile_ocn ) then
2167 
2168  do n = 1, nfiles_ocn
2169  call file_open( &
2170  basename, & ! [IN]
2171  fids_ocn(n), & ! [OUT]
2172  aggregate=.false., & ! [IN]
2173  allnodes=.false., & ! [IN]
2174  rankid=tile_id_ocn(n) ) ! [IN]
2175  end do
2176 
2177  fid_ocn = fids_ocn(1)
2178 
2179  else
2180 
2181  call file_open(basename, fid_ocn, postfix="")
2182 
2183  end if
2184 
2185  return
2186  end subroutine parentoceanopennetcdf
2187 
2188  !-----------------------------------------------------------------------------
2190  subroutine parentoceanfinalizenetcdf
2191  implicit none
2192  !---------------------------------------------------------------------------
2193 
2194  log_info("ParentOceanFinalizeNetCDF",*) 'Real Case/Ocean Finalize'
2195 
2196  if ( allocated(fids_ocn) ) deallocate( fids_ocn )
2197  if ( allocated(tile_id_ocn) ) deallocate( tile_id_ocn )
2198 
2199  call vars_ocean%destroy()
2200 
2201  scale_domid_ocn = -1
2202  nfiles_ocn = 0
2203  fid_ocn = -1
2204 
2205  return
2206  end subroutine parentoceanfinalizenetcdf
2207 
2208  !-----------------------------------------------------------------------------
2209  subroutine parentoceaninputnetcdf( &
2210  IA_org, IS_org, IE_org, &
2211  JA_org, JS_org, JE_org, &
2212  tw_org, &
2213  sst_org, &
2214  albw_org, &
2215  z0w_org, &
2216  omask_org, &
2217  odims, &
2218  it )
2219  use scale_const, only: &
2220  d2r => const_d2r, &
2221  undef => const_undef
2222  implicit none
2223  integer, intent(in) :: ia_org, is_org, ie_org
2224  integer, intent(in) :: ja_org, js_org, je_org
2225 
2226  real(rp), intent(out) :: tw_org(ia_org,ja_org)
2227  real(rp), intent(out) :: sst_org(ia_org,ja_org)
2228  real(rp), intent(out) :: albw_org(ia_org,ja_org,n_rad_dir,n_rad_rgn)
2229  real(rp), intent(out) :: z0w_org(ia_org,ja_org)
2230  real(rp), intent(inout) :: omask_org(ia_org,ja_org)
2231 
2232  integer, intent(in) :: odims(2)
2233  integer, intent(in) :: it
2234 
2235  logical :: exist
2236  integer :: i, j
2237  !---------------------------------------------------------------------------
2238 
2239  if ( first_ocn ) then
2240  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2241  omask_org(:,:), vars_ocean%get("lsmask"), &
2242  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2243  end if
2244 
2245  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2246  tw_org(:,:), vars_ocean%get("OCEAN_SFC_TEMP"), &
2247  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2248  sst_org(:,:) = tw_org(:,:)
2249 
2250  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2251  z0w_org(:,:), vars_ocean%get("OCEAN_SFC_Z0M"), &
2252  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2253  exist = exist )
2254  if ( .not. exist ) then
2255  !$omp parallel do
2256  do j = 1, ja_org
2257  do i = 1, ia_org
2258  z0w_org(:,:) = undef
2259  end do
2260  end do
2261  end if
2262 
2263  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2264  albw_org(:,:,i_r_direct,i_r_vis), vars_ocean%get("OCEAN_SFC_ALB_VIS_dir"), &
2265  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn )
2266  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2267  albw_org(:,:,i_r_diffuse,i_r_vis), vars_ocean%get("OCEAN_SFC_ALB_VIS_dif"), &
2268  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2269  exist = exist )
2270  if ( .not. exist ) then
2271  !$omp parallel do
2272  do j = 1, ja_org
2273  do i = 1, ia_org
2274  albw_org(i,j,i_r_diffuse,i_r_vis) = albw_org(i,j,i_r_direct,i_r_vis)
2275  end do
2276  end do
2277  end if
2278  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2279  albw_org(:,:,i_r_direct,i_r_nir), vars_ocean%get("OCEAN_SFC_ALB_NIR_dir"), &
2280  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2281  exist = exist )
2282  if ( .not. exist ) then
2283  !$omp parallel do
2284  do j = 1, ja_org
2285  do i = 1, ia_org
2286  albw_org(i,j,i_r_direct,i_r_nir) = albw_org(i,j,i_r_direct,i_r_vis)
2287  end do
2288  end do
2289  end if
2290  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2291  albw_org(:,:,i_r_diffuse,i_r_nir), vars_ocean%get("OCEAN_SFC_ALB_NIR_dif"), &
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_diffuse,i_r_nir) = albw_org(i,j,i_r_diffuse,i_r_vis)
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_ir), vars_ocean%get("OCEAN_SFC_ALB_IR_dif"), &
2304  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2305  exist = exist )
2306  if ( .not. exist ) then
2307  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2308  albw_org(:,:,i_r_diffuse,i_r_ir), vars_ocean%get("OCEAN_SFC_EMIS_IR_dif"), &
2309  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2310  exist = exist )
2311  if ( .not. exist ) then
2312  log_error("ParentOceanInputNetCDF",*) '"OCEAN_SFC_ALB_IR_dif" or "OCEAN_SFC_EMIS_IR_dif" is necessary'
2313  call prc_abort
2314  end if
2315  !$omp parallel do
2316  do j = 1, ja_org
2317  do i = 1, ia_org
2318  albw_org(i,j,i_r_diffuse,i_r_ir) = 1.0_rp - albw_org(i,j,i_r_diffuse,i_r_ir)
2319  end do
2320  end do
2321  end if
2322  call read2d( ia_org, is_org, ie_org, ja_org, js_org, je_org, &
2323  albw_org(:,:,i_r_direct,i_r_ir), vars_ocean%get("OCEAN_SFC_ALB_IR_dir"), &
2324  it, nfiles_ocn, fid_ocn, fids_ocn, scale_tile_ocn, scale_domid_ocn, &
2325  exist = exist )
2326  if ( .not. exist ) then
2327  !$omp parallel do
2328  do j = 1, ja_org
2329  do i = 1, ia_org
2330  albw_org(i,j,i_r_direct,i_r_ir) = albw_org(i,j,i_r_diffuse,i_r_ir)
2331  end do
2332  end do
2333  end if
2334 
2335  first_ocn = .false.
2336 
2337  return
2338  end subroutine parentoceaninputnetcdf
2339 
2340  ! private
2341 
2342  subroutine check_filetype(fid, FILE_TYPE, basename_org, SCALE_tile, subname)
2343  use scale_file, only: &
2344  file_open, &
2345  file_get_attribute
2346  integer, intent(out) :: fid
2347  character(len=*), intent(inout) :: FILE_TYPE
2348  logical, intent(out) :: SCALE_tile
2349  character(len=*), intent(in) :: basename_org
2350  character(len=*), intent(in) :: subname
2351 
2352  character(len=FILE_HLONG) :: fname
2353  character(len=32) :: att
2354  logical :: exist
2355  integer :: i
2356 
2357  fname = basename_org
2358  inquire(file=fname, exist=exist)
2359  if ( .not. exist ) then
2360  fname = trim(basename_org)//".nc"
2361  inquire(file=fname, exist=exist)
2362  end if
2363  if ( .not. exist ) then
2364  fname = trim(basename_org)//".pe000000.nc"
2365  inquire(file=fname, exist=exist)
2366  end if
2367  if ( .not. exist ) then
2368  log_error(subname,*) "file is not found: ", trim(basename_org)
2369  call prc_abort
2370  end if
2371  call file_open(fname, fid, postfix="", allnodes=.false.)
2372  if ( file_type == "AUTO" ) then
2373  call file_get_attribute( &
2374  fid, "global", "source", &
2375  att, &
2376  existed = exist )
2377  if ( exist .and. att(:8)=="SCALE-RM" ) then
2378  file_type = "SCALE-RM"
2379  log_info(subname,*) 'FILE-TYPE SCALE-RM was detected'
2380  else
2381  call file_get_attribute( &
2382  fid, "global", "TITLE", &
2383  att, &
2384  existed = exist )
2385  if ( exist .and. index(att, "WRF") > 0 ) then
2386  file_type = "WRFARW"
2387  log_info(subname,*) 'FILE-TYPE WRF was detected'
2388  else
2389  file_type = "NAMELIST"
2390  end if
2391  end if
2392  end if
2393 
2394  scale_tile = .false.
2395  if ( file_type == "SCALE-RM" ) then
2396  call file_get_attribute( &
2397  fid, "global", "scale_cartesC_prc_num_x", &
2398  i, &
2399  existed = exist )
2400  if ( exist .and. i > 1 ) then
2401  scale_tile = .true.
2402  log_info(subname,*) 'Multi files was detected'
2403  else
2404  call file_get_attribute( &
2405  fid, "global", "scale_cartesC_prc_num_y", &
2406  i, &
2407  existed = exist )
2408  if ( exist .and. i > 1 ) then
2409  scale_tile = .true.
2410  log_info(subname,*) 'Multi files was detected'
2411  end if
2412  end if
2413  end if
2414 
2415  return
2416  end subroutine check_filetype
2417 
2418  subroutine read3d( &
2419  KA_org, KS_org, KE_org, &
2420  IA_org, IS_org, IE_org, &
2421  JA_org, JS_org, JE_org, &
2422  val, &
2423  var, &
2424  it, &
2425  nfiles, fid, fids, &
2426  scale_tile, scale_domid, &
2427  exist )
2428  use scale_file, only: &
2429  file_get_datainfo, &
2430  file_get_shape, &
2431  file_read
2432  use scale_comm_cartesc_nest, only: &
2434  integer, intent(in) :: KA_org, KS_org, KE_org
2435  integer, intent(in) :: IA_org, IS_org, IE_org
2436  integer, intent(in) :: JA_org, JS_org, JE_org
2437 
2438  real(RP), intent(out), target :: val(KA_org,IA_org,JA_org)
2439 
2440  class(*), pointer, intent(in) :: var
2441  integer, intent(in) :: it
2442  integer, intent(in) :: nfiles
2443  integer, intent(in) :: fid, fids(nfiles)
2444  logical, intent(in) :: scale_tile
2445  integer, intent(in) :: scale_domid
2446 
2447  logical, intent(out), optional :: exist
2448 
2449  real(RP), allocatable :: buf3d(:,:,:)
2450  real(RP), pointer :: work(:,:,:)
2451  real(RP), allocatable, target :: work_t(:,:,:)
2452 
2453  integer :: dims(3)
2454  integer :: tilei, tilej
2455  integer :: kmax
2456  integer :: cxs, cxe, cys, cye
2457  integer :: pxs, pxe, pys, pye
2458  logical :: has_tdim
2459  logical :: transpose
2460  logical :: exist_
2461  integer :: i0, i1, j0, j1
2462  integer :: kst, ist, jst
2463  integer :: it_
2464  integer :: k, i, j, n
2465 
2466  if ( .not. associated(var) ) then
2467  if ( present(exist) ) then
2468  exist = .false.
2469  else
2470  log_error("read3d",*) 'data is not found '
2471  call prc_abort
2472  end if
2473  return
2474  end if
2475 
2476  select type(var)
2477  type is (vinfo)
2478  if ( var%name == "" ) then
2479  if ( present(exist) ) then
2480  exist = .false.
2481  else
2482  log_error("read3d",*) 'data is not found '
2483  call prc_abort
2484  end if
2485  return
2486  end if
2487 
2488  call file_get_datainfo( fid, var%name, has_tdim=has_tdim, existed=exist_ )
2489  if ( .not. exist_ ) then
2490  if ( present(exist) ) then
2491  exist = .false.
2492  return
2493  else
2494  log_error("read3d",*) 'data is not found: ', trim(var%name)
2495  call prc_abort
2496  end if
2497  end if
2498 
2499  if ( has_tdim ) then
2500  it_ = it
2501  else
2502  it_ = 1
2503  end if
2504 
2505  kmax = ke_org - ks_org + 1
2506 
2507  if ( var%zstg ) then
2508  kst = 1
2509  else
2510  kst = 0
2511  end if
2512  if ( var%xstg ) then
2513  ist = 1
2514  else
2515  ist = 0
2516  end if
2517  if ( var%ystg ) then
2518  jst = 1
2519  else
2520  jst = 0
2521  end if
2522 
2523  call file_get_shape( fid, var%name, dims(:) )
2524  transpose = dims(1) .ne. kmax+kst
2525 
2526  if ( scale_tile ) then
2527  if ( var%xstg .or. var%ystg ) then
2528  allocate( work_t(ka_org,ia_org+ist,ja_org+jst) )
2529  work => work_t
2530  else
2531  work => val
2532  end if
2533  do n = 1, nfiles
2535  tilei, tilej, &
2536  cxs, cxe, cys, cye, &
2537  pxs, pxe, pys, pye, &
2538  scale_domid, n, &
2539  xstg = var%xstg, &
2540  ystg = var%ystg )
2541  i0 = max(is_org - cxs, 0)
2542  i1 = max(cxe - ie_org - ist, 0)
2543  j0 = max(js_org - cys, 0)
2544  j1 = max(cye - je_org - jst, 0)
2545  if ( transpose ) then
2546  allocate( buf3d(pxs+i0:pxe-i1,pys+j0:pye-j1,ks_org:ke_org+kst) )
2547  call file_read( fids(n), var%name, buf3d(:,:,:), &
2548  step=it_, start=(/pxs+i0,pys+j0,1/), count=(/pxe-pxs+1-i1-i0,pye-pys+1-j1-j0,kmax+kst/) )
2549  if ( var%zstg ) then
2550  !$omp parallel do
2551  do j = j0, pye-pys-j1
2552  do i = i0, pxe-pxs-i1
2553  do k = ks_org, ke_org
2554  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
2555  end do
2556  end do
2557  end do
2558  else
2559  !$omp parallel do
2560  do j = j0, pye-pys-j1
2561  do i = i0, pxe-pxs-i1
2562  do k = ks_org, ke_org
2563  work(k,cxs+i-is_org+1,cys+j-js_org+1) = buf3d(pxs+i,pys+j,k) * var%fact + var%offset
2564  end do
2565  end do
2566  end do
2567  end if
2568  deallocate( buf3d )
2569  if ( var%xstg .and. cxs==2 .and. is_org==1 ) then ! tentative
2570  !$omp parallel do
2571  do j = j0, pye-pys-j1
2572  do k = ks_org, ke_org
2573  work(k,1,cys+j-js_org+1) = work(k,2,cys+j-js_org+1)
2574  end do
2575  end do
2576  end if
2577  if ( var%ystg .and. cys==2 .and. js_org==1 ) then ! tentative
2578  !$omp parallel do
2579  do i = i0, pxe-pxs-i1
2580  do k = ks_org, ke_org
2581  work(k,cxs+i-is_org+1,1) = work(k,cxs+i-is_org+1,2)
2582  end do
2583  end do
2584  end if
2585  else
2586  allocate( buf3d(ks_org:ke_org+kst,pxs+i0:pxe-i1,pys+j0:pye-j1) )
2587  call file_read( fids(n), var%name, buf3d(:,:,:), &
2588  step=it_, start=(/1,pxs+i0,pys+j0/), count=(/kmax+kst,pxe-pxs+1-i1-i0,pye-pys+1-j1-j0/) )
2589  if ( var%zstg ) then
2590  !$omp parallel do
2591  do j = j0, pye-pys-j1
2592  do i = i0, pxe-pxs-i1
2593  do k = ks_org, ke_org
2594  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
2595  end do
2596  end do
2597  end do
2598  else
2599  !$omp parallel do
2600  do j = j0, pye-pys-j1
2601  do i = i0, pxe-pxs-i1
2602  do k = ks_org, ke_org
2603  work(k,cxs+i-is_org+1,cys+j-js_org+1) = buf3d(k,pxs+i,pys+j) * var%fact + var%offset
2604  end do
2605  end do
2606  end do
2607  end if
2608  deallocate( buf3d )
2609  end if
2610  end do
2611  if ( var%xstg ) then
2612  !$omp parallel do collapse(2)
2613  do j = 1, ja_org
2614  do i = 1, ia_org
2615  do k = ks_org, ke_org
2616  val(k,i,j) = ( work(k,i,j) + work(k,i+1,j) ) * 0.5_rp
2617  end do
2618  end do
2619  end do
2620  else if ( var%ystg ) then
2621  !$omp parallel do collapse(2)
2622  do j = 1, ja_org
2623  do i = 1, ia_org
2624  do k = ks_org, ke_org
2625  val(k,i,j) = ( work(k,i,j) + work(k,i,j+1) ) * 0.5_rp
2626  end do
2627  end do
2628  end do
2629  end if
2630  if ( var%xstg .or. var%ystg ) then
2631  deallocate( work_t )
2632  end if
2633  nullify( work )
2634  else
2635  if ( transpose ) then
2636  allocate( buf3d(is_org:ie_org+ist,js_org:je_org+jst,ks_org:ke_org+kst) )
2637  call file_read( &
2638  fid, var%name, &
2639  buf3d(:,:,:), &
2640  step=it_, &
2641  start=(/is_org,js_org,1/), &
2642  count=(/ia_org+ist,ja_org+jst,kmax+kst/))
2643  if ( var%zstg ) then
2644  !$omp parallel do collapse(2)
2645  do j = 1, ja_org
2646  do i = 1, ia_org
2647  do k = ks_org, ke_org
2648  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
2649  end do
2650  end do
2651  end do
2652  else if ( var%xstg ) then
2653  !$omp parallel do collapse(2)
2654  do j = 1, ja_org
2655  do i = 1, ia_org
2656  do k = ks_org, ke_org
2657  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
2658  end do
2659  end do
2660  end do
2661  else if ( var%ystg ) then
2662  !$omp parallel do collapse(2)
2663  do j = 1, ja_org
2664  do i = 1, ia_org
2665  do k = ks_org, ke_org
2666  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
2667  end do
2668  end do
2669  end do
2670  else
2671  !$omp parallel do collapse(2)
2672  do j = 1, ja_org
2673  do i = 1, ia_org
2674  do k = ks_org, ke_org
2675  val(k,i,j) = buf3d(i+is_org-1,j+js_org-1,k) * var%fact + var%offset
2676  end do
2677  end do
2678  end do
2679  end if
2680  deallocate( buf3d )
2681  else
2682  allocate( buf3d(ks_org:ke_org+kst,is_org:ie_org+ist,js_org:je_org+jst) )
2683  call file_read( &
2684  fid, var%name, &
2685  buf3d(:,:,:), &
2686  step=it_, &
2687  start=(/1,is_org,js_org/), &
2688  count=(/kmax+kst,ia_org+ist,ja_org+jst/) )
2689  if ( var%zstg ) then
2690  !$omp parallel do collapse(2)
2691  do j = 1, ja_org
2692  do i = 1, ia_org
2693  do k = ks_org, ke_org
2694  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
2695  end do
2696  end do
2697  end do
2698  else if ( var%xstg ) then
2699  !$omp parallel do collapse(2)
2700  do j = 1, ja_org
2701  do i = 1, ia_org
2702  do k = ks_org, ke_org
2703  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
2704  end do
2705  end do
2706  end do
2707  else if ( var%ystg ) then
2708  !$omp parallel do collapse(2)
2709  do j = 1, ja_org
2710  do i = 1, ia_org
2711  do k = ks_org, ke_org
2712  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
2713  end do
2714  end do
2715  end do
2716  else
2717  !$omp parallel do collapse(2)
2718  do j = 1, ja_org
2719  do i = 1, ia_org
2720  do k = ks_org, ke_org
2721  val(k,i,j) = buf3d(k,i+is_org-1,j+js_org-1) * var%fact + var%offset
2722  end do
2723  end do
2724  end do
2725  end if
2726  deallocate( buf3d )
2727  end if
2728  end if
2729 
2730  if ( present(exist) ) exist = .true.
2731  end select
2732 
2733  return
2734  end subroutine read3d
2735 
2736  subroutine read2d( &
2737  IA_org, IS_org, IE_org, &
2738  JA_org, JS_org, JE_org, &
2739  val, &
2740  var, &
2741  it, &
2742  nfiles, fid, fids, &
2743  scale_tile, scale_domid, &
2744  exist )
2745  use scale_file, only: &
2746  file_get_datainfo, &
2747  file_get_shape, &
2748  file_read
2749  use scale_comm_cartesc_nest, only: &
2751  integer, intent(in) :: IA_org, IS_org, IE_org
2752  integer, intent(in) :: JA_org, JS_org, JE_org
2753 
2754  real(RP), intent(out), target :: val(IA_org,JA_org)
2755 
2756  class(*), pointer, intent(in) :: var
2757  integer, intent(in) :: it
2758  integer, intent(in) :: nfiles
2759  integer, intent(in) :: fid, fids(nfiles)
2760  logical, intent(in) :: scale_tile
2761  integer, intent(in) :: scale_domid
2762 
2763  logical, intent(out), optional :: exist
2764 
2765  real(RP), allocatable :: buf2d(:,:)
2766  real(RP), pointer :: work(:,:)
2767  real(RP), allocatable, target :: work_t(:,:)
2768 
2769  integer :: tilei, tilej
2770  integer :: cxs, cxe, cys, cye
2771  integer :: pxs, pxe, pys, pye
2772 
2773  integer :: dims(2)
2774 
2775  integer :: i0, i1, j0, j1
2776  integer :: ist, jst
2777 
2778  logical :: has_tdim
2779  logical :: exist_
2780  integer :: it_
2781  integer :: n, i, j
2782 
2783  if ( .not. associated(var) ) then
2784  if ( present(exist) ) then
2785  exist = .false.
2786  else
2787  log_error("read3d",*) 'data is not found '
2788  call prc_abort
2789  end if
2790  return
2791  end if
2792 
2793  select type(var)
2794  type is (vinfo)
2795  if ( var%name == "" ) then
2796  if ( present(exist) ) then
2797  exist = .false.
2798  else
2799  log_error("read2d",*) 'data is not found '
2800  call prc_abort
2801  end if
2802  return
2803  end if
2804 
2805  call file_get_datainfo( fid, var%name, has_tdim=has_tdim, existed=exist_ )
2806  if ( .not. exist_ ) then
2807  if ( present(exist) ) then
2808  exist = .false.
2809  return
2810  else
2811  log_error("read2d",*) 'data is not found: ', trim(var%name)
2812  call prc_abort
2813  end if
2814  end if
2815 
2816  if ( has_tdim ) then
2817  it_ = it
2818  else
2819  it_ = 1
2820  end if
2821 
2822  if ( var%xstg ) then
2823  ist = 1
2824  else
2825  ist = 0
2826  end if
2827  if ( var%ystg ) then
2828  jst = 1
2829  else
2830  jst = 0
2831  end if
2832 
2833  if ( scale_domid > 0 ) then
2834  if ( var%xstg .or. var%ystg ) then
2835  allocate( work_t(ia_org+ist,ja_org+jst) )
2836  work => work_t
2837  else
2838  work => val
2839  end if
2840  do n = 1, nfiles
2842  tilei, tilej, &
2843  cxs, cxe, cys, cye, &
2844  pxs, pxe, pys, pye, &
2845  scale_domid, n, &
2846  xstg = var%xstg, &
2847  ystg = var%ystg )
2848  i0 = max(is_org - cxs, 0)
2849  i1 = max(cxe - ie_org - ist, 0)
2850  j0 = max(js_org - cys, 0)
2851  j1 = max(cye - je_org - jst, 0)
2852  allocate( buf2d(pxs+i0:pxe-i1,pys+j0:pye-j1) )
2853  call file_read( fids(n), var%name, buf2d(:,:), &
2854  step=it_, start=(/pxs+i0,pys+j0/), count=(/pxe-pxs+1-i1-i0,pye-pys+1-j1-j0/) )
2855  !$omp parallel do
2856  do j = j0, pye-pys-j1
2857  do i = i0, pxe-pxs-i1
2858  work(cxs+i-is_org+1,cys+j-js_org+1) = buf2d(pxs+i,pys+j) * var%fact + var%offset
2859  end do
2860  end do
2861  deallocate( buf2d )
2862  if ( var%xstg .and. cxs==2 .and. is_org==1 ) then ! tentative
2863  !$omp parallel do
2864  do j = j0, pye-pys-j1
2865  work(1,cys+j-js_org+1) = work(2,cys+j-js_org+1)
2866  end do
2867  end if
2868  if ( var%ystg .and. cys==2 .and. js_org==1 ) then ! tentative
2869  !$omp parallel do
2870  do i = i0, pxe-pxs-i1
2871  work(cxs+i-is_org+1,1) = work(cxs+i-is_org+1,2)
2872  end do
2873  end if
2874  end do
2875  if ( var%xstg ) then
2876  !$omp parallel do
2877  do j = 1, ja_org
2878  do i = 1, ia_org
2879  val(i,j) = ( work(i,j) + work(i+1,j) ) * 0.5_rp
2880  end do
2881  end do
2882  else if ( var%ystg ) then
2883  !$omp parallel do
2884  do j = 1, ja_org
2885  do i = 1, ia_org
2886  val(i,j) = ( work(i,j) + work(i,j+1) ) * 0.5_rp
2887  end do
2888  end do
2889  end if
2890  if ( var%xstg .or. var%ystg ) then
2891  deallocate( work_t )
2892  end if
2893  nullify( work )
2894  else
2895  if ( var%xstg .or. var%ystg ) then
2896  allocate( work_t(is_org:ie_org+ist,js_org:je_org+jst) )
2897  work => work_t
2898  else
2899  work => val
2900  end if
2901  call file_read( &
2902  fid, var%name, &
2903  work(:,:), &
2904  step=it_, &
2905  start=(/is_org,js_org/), &
2906  count=(/ia_org+ist,ja_org+jst/) )
2907  if ( var%xstg ) then
2908  !$omp parallel do
2909  do j = 1, ja_org
2910  do i = 1, ia_org
2911  val(i,j) = ( work(i,j) + work(i+1,j) ) * 0.5_rp * var%fact + var%offset
2912  end do
2913  end do
2914  else if ( var%ystg ) then
2915  !$omp parallel do
2916  do j = 1, ja_org
2917  do i = 1, ia_org
2918  val(i,j) = ( work(i,j) + work(i,j+1) ) * 0.5_rp * var%fact + var%offset
2919  end do
2920  end do
2921  else if ( var%fact .ne. 1.0_rp .or. var%offset .ne. 0.0_rp ) then
2922  !$omp parallel do
2923  do j = 1, ja_org
2924  do i = 1, ia_org
2925  val(i,j) = val(i,j) * var%fact + var%offset
2926  end do
2927  end do
2928  end if
2929  if ( var%xstg .or. var%ystg ) then
2930  deallocate( work_t )
2931  end if
2932  nullify( work )
2933  end if
2934 
2935  if ( present(exist) ) exist = .true.
2936  end select
2937 
2938  return
2939  end subroutine read2d
2940 
2941  subroutine read1d( &
2942  KA_org, &
2943  val, &
2944  var, &
2945  it, &
2946  fid, &
2947  exist )
2948  use scale_file, only: &
2949  file_get_datainfo, &
2950  file_read
2951  integer, intent(in) :: KA_org
2952 
2953  real(RP), intent(out) :: val(KA_org)
2954 
2955  class(*), pointer, intent(in) :: var
2956  integer, intent(in) :: it
2957  integer, intent(in) :: fid
2958 
2959  logical, intent(out), optional :: exist
2960 
2961  logical :: has_tdim
2962  logical :: exist_
2963 
2964  integer :: it_
2965 
2966  if ( .not. associated(var) ) then
2967  if ( present(exist) ) then
2968  exist = .false.
2969  else
2970  log_error("read3d",*) 'data is not found '
2971  call prc_abort
2972  end if
2973  return
2974  end if
2975 
2976  select type(var)
2977  type is ( vinfo )
2978  if ( var%name == "" ) then
2979  if ( present(exist) ) then
2980  exist = .false.
2981  else
2982  log_error("read1d",*) 'data is not found '
2983  call prc_abort
2984  end if
2985  return
2986  end if
2987 
2988  call file_get_datainfo( fid, var%name, has_tdim=has_tdim, existed=exist_ )
2989  if ( .not. exist_ ) then
2990  if ( present(exist) ) then
2991  exist = .false.
2992  return
2993  else
2994  log_error("read1d",*) 'data is not found: ', trim(var%name)
2995  call prc_abort
2996  end if
2997  end if
2998 
2999  if ( has_tdim ) then
3000  it_ = it
3001  else
3002  it_ = 1
3003  end if
3004 
3005  call file_read( fid, var%name, val(:), step=it_ )
3006  if ( var%fact .ne. 1.0_rp .or. var%offset .ne. 0.0_rp ) then
3007  val(:) = val(:) * var%fact + var%offset
3008  end if
3009 
3010  if ( present(exist) ) exist = .true.
3011  end select
3012 
3013  return
3014  end subroutine read1d
3015 
3016 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:2150
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:1649
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:349
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:2191
mod_realinput_netcdf::parentlandfinalizenetcdf
subroutine, public parentlandfinalizenetcdf
Land Finalize.
Definition: mod_realinput_netcdf.F90:1690
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:535
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:1724
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
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_mptype, sfc_diagnoses, update_coord, dims, it)
Definition: mod_realinput_netcdf.F90:727
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:90
scale_file::file_get_dimlength
subroutine, public file_get_dimlength(fid, dimname, len, error)
get length of dimension
Definition: scale_file.F90:632
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:2948
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:2343
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:2428
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:1879
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::parentatmossetupnetcdf
subroutine, public parentatmossetupnetcdf(dims, timelen, mixing_ratio, update_coord, mapping_info, qtrc_flag, lon_all, lat_all, basename_org, basename_num, same_mptype, pt_dry, serial, do_read)
Atmos Setup.
Definition: mod_realinput_netcdf.F90:150
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:1362
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:2745
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:2219
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:91
mod_realinput_netcdf::parentatmosopennetcdf
subroutine, public parentatmosopennetcdf(basename_org, basename_num)
Atmos Open.
Definition: mod_realinput_netcdf.F90:642