SCALE-RM
mod_realinput_grads.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_stdio
20  use scale_tracer
21  use scale_process, only: &
22  myrank => prc_myrank, &
24  !-----------------------------------------------------------------------------
25  implicit none
26  private
27  !-----------------------------------------------------------------------------
28  !
29  !++ Public procedure
30  !
31  public :: parentatomsetupgrads
32  public :: parentatomopengrads
33  public :: parentatominputgrads
34  public :: parentlandsetupgrads
35  public :: parentlandinputgrads
36  public :: parentoceansetupgrads
37  public :: parentoceanopengrads
38  public :: parentoceaninputgrads
39 
40  !-----------------------------------------------------------------------------
41  !
42  !++ Public parameters & variables
43  !
44  !-----------------------------------------------------------------------------
45  !
46  !++ Private procedure
47  !
48  !-----------------------------------------------------------------------------
49  !
50  !++ Private parameters & variables
51  !
52  integer, parameter :: grads_vars_limit = 1000
53  integer, parameter :: num_item_list = 17
54  integer, parameter :: num_item_list_atom = 17
55  integer, parameter :: num_item_list_land = 11
56  integer, parameter :: num_item_list_ocean = 9
57  logical :: data_available(num_item_list_atom,3) ! 1:atom, 2:land, 3:ocean
58  character(len=H_SHORT) :: item_list_atom (num_item_list_atom)
59  character(len=H_SHORT) :: item_list_land (num_item_list_land)
60  character(len=H_SHORT) :: item_list_ocean(num_item_list_ocean)
61  data item_list_atom /'lon','lat','plev','U','V','T','HGT','QV','RH','MSLP','PSFC','U10','V10','T2','Q2','RH2','TOPO' /
62  data item_list_land /'lsmask','lon','lat','lon_sfc','lat_sfc','llev', &
63  'STEMP','SMOISVC','SMOISDS','SKINT','TOPO' /
64  data item_list_ocean /'lsmask','lon','lat','lon_sfc','lat_sfc','lon_sst','lat_sst','SKINT','SST'/
65 
66  integer, parameter :: ia_lon = 1
67  integer, parameter :: ia_lat = 2
68  integer, parameter :: ia_p = 3 ! Pressure (Pa)
69  integer, parameter :: ia_u = 4
70  integer, parameter :: ia_v = 5
71  integer, parameter :: ia_t = 6
72  integer, parameter :: ia_hgt = 7 ! Geopotential height (m)
73  integer, parameter :: ia_qv = 8
74  integer, parameter :: ia_rh = 9 ! Percentage (%)
75  integer, parameter :: ia_slp = 10 ! Sea level pressure (Pa)
76  integer, parameter :: ia_ps = 11 ! Surface pressure (Pa)
77  integer, parameter :: ia_u10 = 12
78  integer, parameter :: ia_v10 = 13
79  integer, parameter :: ia_t2 = 14
80  integer, parameter :: ia_q2 = 15
81  integer, parameter :: ia_rh2 = 16 ! Percentage (%)
82  integer, parameter :: ia_topo = 17
83 
84  integer, parameter :: il_lsmask = 1
85  integer, parameter :: il_lon = 2
86  integer, parameter :: il_lat = 3
87  integer, parameter :: il_lon_sfc = 4
88  integer, parameter :: il_lat_sfc = 5
89  integer, parameter :: il_lz = 6 ! Level(depth) of stemp & smois (m)
90  integer, parameter :: il_stemp = 7
91  integer, parameter :: il_smoisvc = 8 ! soil moisture (vormetric water content)
92  integer, parameter :: il_smoisds = 9 ! soil moisture (degree of saturation)
93  integer, parameter :: il_skint = 10
94  integer, parameter :: il_topo = 11
95 
96  integer, parameter :: io_lsmask = 1
97  integer, parameter :: io_lon = 2
98  integer, parameter :: io_lat = 3
99  integer, parameter :: io_lon_sfc = 4
100  integer, parameter :: io_lat_sfc = 5
101  integer, parameter :: io_lon_sst = 6
102  integer, parameter :: io_lat_sst = 7
103  integer, parameter :: io_skint = 8
104  integer, parameter :: io_sst = 9
105 
106 
107  integer, parameter :: lvars_limit = 1000 ! limit of values for levels data
108  real(RP), parameter :: large_number_one = 9.999e+15_rp
109 
110 
111  character(len=H_SHORT) :: upper_qv_type = "ZERO"
114 
115  character(len=H_SHORT) :: grads_item (num_item_list,3)
116  character(len=H_LONG) :: grads_dtype (num_item_list,3)
117  character(len=H_LONG) :: grads_fname (num_item_list,3)
118  character(len=H_SHORT) :: grads_fendian (num_item_list,3)
119  character(len=H_SHORT) :: grads_yrev (num_item_list,3)
120  real(RP) :: grads_swpoint (num_item_list,3)
121  real(RP) :: grads_dd (num_item_list,3)
122  integer :: grads_lnum (num_item_list,3)
123  real(RP) :: grads_lvars (lvars_limit,num_item_list,3)
124  integer :: grads_startrec(num_item_list,3)
125  integer :: grads_totalrec(num_item_list,3)
126  integer :: grads_knum (num_item_list,3)
127  real(SP) :: grads_missval (num_item_list,3)
128 
129  real(SP), allocatable :: gdata2d(:,:)
130  real(SP), allocatable :: gdata3d(:,:,:)
131  real(SP), allocatable :: gland2d(:,:)
132  real(SP), allocatable :: gland3d(:,:,:)
133  real(SP), allocatable :: gsst2d (:,:)
134 
135  integer :: io_fid_grads_nml = -1
136  integer :: io_fid_grads_data = -1
137 
138 
139  ! atmos data
140  integer :: outer_nx = -1
141  integer :: outer_ny = -1
142  integer :: outer_nz = -1 ! number of atmos layers
143  integer :: outer_nl = -1 ! number of land layers
144  ! surface data
145  integer :: outer_nx_sfc = -1
146  integer :: outer_ny_sfc = -1
147  ! sst data
148  integer :: outer_nx_sst = -1
149  integer :: outer_ny_sst = -1
150 
151  namelist / nml_grads_grid / &
152  outer_nx, &
153  outer_ny, &
154  outer_nz, &
155  outer_nl, &
156  outer_nx_sfc, &
157  outer_ny_sfc, &
158  outer_nx_sst, &
159  outer_ny_sst
160 
161  character(len=H_SHORT) :: item ! up to 16 characters
162  integer :: knum ! optional: vertical level
163  character(len=H_SHORT) :: dtype ! 'linear','levels','map'
164  character(len=H_LONG) :: fname ! head of file name
165  real(RP) :: swpoint ! start point (south-west point) for "linear"
166  real(RP) :: dd ! dlon,dlat for "linear"
167  integer :: lnum ! number of data
168  real(RP) :: lvars(lvars_limit) = large_number_one ! values for "levels"
169  integer :: startrec ! record position
170  integer :: totalrec ! total record number per one time
171  real(SP) :: missval ! missing value
172  character(len=H_SHORT) :: fendian='big' ! option for "map"
173  character(len=H_SHORT) :: yrev='off' ! option for "map", if yrev=on, order of data is NW to SE.
174 
175 
176  !-----------------------------------------------------------------------------
177 contains
178  !-----------------------------------------------------------------------------
180  subroutine parentatomsetupgrads( &
181  dims, & ! (out)
182  basename ) ! (in)
183  implicit none
184 
185  integer, intent(out) :: dims(6)
186  character(len=H_LONG), intent(in) :: basename
187 
188 
189  namelist / param_mkinit_real_grads / &
190  upper_qv_type
191 
192  integer :: ielem
193  integer :: k, n
194 
195  integer :: ierr
196  !---------------------------------------------------------------------------
197 
198  if( io_l ) write(io_fid_log,*) '+++ Real Case/Atom Input File Type: GrADS format'
199 
200  !--- read namelist
201  rewind(io_fid_conf)
202  read(io_fid_conf,nml=param_mkinit_real_grads,iostat=ierr)
203 
204  if( ierr > 0 ) then
205  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_GrADS. Check!'
206  call prc_mpistop
207  endif
208  if( io_lnml ) write(io_fid_log,nml=param_mkinit_real_grads)
209 
210 
211  if ( len_trim(basename) == 0 ) then
212  if( io_l ) write(io_fid_log,*) &
213  'xxx "BASENAME_ORG" is not specified in "PARAM_MKINIT_REAL_ATMOS"!',trim(basename)
214  call prc_mpistop
215  endif
216 
217  !--- read namelist
218  io_fid_grads_nml = io_get_available_fid()
219  open( io_fid_grads_nml, &
220  file = trim(basename), &
221  form = 'formatted', &
222  status = 'old', &
223  action = 'read', &
224  iostat = ierr )
225  if ( ierr /= 0 ) then
226  if( io_l ) write(io_fid_log,*) 'xxx Input file is not found! ', trim(basename)
227  call prc_mpistop
228  endif
229 
230  read(io_fid_grads_nml,nml=nml_grads_grid,iostat=ierr)
231  if( ierr /= 0 ) then !--- missing or fatal error
232  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in nml_grads_grid in ', trim(basename),'. Check!'
233  call prc_mpistop
234  endif
235  if( io_lnml ) write(io_fid_log,nml=nml_grads_grid)
236 
237  ! full level
238  dims(1) = outer_nz ! bottom_top
239  dims(2) = outer_nx ! west_east
240  dims(3) = outer_ny ! south_north
241  ! half level
242  dims(4) = outer_nz ! bottom_top_stag
243  dims(5) = outer_nx ! west_east for 2dim data
244  dims(6) = outer_ny ! south_north for 2dim data
245 
246  allocate( gdata2d( dims(2), dims(3) ) )
247  allocate( gdata3d( dims(2), dims(3), dims(1) ) )
248 
249  call read_namelist( &
250  grads_item(:,1), & ! (out)
251  grads_fname(:,1), & ! (out)
252  grads_dtype(:,1), & ! (out)
253  grads_swpoint(:,1), & ! (out)
254  grads_dd(:,1), & ! (out)
255  grads_lnum(:,1), & ! (out)
256  grads_lvars(:,:,1), & ! (out)
257  grads_startrec(:,1), & ! (out)
258  grads_totalrec(:,1), & ! (out)
259  grads_knum(:,1), & ! (out)
260  grads_yrev(:,1), & ! (out)
261  grads_fendian(:,1), & ! (out)
262  grads_missval(:,1), & ! (out)
263  data_available(:,1), & ! (out)
264  item_list_atom, & ! (in)
265  num_item_list_atom, & ! (in)
266  basename, & ! (in)
267  io_fid_grads_nml ) ! (in)
268 
269  close( io_fid_grads_nml )
270 
271  do ielem = 1, num_item_list_atom
272  item = item_list_atom(ielem)
273  !--- check data
274  select case(trim(item))
275  case('QV')
276  if (.not. data_available(ia_qv,1)) then
277  if (.not.data_available(ia_rh,1)) then
278  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : QV and RH'
279  call prc_mpistop
280  else ! will read RH
281  cycle
282  endif
283  endif
284  case('RH')
285  if (.not. data_available(ia_qv,1))then
286  if(data_available(ia_rh,1)) then
287  if ((.not. data_available(ia_t,1)).or.(.not. data_available(ia_p,1))) then
288  if( io_l ) write(io_fid_log,*) 'xxx Temperature and pressure are required to convert from RH to QV ! '
289  call prc_mpistop
290  else
291  cycle ! read RH and estimate QV
292  endif
293  else
294  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : QV and RH'
295  call prc_mpistop
296  endif
297  endif
298  case('MSLP','PSFC','U10','V10','T2')
299  if (.not. data_available(ielem,1)) then
300  if (io_l) write(io_fid_log,*) 'warning: ',trim(item),' is not found & will be estimated.'
301  cycle
302  endif
303  case('Q2')
304  if ( .not. data_available(ia_q2,1) ) then
305  if (io_l) write(io_fid_log,*) 'warning: Q2 is not found & will be estimated.'
306  cycle
307  endif
308  case('RH2')
309  if ( data_available(ia_q2,1) ) then
310  cycle
311  else
312  if ( data_available(ia_rh2,1) ) then
313  if ((.not. data_available(ia_t2,1)).or.(.not. data_available(ia_ps,1))) then
314  if (io_l) write(io_fid_log,*) 'warning: T2 and PSFC are required to convert from RH2 to Q2 !'
315  if (io_l) write(io_fid_log,*) ' Q2 will be copied from data at above level.'
316  data_available(ia_rh2,1) = .false.
317  cycle
318  endif
319  else
320  if (io_l) write(io_fid_log,*) 'warning: Q2 and RH2 are not found, Q2 will be estimated.'
321  cycle
322  endif
323  endif
324  case('TOPO')
325  if ( .not. data_available(ielem,1) ) then
326  if (io_l) write(io_fid_log,*) 'warning: ',trim(item),' is not found & not used.'
327  cycle
328  endif
329  case default ! lon, lat, plev, U, V, T, HGT
330  if ( .not. data_available(ielem,1) ) then
331  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : ',trim(item_list_atom(ielem))
332  call prc_mpistop
333  endif
334  end select
335 
336  end do
337 
338  return
339  end subroutine parentatomsetupgrads
340 
341  !-----------------------------------------------------------------------------
342  subroutine parentatomopengrads
343  implicit none
344 
345  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[AtomOpenGrADS]'
346 
347  return
348  end subroutine parentatomopengrads
349 
350  !-----------------------------------------------------------------------------
351  subroutine parentatominputgrads( &
352  velz_org, &
353  velx_org, &
354  vely_org, &
355  pres_org, &
356  temp_org, &
357  qtrc_org, &
358  lon_org, &
359  lat_org, &
360  cz_org, &
361  basename_num, &
362  dims, &
363  nt )
364  use scale_const, only: &
365  d2r => const_d2r, &
366  eps => const_eps, &
367  epsvap => const_epsvap, &
368  grav => const_grav, &
369  laps => const_laps, &
370  p00 => const_pre00, &
371  rdry => const_rdry, &
372  cpdry => const_cpdry
373  use scale_atmos_saturation, only: &
374  psat => atmos_saturation_psat_liq
375  implicit none
376 
377 
378  real(RP), intent(out) :: velz_org(:,:,:)
379  real(RP), intent(out) :: velx_org(:,:,:)
380  real(RP), intent(out) :: vely_org(:,:,:)
381  real(RP), intent(out) :: pres_org(:,:,:)
382  real(RP), intent(out) :: temp_org(:,:,:)
383  real(RP), intent(out) :: qtrc_org(:,:,:,:)
384  real(RP), intent(out) :: lon_org(:,:)
385  real(RP), intent(out) :: lat_org(:,:)
386  real(RP), intent(out) :: cz_org(:,:,:)
387  character(len=*), intent(in) :: basename_num
388  integer, intent(in) :: dims(6)
389  integer, intent(in) :: nt
390 
391  real(RP) :: rhprs_org(dims(1)+2,dims(2),dims(3))
392  real(RP) :: pott(dims(2),dims(3))
393 
394  real(RP) :: RovCP
395  real(RP) :: CPovR
396 
397  ! data
398  character(len=H_LONG) :: gfile
399 
400  integer :: QA_outer = 1
401  real(RP) :: p_sat, qm, rhsfc
402 
403  integer :: i, j, k, ielem
404 
405  !---------------------------------------------------------------------------
406 
407  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[AtomInputGrADS]'
408 
409  qtrc_org = 0.0_rp
410 
411  !--- read grads data
412  loop_inputatomgrads : do ielem = 1, num_item_list_atom
413 
414  if ( .not. data_available(ielem,1) ) cycle
415 
416  item = grads_item(ielem,1)
417  dtype = grads_dtype(ielem,1)
418  fname = grads_fname(ielem,1)
419  lnum = grads_lnum(ielem,1)
420 
421  if ( dims(1) < grads_knum(ielem,1) ) then
422  write(*,*) 'xxx "knum" must be less than or equal to outer_nz. knum:',knum,'> outer_nz:',dims(1),trim(item)
423  call prc_mpistop
424  else if ( grads_knum(ielem,1) > 0 )then
425  knum = grads_knum(ielem,1) ! not missing
426  else
427  knum = dims(1)
428  endif
429 
430  select case (trim(dtype))
431  case("linear")
432  swpoint = grads_swpoint(ielem,1)
433  dd = grads_dd(ielem,1)
434  if( (abs(swpoint-large_number_one)<eps).or.(abs(dd-large_number_one)<eps) )then
435  write(*,*) 'xxx "swpoint" is required in grads namelist! ',swpoint
436  write(*,*) 'xxx "dd" is required in grads namelist! ',dd
437  call prc_mpistop
438  endif
439  case("levels")
440  if ( lnum < 0 )then
441  write(*,*) 'xxx "lnum" is required in grads namelist for levels data! '
442  call prc_mpistop
443  endif
444  do k=1, lnum
445  lvars(k)=grads_lvars(k,ielem,1)
446  enddo
447  if(abs(lvars(1)-large_number_one)<eps)then
448  write(*,*) 'xxx "lvars" must be specified in grads namelist for levels data! '
449  call prc_mpistop
450  endif
451  case("map")
452  startrec = grads_startrec(ielem,1)
453  totalrec = grads_totalrec(ielem,1)
454  fendian = grads_fendian(ielem,1)
455  yrev = grads_yrev(ielem,1)
456  if( (startrec<0).or.(totalrec<0) )then
457  write(*,*) 'xxx "startrec" is required in grads namelist! ',startrec
458  write(*,*) 'xxx "totalrec" is required in grads namelist! ',totalrec
459  call prc_mpistop
460  endif
461  ! get file_id
462  if(io_fid_grads_data < 0)then
463  io_fid_grads_data = io_get_available_fid()
464  endif
465  gfile=trim(fname)//trim(basename_num)//'.grd'
466  if( len_trim(fname)==0 )then
467  write(*,*) 'xxx "fname" is required in grads namelist for map data! ',trim(fname)
468  call prc_mpistop
469  endif
470  end select
471 
472  ! read data
473  select case (trim(item))
474  case("lon")
475  if ( trim(dtype) == "linear" ) then
476  do j = 1, dims(3)
477  do i = 1, dims(2)
478  lon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * D2R
479  enddo
480  enddo
481  else if ( trim(dtype) == "map" ) then
482  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,1,item,startrec,totalrec,yrev,gdata2d)
483  lon_org(:,:) = real(gdata2D(:,:), kind=RP) * D2R
484  endif
485  case("lat")
486  if ( trim(dtype) == "linear" ) then
487  do j = 1, dims(3)
488  do i = 1, dims(2)
489  lat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * D2R
490  enddo
491  enddo
492  else if ( trim(dtype) == "map" ) then
493  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,1,item,startrec,totalrec,yrev,gdata2d)
494  lat_org(:,:) = real(gdata2D(:,:), kind=RP) * D2R
495  endif
496  case("plev")
497  if(dims(1)/=knum)then
498  write(*,*) 'xxx "knum" must be equal to outer_nz for plev. knum:',knum,'> outer_nz:',dims(1)
499  call prc_mpistop
500  endif
501  if ( trim(dtype) == "levels" ) then
502  if(dims(1)/=lnum)then
503  write(*,*) 'xxx lnum must be same as the outer_nz for plev! ',dims(1),lnum
504  call prc_mpistop
505  endif
506  do j = 1, dims(3)
507  do i = 1, dims(2)
508  do k = 1, dims(1)
509  pres_org(k+2,i,j) = real(lvars(k), kind=rp)
510  enddo
511  enddo
512  enddo
513  else if ( trim(dtype) == "map" ) then
514  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),dims(1),nt,item,startrec,totalrec,yrev,gdata3d)
515  do j = 1, dims(3)
516  do i = 1, dims(2)
517  do k = 1, dims(1)
518  pres_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
519  enddo
520  enddo
521  enddo
522  endif
523  case('U')
524  if ( trim(dtype) == "map" ) then
525  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
526  do j = 1, dims(3)
527  do i = 1, dims(2)
528  velx_org(1:2,i,j) = 0.0_rp
529  do k = 1, knum
530  velx_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
531  enddo
532  if(dims(1)>knum)then
533  do k = knum+1, dims(1)
534  velx_org(k+2,i,j) = velx_org(knum+2,i,j)
535  enddo
536  endif
537  enddo
538  enddo
539  endif
540  case('V')
541  if ( trim(dtype) == "map" ) then
542  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
543  do j = 1, dims(3)
544  do i = 1, dims(2)
545  vely_org(1:2,i,j) = 0.0_rp
546  do k = 1, knum
547  vely_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
548  enddo
549  if(dims(1)>knum)then
550  do k = knum+1, dims(1)
551  vely_org(k+2,i,j) = vely_org(knum+2,i,j)
552  enddo
553  endif
554  enddo
555  enddo
556  endif
557  case('T')
558  if ( trim(dtype) == "map" ) then
559  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
560  do j = 1, dims(3)
561  do i = 1, dims(2)
562  do k = 1, knum
563  temp_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
564  enddo
565  if(dims(1)>knum)then
566  do k = knum+1, dims(1)
567  temp_org(k+2,i,j) = temp_org(knum+2,i,j)
568  enddo
569  endif
570  enddo
571  enddo
572  endif
573  case('HGT')
574  if(dims(1)/=knum)then
575  write(*,*) 'xxx The number of levels for HGT must be same as plevs! knum:',knum,'> outer_nz:',dims(1)
576  call prc_mpistop
577  endif
578  if ( trim(dtype) == "map" ) then
579  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),dims(1),nt,item,startrec,totalrec,yrev,gdata3d)
580  do j = 1, dims(3)
581  do i = 1, dims(2)
582  do k = 1, dims(1)
583  cz_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
584  enddo
585  cz_org(1,i,j) = 0.0_rp
586  enddo
587  enddo
588  endif
589  case('QV')
590  if ( trim(dtype) == "map" ) then
591  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
592  do j = 1, dims(3)
593  do i = 1, dims(2)
594  do k = 1, knum
595  qtrc_org(k+2,i,j,qa_outer) = real(gdata3D(i,j,k), kind=rp)
596  enddo
597  qtrc_org(1:2,i,j,qa_outer) = qtrc_org(3,i,j,qa_outer)
598  enddo
599  enddo
600  if( dims(1)>knum ) then
601  select case ( upper_qv_type )
602  case ("COPY")
603  do j = 1, dims(3)
604  do i = 1, dims(2)
605  do k = knum+1, dims(1)
606  qtrc_org(k+2,i,j,qa_outer) = qtrc_org(knum+2,i,j,qa_outer)
607  enddo
608  enddo
609  enddo
610  case ("ZERO")
611  ! do nothing
612  case default
613  write(*,*) 'xxx upper_qv_type in PARAM_MKINIT_REAL_GrADS is invalid! ', upper_qv_type
614  call prc_mpistop
615  end select
616  endif
617  endif
618  case('RH')
619  if (data_available(ia_qv,1)) cycle ! use QV
620  if ( trim(dtype) == "map" ) then
621  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
622  do j = 1, dims(3)
623  do i = 1, dims(2)
624  do k = 1, knum
625  rhprs_org(k+2,i,j) = real(gdata3D(i,j,k), kind=RP) / 100.0_RP ! relative humidity
626  call psat( p_sat, temp_org(k+2,i,j) ) ! satulation pressure
627  qm = epsvap * rhprs_org(k+2,i,j) * p_sat &
628  / ( pres_org(k+2,i,j) - rhprs_org(k+2,i,j) * p_sat ) ! mixing ratio
629  qtrc_org(k+2,i,j,qa_outer) = qm / ( 1.0_rp + qm ) ! specific humidity
630  enddo
631  qtrc_org(1:2,i,j,qa_outer) = qtrc_org(3,i,j,qa_outer)
632  enddo
633  enddo
634  if( dims(3)>knum ) then
635  select case ( upper_qv_type )
636  case ("COPY")
637  do j = 1, dims(3)
638  do i = 1, dims(2)
639  do k = knum+1, dims(1)
640  rhprs_org(k+2,i,j) = rhprs_org(knum+2,i,j) ! relative humidity
641  call psat( p_sat, temp_org(k+2,i,j) ) ! satulated specific humidity
642  qm = epsvap * rhprs_org(k+2,i,j) * p_sat &
643  / ( pres_org(k+2,i,j) - rhprs_org(k+2,i,j) * p_sat ) ! mixing ratio
644  qtrc_org(k+2,i,j,qa_outer) = qm / ( 1.0_rp + qm ) ! specific humidity
645  qtrc_org(k+2,i,j,qa_outer) = min(qtrc_org(k+2,i,j,qa_outer),qtrc_org(k+1,i,j,qa_outer))
646  enddo
647  enddo
648  enddo
649  case ("ZERO")
650  ! do nothing
651  case default
652  write(*,*) 'xxx upper_qv_type in PARAM_MKINIT_REAL_GrADS is invalid! ', upper_qv_type
653  call prc_mpistop
654  end select
655  endif
656  endif
657  case('MSLP')
658  if ( trim(dtype) == "map" ) then
659  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
660  do j = 1, dims(3)
661  do i = 1, dims(2)
662  pres_org(1,i,j) = real(gdata2D(i,j), kind=rp)
663  enddo
664  enddo
665  endif
666  case('PSFC')
667  if ( trim(dtype) == "map" ) then
668  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
669  do j = 1, dims(3)
670  do i = 1, dims(2)
671  pres_org(2,i,j) = real(gdata2D(i,j), kind=rp)
672  enddo
673  enddo
674  endif
675  case('U10')
676  if ( trim(dtype) == "map" ) then
677  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
678  do j = 1, dims(3)
679  do i = 1, dims(2)
680  velx_org(2,i,j) = real(gdata2D(i,j), kind=rp)
681  enddo
682  enddo
683  endif
684  case('V10')
685  if ( trim(dtype) == "map" ) then
686  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
687  do j = 1, dims(3)
688  do i = 1, dims(2)
689  vely_org(2,i,j) = real(gdata2D(i,j), kind=rp)
690  enddo
691  enddo
692  endif
693  case('T2')
694  if ( trim(dtype) == "map" ) then
695  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
696  do j = 1, dims(3)
697  do i = 1, dims(2)
698  temp_org(2,i,j) = real(gdata2D(i,j), kind=rp)
699  enddo
700  enddo
701  endif
702  case('Q2')
703  if ( trim(dtype) == "map" ) then
704  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
705  do j = 1, dims(3)
706  do i = 1, dims(2)
707  qtrc_org(2,i,j,qa_outer) = real(gdata2D(i,j), kind=rp)
708  enddo
709  enddo
710  endif
711  case('RH2')
712  if (data_available(ia_q2,1)) cycle ! use QV
713  if ( trim(dtype) == "map" ) then
714  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
715  do j = 1, dims(3)
716  do i = 1, dims(2)
717  rhsfc = real(gdata2D(i,j), kind=RP) / 100.0_RP ! relative humidity
718  call psat( p_sat, temp_org(2,i,j) ) ! satulation pressure
719  qm = epsvap * rhsfc * p_sat &
720  / ( pres_org(2,i,j) - rhsfc * p_sat ) ! mixing ratio
721  qtrc_org(2,i,j,qa_outer) = qm / ( 1.0_rp + qm ) ! specific humidity
722  enddo
723  enddo
724  endif
725  case('TOPO')
726  if ( trim(dtype) == "map" ) then
727  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
728  do j = 1, dims(3)
729  do i = 1, dims(2)
730  cz_org(2,i,j) = real(gdata2D(i,j), kind=rp)
731  enddo
732  enddo
733  endif
734  end select
735  enddo loop_inputatomgrads
736 
737 
738  rovcp = rdry / cpdry
739  cpovr = cpdry / rdry
740 
741  if ( data_available(ia_t2,1) .and. data_available(ia_ps,1) ) then
742  do j = 1, dims(3)
743  do i = 1, dims(2)
744  pott(i,j) = temp_org(2,i,j) * (p00/pres_org(2,i,j))**rovcp
745  end do
746  end do
747  else
748  do j = 1, dims(3)
749  do i = 1, dims(2)
750  pott(i,j) = temp_org(3,i,j) * (p00/pres_org(3,i,j))**rovcp
751  end do
752  end do
753  end if
754 
755  if ( .not. data_available(ia_t2,1) ) then
756  if ( data_available(ia_ps,1) ) then
757  do j = 1, dims(3)
758  do i = 1, dims(2)
759  temp_org(2,i,j) = pott(i,j) * (pres_org(2,i,j)/p00)**rovcp
760  end do
761  end do
762  else
763  if ( data_available(ia_topo,1) ) then
764  do j = 1, dims(3)
765  do i = 1, dims(2)
766  temp_org(2,i,j) = temp_org(3,i,j) &
767  + laps * (cz_org(3,i,j)-cz_org(2,i,j))
768  end do
769  end do
770  else
771  do j = 1, dims(3)
772  do i = 1, dims(2)
773  temp_org(2,i,j) = temp_org(3,i,j)
774  end do
775  end do
776  end if
777  end if
778  end if
779 
780  if ( .not. data_available(ia_ps,1) ) then
781  do j = 1, dims(3)
782  do i = 1, dims(2)
783  pres_org(2,i,j) = p00 * (temp_org(2,i,j)/pott(i,j))**cpovr
784  end do
785  end do
786  end if
787 
788  if ( data_available(ia_slp,1) ) then
789  do j = 1, dims(3)
790  do i = 1, dims(2)
791  temp_org(1,i,j) = pott(i,j) * (pres_org(1,i,j)/p00)**rovcp
792  end do
793  end do
794  else
795  if ( data_available(ia_t2,1) .and. data_available(ia_topo,1) ) then
796  do j = 1, dims(3)
797  do i = 1, dims(2)
798  temp_org(1,i,j) = temp_org(2,i,j) + laps * cz_org(2,i,j)
799  end do
800  end do
801  else
802  do j = 1, dims(3)
803  do i = 1, dims(2)
804  temp_org(1,i,j) = temp_org(3,i,j) + laps * cz_org(3,i,j)
805  end do
806  end do
807  end if
808  do j = 1, dims(3)
809  do i = 1, dims(2)
810  pres_org(1,i,j) = p00 * (temp_org(1,i,j)/pott(i,j))**cpovr
811  end do
812  end do
813  end if
814 
815  if ( .not. data_available(ia_topo,1) ) then
816  ! guess surface height (elevation)
817  do j = 1, dims(3)
818  do i = 1, dims(2)
819  cz_org(2,i,j) = max( 0.0_rp, &
820  cz_org(3,i,j) &
821  * ( log(pres_org(2,i,j)/pres_org(1,i,j)) ) &
822  / ( log(pres_org(3,i,j)/pres_org(1,i,j)) ) )
823  end do
824  end do
825  end if
826 
827  velz_org = 0.0_rp
828 
829 
830  ! check verticaly extrapolated data in outer model
831  do j = 1, dims(3)
832  do i = 1, dims(2)
833  do k = 3, dims(1)+2
834 
835  if( pres_org(k,i,j)>pres_org(2,i,j) )then ! if Pressure is larger than Surface pressure
836  !velz_org(k,i,j)=velz_org(2,i,j)
837  velx_org(k,i,j)=velx_org(2,i,j)
838  vely_org(k,i,j)=vely_org(2,i,j)
839  temp_org(k,i,j)=temp_org(2,i,j)
840  qtrc_org(k,i,j,:)=qtrc_org(2,i,j,:)
841  cz_org(k,i,j)=cz_org(2,i,j)
842  end if
843 
844  enddo
845  enddo
846  enddo
847 
848 
849  !do it = 1, nt
850  ! k=1 ; j=int(dims(3)/2) ; i=int(dims(2)/2) ; iq = 1
851  ! write(*,*) "read 3D grads data",i,j,k
852  ! write(*,*) "lon_org ",lon_org (i,j)/D2R
853  ! write(*,*) "lat_org ",lat_org (i,j)/D2R
854  ! write(*,*) "pres_org ",pres_org (k,i,j)
855  ! write(*,*) "usfc_org ",usfc_org (i,j)
856  ! write(*,*) "vsfc_org ",vsfc_org (i,j)
857  ! write(*,*) "tsfc_org ",tsfc_org (i,j)
858  ! write(*,*) "qsfc_org ",qsfc_org (i,j,iq)
859  ! write(*,*) "rhsfc_org ",rhsfc_org (i,j)
860  ! write(*,*) "velx_org ",velx_org (k,i,j)
861  ! write(*,*) "vely_org ",vely_org (k,i,j)
862  ! write(*,*) "temp_org ",temp_org (k,i,j)
863  ! write(*,*) "hgt_org ",hgt_org (k,i,j)
864  ! write(*,*) "qtrc_org ",qtrc_org (k,i,j,iq)
865  ! write(*,*) "rhprs_org ",rhprs_org (k,i,j)
866  !enddo
867 
868  return
869  end subroutine parentatominputgrads
870 
871  !-----------------------------------------------------------------------------
873  subroutine parentlandsetupgrads( &
874  ldims, & ! (out)
875  use_waterratio, & ! (out)
876  use_file_landwater, & ! (in)
877  basename )
878  implicit none
879 
880  integer, intent(out) :: ldims(3)
881  logical, intent(out) :: use_waterratio
882  logical, intent(in) :: use_file_landwater ! use landwater data from files
883  character(len=*), intent(in) :: basename
884 
885  integer :: ielem
886  integer :: k, n
887 
888  integer :: ierr
889  !---------------------------------------------------------------------------
890 
891  if( io_l ) write(io_fid_log,*) '+++ Real Case/Land Input File Type: GrADS format'
892 
893  !--- initialization
894  use_waterratio = .false.
895 
896  if ( len_trim(basename) == 0 ) then
897  if( io_l ) write(io_fid_log,*) &
898  'xxx "BASEMAAME" is not specified in "PARAM_MKINIT_REAL_ATOMS"!',trim(basename)
899  call prc_mpistop
900  endif
901 
902  !--- read namelist
903  io_fid_grads_nml = io_get_available_fid()
904  open( io_fid_grads_nml, &
905  file = trim(basename), &
906  form = 'formatted', &
907  status = 'old', &
908  action = 'read', &
909  iostat = ierr )
910  if ( ierr /= 0 ) then
911  if( io_l ) write(io_fid_log,*) 'xxx Input file is not found! ', trim(basename)
912  call prc_mpistop
913  endif
914 
915  read(io_fid_grads_nml,nml=nml_grads_grid,iostat=ierr)
916  if( ierr /= 0 ) then !--- missing or fatal error
917  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in nml_grads_grid in ', trim(basename),'. Check!'
918  call prc_mpistop
919  endif
920  if( io_lnml ) write(io_fid_log,nml=nml_grads_grid)
921 
922  ! land
923  ldims(1) = outer_nl ! soil_layers_stag
924  if(outer_nx_sfc > 0)then
925  ldims(2) = outer_nx_sfc
926  else
927  ldims(2) = outer_nx
928  outer_nx_sfc = outer_nx
929  endif
930  if(outer_ny_sfc > 0)then
931  ldims(3) = outer_ny_sfc
932  else
933  ldims(3) = outer_ny
934  outer_ny_sfc = outer_ny
935  endif
936 
937  allocate( gland2d( ldims(2), ldims(3) ) )
938  allocate( gland3d( ldims(2), ldims(3), ldims(1) ) )
939 
940  call read_namelist( &
941  grads_item(:,2), & ! (out)
942  grads_fname(:,2), & ! (out)
943  grads_dtype(:,2), & ! (out)
944  grads_swpoint(:,2), & ! (out)
945  grads_dd(:,2), & ! (out)
946  grads_lnum(:,2), & ! (out)
947  grads_lvars(:,:,2), & ! (out)
948  grads_startrec(:,2), & ! (out)
949  grads_totalrec(:,2), & ! (out)
950  grads_knum(:,2), & ! (out)
951  grads_yrev(:,2), & ! (out)
952  grads_fendian(:,2), & ! (out)
953  grads_missval(:,2), & ! (out)
954  data_available(:,2), & ! (out)
955  item_list_land, & ! (in)
956  num_item_list_land, & ! (in)
957  basename, & ! (in)
958  io_fid_grads_nml ) ! (in)
959 
960  close( io_fid_grads_nml )
961 
962  do ielem = 1, num_item_list_land
963  item = item_list_land(ielem)
964  !--- check data
965  select case(trim(item))
966  case('TOPO','lsmask')
967  if ( .not. data_available(ielem,2) ) then
968  if (io_l) write(io_fid_log,*) 'warning: ',trim(item),' is not found & not used.'
969  cycle
970  endif
971  case('lon', 'lat', 'lon_sfc', 'lat_sfc')
972  cycle
973  case('SMOISVC', 'SMOISDS')
974  if ( use_file_landwater ) then
975  if (.not. data_available(il_smoisvc,2) .and. .not. data_available(il_smoisds,2)) then
976  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : ',trim(item_list_land(ielem))
977  call prc_mpistop
978  end if
979  use_waterratio = data_available(il_smoisds,2)
980  else
981  cycle
982  end if
983  case default ! llev, SKINT, STEMP
984  if ( .not. data_available(ielem,2) ) then
985  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : ',trim(item_list_land(ielem))
986  call prc_mpistop
987  endif
988  end select
989 
990  end do
991 
992  return
993  end subroutine parentlandsetupgrads
994 
995  subroutine parentlandinputgrads( &
996  tg_org, & ! (out)
997  strg_org, & ! (out)
998  smds_org, & ! (out)
999  lst_org, & ! (out)
1000  llon_org, & ! (out)
1001  llat_org, & ! (out)
1002  lz_org, & ! (out)
1003  topo_org, & ! (out)
1004  lmask_org, & ! (out)
1005  basename_num, & ! (in)
1006  ldims, & ! (in)
1007  use_file_landwater, & ! (in)
1008  nt ) ! (in)
1009  use scale_const, only: &
1010  undef => const_undef, &
1011  d2r => const_d2r, &
1012  tem00 => const_tem00, &
1013  eps => const_eps
1014  implicit none
1015 
1016  real(RP), intent(out) :: tg_org (:,:,:)
1017  real(RP), intent(out) :: strg_org (:,:,:)
1018  real(RP), intent(out) :: smds_org (:,:,:)
1019  real(RP), intent(out) :: lst_org (:,:)
1020  real(RP), intent(out) :: llon_org (:,:)
1021  real(RP), intent(out) :: llat_org (:,:)
1022  real(RP), intent(out) :: lz_org (:)
1023  real(RP), intent(out) :: topo_org(:,:)
1024  real(RP), intent(out) :: lmask_org(:,:)
1025 
1026  character(len=*), intent(in) :: basename_num
1027  integer, intent(in) :: ldims(3)
1028  logical, intent(in) :: use_file_landwater ! use land water data from files
1029  integer, intent(in) :: nt
1030  ! ----------------------------------------------------------------
1031 
1033  character(len=H_LONG) :: gfile
1034 
1035  real(RP) :: qvsat, qm
1036 
1037  integer :: i, j, k, ielem, n
1038  integer :: ierr
1039 
1040  !---------------------------------------------------------------------------
1041 
1042  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[LandInputGrADS]'
1043 
1044 
1045  loop_inputlandgrads : do ielem = 1, num_item_list_land
1046 
1047  item = item_list_land(ielem)
1048 
1049  dtype = grads_dtype(ielem,2)
1050  fname = grads_fname(ielem,2)
1051  lnum = grads_lnum(ielem,2)
1052  missval = grads_missval(ielem,2)
1053 
1054  if ( grads_knum(ielem,2) > 0 )then
1055  knum = grads_knum(ielem,2)
1056  else
1057  knum = ldims(1)
1058  endif
1059 
1060  select case (trim(dtype))
1061  case("linear")
1062  swpoint = grads_swpoint(ielem,2)
1063  dd = grads_dd(ielem,2)
1064  if( (abs(swpoint-large_number_one)<eps).or.(abs(dd-large_number_one)<eps) )then
1065  write(*,*) 'xxx "swpoint" is required in grads namelist! ',swpoint
1066  write(*,*) 'xxx "dd" is required in grads namelist! ',dd
1067  call prc_mpistop
1068  endif
1069  case("levels")
1070  if ( lnum < 0 )then
1071  write(*,*) 'xxx "lnum" in grads namelist is required for levels data! '
1072  call prc_mpistop
1073  endif
1074  do k=1, lnum
1075  lvars(k)=grads_lvars(k,ielem,2)
1076  enddo
1077  if(abs(lvars(1)-large_number_one)<eps)then
1078  write(*,*) 'xxx "lvars" must be specified in grads namelist for levels data!',(lvars(k),k=1,lnum)
1079  call prc_mpistop
1080  endif
1081  case("map")
1082  startrec = grads_startrec(ielem,2)
1083  totalrec = grads_totalrec(ielem,2)
1084  fendian = grads_fendian(ielem,2)
1085  yrev = grads_yrev(ielem,2)
1086  if( (startrec<0).or.(totalrec<0) )then
1087  write(*,*) 'xxx "startrec" is required in grads namelist! ',startrec
1088  write(*,*) 'xxx "totalrec" is required in grads namelist! ',totalrec
1089  call prc_mpistop
1090  endif
1091  ! get file_io
1092  if(io_fid_grads_data < 0)then
1093  io_fid_grads_data = io_get_available_fid()
1094  endif
1095  gfile=trim(fname)//trim(basename_num)//'.grd'
1096  if( len_trim(fname)==0 )then
1097  write(*,*) 'xxx "fname" is required in grads namelist for map data! ',trim(fname)
1098  call prc_mpistop
1099  endif
1100  end select
1101 
1102  ! read data
1103  select case (trim(item))
1104  case("lsmask")
1105  if ( data_available(il_lsmask,2) ) then
1106  if ( trim(dtype) == "map" ) then
1107  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1108  lmask_org(:,:) = real(gland2D(:,:), kind=rp)
1109  endif
1110  else
1111  lmask_org = undef
1112  end if
1113  case("lon")
1114  if ( .not. data_available(il_lon_sfc,2) ) then
1115  if ( ldims(2).ne.outer_nx .or. ldims(3).ne.outer_ny ) then
1116  write(*,*) 'xxx namelist of "lon_sfc" is not found in grads namelist!'
1117  write(*,*) 'xxx dimension is different: outer_nx and outer_nx_sfc! ', outer_nx, ldims(2)
1118  write(*,*) ' : outer_ny and outer_ny_sfc! ', outer_ny, ldims(3)
1119  call prc_mpistop
1120  end if
1121  if ( trim(dtype) == "linear" ) then
1122  do j = 1, ldims(3)
1123  do i = 1, ldims(2)
1124  llon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * D2R
1125  enddo
1126  enddo
1127  else if ( trim(dtype) == "map" ) then
1128  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1129  llon_org(:,:) = real(gland2D(:,:), kind=RP) * D2R
1130  endif
1131  end if
1132  case("lon_sfc")
1133  if ( .not. data_available(il_lon_sfc,2) ) cycle
1134  if ( trim(dtype) == "linear" ) then
1135  do j = 1, ldims(3)
1136  do i = 1, ldims(2)
1137  llon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * D2R
1138  enddo
1139  enddo
1140  else if ( trim(dtype) == "map" ) then
1141  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1142  llon_org(:,:) = real(gland2D(:,:), kind=RP) * D2R
1143  endif
1144  case("lat")
1145  if ( .not. data_available(il_lat_sfc,2) ) then
1146  if ( ldims(2).ne.outer_nx .or. ldims(3).ne.outer_ny ) then
1147  write(*,*) 'xxx namelist of "lat_sfc" is not found in grads namelist!'
1148  write(*,*) 'xxx dimension is different: outer_nx and outer_nx_sfc! ', outer_nx, ldims(2)
1149  write(*,*) ' : outer_ny and outer_ny_sfc! ', outer_nx, ldims(3)
1150  call prc_mpistop
1151  end if
1152  if ( trim(dtype) == "linear" ) then
1153  do j = 1, ldims(3)
1154  do i = 1, ldims(2)
1155  llat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * D2R
1156  enddo
1157  enddo
1158  else if ( trim(dtype) == "map" ) then
1159  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1160  llat_org(:,:) = real(gland2D(:,:), kind=RP) * D2R
1161  endif
1162  end if
1163  case("lat_sfc")
1164  if ( .not. data_available(il_lat_sfc,2) ) cycle
1165  if ( trim(dtype) == "linear" ) then
1166  do j = 1, ldims(3)
1167  do i = 1, ldims(2)
1168  llat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * D2R
1169  enddo
1170  enddo
1171  else if ( trim(dtype) == "map" ) then
1172  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1173  llat_org(:,:) = real(gland2D(:,:), kind=RP) * D2R
1174  endif
1175  case("llev")
1176  if(ldims(1)/=knum)then
1177  write(*,*) 'xxx "knum" must be equal to outer_nl for llev. knum:',knum,'> outer_nl:',ldims(1)
1178  call prc_mpistop
1179  endif
1180  if ( trim(dtype) == "levels" ) then
1181  if(ldims(1)/=lnum)then
1182  write(*,*) 'xxx lnum must be same as the outer_nl for llev! ',ldims(1),lnum
1183  call prc_mpistop
1184  endif
1185  do k = 1, ldims(1)
1186  lz_org(k) = real(lvars(k), kind=rp)
1187  enddo
1188 ! else if ( trim(dtype) == "map" ) then
1189 ! call read_grads_file_3d(io_fid_grads_data,gfile,ldims(2),ldims(3),ldims(1),nt,item,startrec,totalrec,yrev,gland)
1190 ! do j = 1, ldims(3)
1191 ! do i = 1, ldims(2)
1192 ! do k = 1, ldims(1)
1193 ! lz_org(k,i,j) = real(gland(i,j,k), kind=RP)
1194 ! enddo
1195 ! enddo
1196 ! enddo
1197  endif
1198  case('STEMP')
1199  if(ldims(1)/=knum)then
1200  write(*,*) 'xxx The number of levels for STEMP must be same as llevs! ',ldims(1),knum
1201  call prc_mpistop
1202  endif
1203  if ( trim(dtype) == "map" ) then
1204  call read_grads_file_3d(io_fid_grads_data,gfile,ldims(2),ldims(3),ldims(1),nt,item,startrec,totalrec,yrev,gland3d)
1205  do j = 1, ldims(3)
1206  do i = 1, ldims(2)
1207  do k = 1, ldims(1)
1208  if ( abs(gland3d(i,j,k)-missval) < eps ) then
1209  tg_org(k,i,j) = undef
1210  else
1211  tg_org(k,i,j) = real(gland3D(i,j,k), kind=rp)
1212  end if
1213  enddo
1214  enddo
1215  enddo
1216  endif
1217  case('SMOISVC')
1218  if ( use_file_landwater ) then
1219  if(ldims(1)/=knum)then
1220  write(*,*) 'xxx The number of levels for SMOISVC must be same as llevs! ',ldims(1),knum
1221  call prc_mpistop
1222  endif
1223  if ( trim(dtype) == "map" ) then
1224  call read_grads_file_3d(io_fid_grads_data,gfile,ldims(2),ldims(3),ldims(1),nt,item,startrec,totalrec,yrev,gland3d)
1225  do j = 1, ldims(3)
1226  do i = 1, ldims(2)
1227  do k = 1, ldims(1)
1228  if ( abs(gland3d(i,j,k)-missval) < eps ) then
1229  strg_org(k,i,j) = undef
1230  else
1231  strg_org(k,i,j) = real(gland3D(i,j,k), kind=rp)
1232  end if
1233  enddo
1234  enddo
1235  enddo
1236  endif
1237  endif
1238  case('SMOISDS')
1239  if ( use_file_landwater ) then
1240  if(ldims(1)/=knum)then
1241  write(*,*) 'xxx The number of levels for SMOISDS must be same as llevs! ',ldims(1),knum
1242  call prc_mpistop
1243  endif
1244  if ( trim(dtype) == "map" ) then
1245  call read_grads_file_3d(io_fid_grads_data,gfile,ldims(2),ldims(3),ldims(1),nt,item,startrec,totalrec,yrev,gland3d)
1246  do j = 1, ldims(3)
1247  do i = 1, ldims(2)
1248  do k = 1, ldims(1)
1249  if ( abs(gland3d(i,j,k)-missval) < eps ) then
1250  smds_org(k,i,j) = undef
1251  else
1252  smds_org(k,i,j) = real(gland3D(i,j,k), kind=rp)
1253  end if
1254  enddo
1255  enddo
1256  enddo
1257  endif
1258  endif
1259  case('SKINT')
1260  if ( trim(dtype) == "map" ) then
1261  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,nt,item,startrec,totalrec,yrev,gland2d)
1262  do j = 1, ldims(3)
1263  do i = 1, ldims(2)
1264  if ( abs(gland2d(i,j)-missval) < eps ) then
1265  lst_org(i,j) = undef
1266  else
1267  lst_org(i,j) = real(gland2D(i,j), kind=rp)
1268  end if
1269  enddo
1270  enddo
1271  endif
1272  case('TOPO')
1273  if ( data_available(il_topo,2) ) then
1274  if ( trim(dtype) == "map" ) then
1275  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,nt,item,startrec,totalrec,yrev,gland2d)
1276  do j = 1, ldims(3)
1277  do i = 1, ldims(2)
1278  if ( abs(gland2d(i,j)-missval) < eps ) then
1279  topo_org(i,j) = undef
1280  else
1281  topo_org(i,j) = real(gland2D(i,j), kind=rp)
1282  end if
1283  enddo
1284  enddo
1285  endif
1286  else
1287  topo_org = undef
1288  endif
1289  end select
1290  enddo loop_inputlandgrads
1291 
1292  !do it = 1, nt
1293  ! i=int(ldims(2)/2) ; j=int(ldims(3)/2)
1294  ! write(*,*) "read 2D grads data",ldims(2),ldims(3),i,j,it
1295  ! write(*,*) "lon_org ",lon_org (i,j)
1296  ! write(*,*) "lat_org ",lat_org (i,j)
1297  ! write(*,*) "lst_org ",lst_org(i,j)
1298  ! do k=1,dims(7)
1299  ! write(*,*) "tg_org ",tg_org (k,i,j)," k= ",k
1300  ! write(*,*) "strg_org ",strg_org (k,i,j)," k= ",k
1301  ! enddo
1302  !enddo
1303 
1304  return
1305  end subroutine parentlandinputgrads
1306 
1307  !-----------------------------------------------------------------------------
1309  subroutine parentoceansetupgrads( &
1310  odims, & ! (out)
1311  timelen, & ! (out)
1312  basename ) ! (in)
1313  implicit none
1314 
1315  integer, intent(out) :: odims(2)
1316  integer, intent(out) :: timelen
1317  character(len=*), intent(in) :: basename
1318 
1319  character(len=H_LONG) :: grads_ctl
1320  integer :: ielem
1321  integer :: n
1322 
1323  integer :: ierr
1324  !---------------------------------------------------------------------------
1325 
1326  if( io_l ) write(io_fid_log,*) '+++ Real Case/Ocean Input File Type: GrADS format'
1327 
1328  !--- read namelist
1329 
1330  if ( len_trim(basename) == 0 ) then
1331  grads_ctl = "namelist.grads_boundary"
1332  else
1333  grads_ctl = basename
1334  endif
1335 
1336  !--- read namelist
1337  io_fid_grads_nml = io_get_available_fid()
1338  open( io_fid_grads_nml, &
1339  file = trim(grads_ctl), &
1340  form = 'formatted', &
1341  status = 'old', &
1342  action = 'read', &
1343  iostat = ierr )
1344  if ( ierr /= 0 ) then
1345  if( io_l ) write(io_fid_log,*) 'xxx Input file is not found! ', trim(grads_ctl)
1346  call prc_mpistop
1347  endif
1348 
1349  read(io_fid_grads_nml,nml=nml_grads_grid,iostat=ierr)
1350  if( ierr /= 0 ) then !--- missing or fatal error
1351  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in nml_grads_grid in ', trim(grads_ctl),'. Check!'
1352  call prc_mpistop
1353  endif
1354  if( io_lnml ) write(io_fid_log,nml=nml_grads_grid)
1355 
1356  timelen = 0 ! will be replaced later
1357 
1358  ! sst
1359  if(outer_nx_sst > 0)then
1360  odims(1) = outer_nx_sst
1361  else if (outer_nx_sfc > 0) then
1362  odims(1) = outer_nx_sfc
1363  outer_nx_sst = outer_nx_sfc
1364  else
1365  odims(1) = outer_nx
1366  outer_nx_sst = outer_nx
1367  endif
1368  if(outer_ny_sst > 0)then
1369  odims(2) = outer_ny_sst
1370  else if(outer_ny_sfc > 0)then
1371  odims(2) = outer_ny_sfc
1372  outer_ny_sst = outer_ny_sfc
1373  else
1374  odims(2) = outer_ny
1375  outer_ny_sst = outer_ny
1376  endif
1377 
1378  allocate( gsst2d( odims(1), odims(2) ) )
1379 
1380 
1381  call read_namelist( &
1382  grads_item(:,3), & ! (out)
1383  grads_fname(:,3), & ! (out)
1384  grads_dtype(:,3), & ! (out)
1385  grads_swpoint(:,3), & ! (out)
1386  grads_dd(:,3), & ! (out)
1387  grads_lnum(:,3), & ! (out)
1388  grads_lvars(:,:,3), & ! (out)
1389  grads_startrec(:,3), & ! (out)
1390  grads_totalrec(:,3), & ! (out)
1391  grads_knum(:,3), & ! (out)
1392  grads_yrev(:,3), & ! (out)
1393  grads_fendian(:,3), & ! (out)
1394  grads_missval(:,3), & ! (out)
1395  data_available(:,3), & ! (out)
1396  item_list_ocean, & ! (in)
1397  num_item_list_ocean, & ! (in)
1398  grads_ctl, & ! (in)
1399  io_fid_grads_nml ) ! (in)
1400 
1401  close( io_fid_grads_nml )
1402 
1403  do ielem = 1, num_item_list_ocean
1404  item = item_list_ocean(ielem)
1405  !--- check data
1406  select case(trim(item))
1407  case('lsmask')
1408  if ( .not. data_available(ielem,3) ) then
1409  if (io_l) write(io_fid_log,*) 'warning: ',trim(item),' is not found & not used.'
1410  cycle
1411  endif
1412  case('lon', 'lat', 'lon_sfc', 'lat_sfc', 'lon_sst', 'lat_sst')
1413  cycle
1414  case('SST')
1415  if (.not. data_available(io_sst,3) .and. .not. data_available(io_skint,3) ) then
1416  if (io_l) write(io_fid_log,*) 'xxx SST and SKINT are found in grads namelist!'
1417  call prc_mpistop
1418  endif
1419  if (.not. data_available(io_sst,3)) then
1420  if (io_l) write(io_fid_log,*) 'warning: SST is found in grads namelist. SKINT is used in place of SST.'
1421  cycle
1422  endif
1423  case('SKINT')
1424  cycle
1425  case default !
1426  if ( .not. data_available(ielem,3) ) then
1427  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : ',trim(item_list_ocean(ielem))
1428  call prc_mpistop
1429  endif
1430  end select
1431 
1432  end do
1433 
1434  return
1435  end subroutine parentoceansetupgrads
1436 
1437  !-----------------------------------------------------------------------------
1438  subroutine parentoceanopengrads
1439  implicit none
1440 
1441  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanOpenGrADS]'
1442 
1443  return
1444  end subroutine parentoceanopengrads
1445 
1446  !-----------------------------------------------------------------------------
1447  subroutine parentoceaninputgrads( &
1448  tw_org, & ! (out)
1449  sst_org, & ! (out)
1450  omask_org, & ! (out)
1451  olon_org, & ! (out)
1452  olat_org, & ! (out)
1453  basename_num, & ! (in)
1454  odims, & ! (in)
1455  nt ) ! (in)
1456  use scale_const, only: &
1457  undef => const_undef, &
1458  d2r => const_d2r, &
1459  tem00 => const_tem00, &
1460  eps => const_eps
1461  implicit none
1462 
1463  real(RP), intent(out) :: tw_org (:,:)
1464  real(RP), intent(out) :: sst_org (:,:)
1465  real(RP), intent(out) :: omask_org(:,:)
1466  real(RP), intent(out) :: olon_org (:,:)
1467  real(RP), intent(out) :: olat_org (:,:)
1468 
1469  character(len=*), intent(in) :: basename_num
1470  integer, intent(in) :: odims(2)
1471  integer, intent(in) :: nt
1472  ! ----------------------------------------------------------------
1473 
1475  character(len=H_LONG) :: gfile
1476 
1477  real(RP) :: qvsat, qm
1478 
1479  integer :: i, j, ielem, n
1480  integer :: ierr
1481 
1482  !---------------------------------------------------------------------------
1483 
1484  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanInputGrADS]'
1485 
1486 
1487  loop_inputoceangrads : do ielem = 1, num_item_list_ocean
1488 
1489  item = item_list_ocean(ielem)
1490 
1491  dtype = grads_dtype(ielem,3)
1492  fname = grads_fname(ielem,3)
1493  lnum = grads_lnum(ielem,3)
1494  missval = grads_missval(ielem,3)
1495 
1496  select case (trim(dtype))
1497  case("linear")
1498  swpoint = grads_swpoint(ielem,3)
1499  dd = grads_dd(ielem,3)
1500  if( (abs(swpoint-large_number_one)<eps).or.(abs(dd-large_number_one)<eps) )then
1501  write(*,*) 'xxx "swpoint" is required in grads namelist! ',swpoint
1502  write(*,*) 'xxx "dd" is required in grads namelist! ',dd
1503  call prc_mpistop
1504  endif
1505  case("levels")
1506  write(*,*) 'xxx "lnum" in grads namelist is invalid for ocean data'
1507  call prc_mpistop
1508  case("map")
1509  startrec = grads_startrec(ielem,3)
1510  totalrec = grads_totalrec(ielem,3)
1511  fendian = grads_fendian(ielem,3)
1512  yrev = grads_yrev(ielem,3)
1513  if( (startrec<0).or.(totalrec<0) )then
1514  write(*,*) 'xxx "startrec" is required in grads namelist! ',startrec
1515  write(*,*) 'xxx "totalrec" is required in grads namelist! ',totalrec
1516  call prc_mpistop
1517  endif
1518  ! get file_io
1519  if(io_fid_grads_data < 0)then
1520  io_fid_grads_data = io_get_available_fid()
1521  endif
1522  gfile=trim(fname)//trim(basename_num)//'.grd'
1523  if( len_trim(fname)==0 )then
1524  write(*,*) 'xxx "fname" is required in grads namelist for map data! ',trim(fname)
1525  call prc_mpistop
1526  endif
1527  end select
1528 
1529  ! read data
1530  select case (trim(item))
1531  case("lsmask")
1532  if ( data_available(io_lsmask,3) ) then
1533  if ( trim(dtype) == "map" ) then
1534  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1535  omask_org(:,:) = real(gsst2D(:,:), kind=rp)
1536  endif
1537  else
1538  omask_org = undef
1539  end if
1540  case("lon")
1541  if ( .not. data_available(io_lon_sst,3) .and. .not. data_available(io_lon_sfc,3) ) then
1542  if ( odims(1).ne.outer_nx .or. odims(2).ne.outer_ny ) then
1543  write(*,*) 'xxx namelist of "lon_sst" is not found in grads namelist!'
1544  write(*,*) 'xxx dimension is different: outer_nx and outer_nx_sst! ', outer_nx, odims(1)
1545  write(*,*) ' : outer_ny and outer_ny_sst! ', outer_ny, odims(2)
1546  call prc_mpistop
1547  end if
1548  if ( trim(dtype) == "linear" ) then
1549  do j = 1, odims(2)
1550  do i = 1, odims(1)
1551  olon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * D2R
1552  enddo
1553  enddo
1554  else if ( trim(dtype) == "map" ) then
1555  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1556  olon_org(:,:) = real(gsst2D(:,:), kind=RP) * D2R
1557  endif
1558  end if
1559  case("lon_sfc")
1560  if ( .not. data_available(io_lon_sst,3) .and. data_available(io_lon_sfc,3) ) then
1561  if ( odims(1).ne.outer_nx_sfc .or. odims(2).ne.outer_ny_sfc ) then
1562  write(*,*) 'xxx namelist of "lon_sst" is not found in grads namelist!'
1563  write(*,*) 'xxx dimension is different: outer_nx_sfc and outer_nx_sst! ', outer_nx_sfc, odims(1)
1564  write(*,*) ' : outer_ny_sfc and outer_ny_sst! ', outer_ny_sfc, odims(2)
1565  call prc_mpistop
1566  end if
1567  if ( trim(dtype) == "linear" ) then
1568  do j = 1, odims(2)
1569  do i = 1, odims(1)
1570  olon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * D2R
1571  enddo
1572  enddo
1573  else if ( trim(dtype) == "map" ) then
1574  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1575  olon_org(:,:) = real(gsst2D(:,:), kind=RP) * D2R
1576  endif
1577  end if
1578  case("lon_sst")
1579  if ( .not. data_available(io_lon_sst,3) ) cycle
1580  if ( trim(dtype) == "linear" ) then
1581  do j = 1, odims(2)
1582  do i = 1, odims(1)
1583  olon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * D2R
1584  enddo
1585  enddo
1586  else if ( trim(dtype) == "map" ) then
1587  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1588  olon_org(:,:) = real(gsst2D(:,:), kind=RP) * D2R
1589  endif
1590  case("lat")
1591  if ( .not. data_available(io_lat_sfc,3) .and. .not. data_available(io_lat_sst,3) ) then
1592  if ( odims(1).ne.outer_nx .or. odims(2).ne.outer_ny ) then
1593  write(*,*) 'xxx namelist of "lat_sst" is not found in grads namelist!'
1594  write(*,*) 'xxx dimension is different: outer_nx and outer_nx_sst! ', outer_nx, odims(1)
1595  write(*,*) ' : outer_ny and outer_ny_sst! ', outer_ny, odims(2)
1596  call prc_mpistop
1597  end if
1598  if ( trim(dtype) == "linear" ) then
1599  do j = 1, odims(2)
1600  do i = 1, odims(1)
1601  olat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * D2R
1602  enddo
1603  enddo
1604  else if ( trim(dtype) == "map" ) then
1605  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1606  olat_org(:,:) = real(gsst2D(:,:), kind=RP) * D2R
1607  endif
1608  end if
1609  case("lat_sfc")
1610  if ( .not. data_available(io_lat_sst,3) .and. data_available(io_lat_sfc,3) ) then
1611  if ( odims(1).ne.outer_nx .or. odims(1).ne.outer_ny ) then
1612  write(*,*) 'xxx namelist of "lat_sst" is not found in grads namelist!'
1613  write(*,*) 'xxx dimension is different: outer_nx_sfc and outer_nx_sst! ', outer_nx_sfc, odims(1)
1614  write(*,*) ' : outer_ny_sfc and outer_ny_sst! ', outer_ny_sfc, odims(2)
1615  call prc_mpistop
1616  end if
1617  if ( trim(dtype) == "linear" ) then
1618  do j = 1, odims(2)
1619  do i = 1, odims(1)
1620  olat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * D2R
1621  enddo
1622  enddo
1623  else if ( trim(dtype) == "map" ) then
1624  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1625  olat_org(:,:) = real(gsst2D(:,:), kind=RP) * D2R
1626  endif
1627  end if
1628  case("lat_sst")
1629  if ( .not. data_available(io_lat_sst,3) ) cycle
1630  if ( trim(dtype) == "linear" ) then
1631  do j = 1, odims(2)
1632  do i = 1, odims(1)
1633  olat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * D2R
1634  enddo
1635  enddo
1636  else if ( trim(dtype) == "map" ) then
1637  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1638  olat_org(:,:) = real(gsst2D(:,:), kind=RP) * D2R
1639  endif
1640  case('SKINT')
1641  if ( .not. data_available(io_sst,3) ) then
1642  if ( odims(1).ne.outer_nx_sfc .or. odims(2).ne.outer_ny_sfc ) then
1643  write(*,*) 'xxx dimsntion is different: outer_nx_sst/outer_nx_sfc and outer_nx_sst! ', odims(1), outer_nx_sfc
1644  write(*,*) ' : outer_ny_sst/outer_ny_sfc and outer_ny_sst! ', odims(2), outer_ny_sfc
1645  call prc_mpistop
1646  end if
1647  if ( trim(dtype) == "map" ) then
1648  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,nt,item,startrec,totalrec,yrev,gsst2d)
1649  do j = 1, odims(2)
1650  do i = 1, odims(1)
1651  if ( abs(gsst2d(i,j)-missval) < eps ) then
1652  sst_org(i,j) = undef
1653  else
1654  sst_org(i,j) = real(gsst2D(i,j), kind=rp)
1655  end if
1656  enddo
1657  enddo
1658  end if
1659  endif
1660  case('SST')
1661  if ( .not. data_available(io_sst,3) ) cycle
1662  if ( trim(dtype) == "map" ) then
1663  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,nt,item,startrec,totalrec,yrev,gsst2d)
1664  do j = 1, odims(2)
1665  do i = 1, odims(1)
1666  if ( abs(gsst2d(i,j)-missval) < eps ) then
1667  sst_org(i,j) = undef
1668  else
1669  sst_org(i,j) = real(gsst2D(i,j), kind=rp)
1670  end if
1671  enddo
1672  enddo
1673  end if
1674  end select
1675  enddo loop_inputoceangrads
1676 
1677  tw_org = sst_org
1678 
1679  !do it = 1, nt
1680  ! i=int(dims(8)/2) ; j=int(dims(9)/2)
1681  ! write(*,*) "read 2D grads data",dims(8),dims(9),i,j,it
1682  ! write(*,*) "lon_org ",lon_org (i,j)
1683  ! write(*,*) "lat_org ",lat_org (i,j)
1684  ! write(*,*) "sst_org ",sst_org (i,j)
1685  ! write(*,*) "lst_org ",lst_org(i,j)
1686  ! do k=1,dims(7)
1687  ! write(*,*) "tg_org ",tg_org (k,i,j)," k= ",k
1688  ! write(*,*) "strg_org ",strg_org (k,i,j)," k= ",k
1689  ! enddo
1690  !enddo
1691 
1692  return
1693  end subroutine parentoceaninputgrads
1694 
1695  subroutine read_namelist( &
1696  grads_item, &
1697  grads_fname, &
1698  grads_dtype, &
1699  grads_swpoint, &
1700  grads_dd, &
1701  grads_lnum, &
1702  grads_lvars, &
1703  grads_startrec, &
1704  grads_totalrec, &
1705  grads_knum, &
1706  grads_yrev, &
1707  grads_fendian, &
1708  grads_missval, &
1709  data_available, &
1710  item_list, &
1711  num_item_list, &
1712  basename, &
1713  io_fid_grads_nml )
1714  implicit none
1715  character(len=H_SHORT), intent(out) :: grads_item (:)
1716  character(len=H_LONG), intent(out) :: grads_fname (:)
1717  character(len=H_LONG), intent(out) :: grads_dtype (:)
1718  real(RP), intent(out) :: grads_swpoint (:)
1719  real(RP), intent(out) :: grads_dd (:)
1720  integer, intent(out) :: grads_lnum (:)
1721  real(RP), intent(out) :: grads_lvars (:,:)
1722  integer, intent(out) :: grads_startrec(:)
1723  integer, intent(out) :: grads_totalrec(:)
1724  integer, intent(out) :: grads_knum (:)
1725  character(len=H_SHORT), intent(out) :: grads_yrev (:)
1726  character(len=H_SHORT), intent(out) :: grads_fendian (:)
1727  real(SP), intent(out) :: grads_missval (:)
1728  logical, intent(out) :: data_available(:)
1729  character(len=H_SHORT), intent(in) :: item_list (:)
1730  integer, intent(in) :: num_item_list
1731  character(len=*), intent(in) :: basename
1732  integer, intent(in) :: io_fid_grads_nml
1733 
1734  integer :: grads_vars_nmax
1735  integer :: k, n, ielem, ierr
1736 
1737  namelist /grdvar/ &
1738  item, & ! necessary
1739  dtype, & ! necessary
1740  fname, & ! necessary except for linear data
1741  swpoint, & ! for linear data
1742  dd, & ! for linear data
1743  lnum, & ! for levels data
1744  lvars, & ! for levels data
1745  startrec, & ! for map data
1746  totalrec, & ! for map data
1747  missval, & ! option
1748  knum, & ! option
1749  yrev, & ! option
1750  fendian ! option
1751 
1752  ! listup variables
1753  if ( io_fid_grads_nml > 0 ) then
1754  rewind( io_fid_grads_nml )
1755  grads_vars_nmax = 0
1756  do n = 1, grads_vars_limit
1757  read(io_fid_grads_nml, nml=grdvar, iostat=ierr)
1758  if( ierr > 0 )then
1759  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in grdvar in ', trim(basename),'. Check!'
1760  call prc_mpistop
1761  else if( ierr < 0 )then
1762  exit
1763  endif
1764  grads_vars_nmax = grads_vars_nmax + 1
1765  enddo
1766  else
1767  if( io_l ) write(io_fid_log,*) 'xxx namelist file is not open! ', trim(basename)
1768  call prc_mpistop
1769  endif
1770 
1771  if ( grads_vars_nmax > grads_vars_limit ) then
1772  if( io_l ) write(io_fid_log,*) &
1773  'xxx The number of grads vars exceeds grads_vars_limit! ',grads_vars_nmax,' >', grads_vars_limit
1774  call prc_mpistop
1775  endif
1776 
1777 
1778  ! check data availability
1779  data_available(:) = .false.
1780  do ielem = 1, num_item_list
1781  if ( io_fid_grads_nml > 0 ) rewind( io_fid_grads_nml )
1782  do n = 1, grads_vars_nmax
1783 
1784  ! set default
1785  item = ''
1786  dtype = ''
1787  fname = ''
1788  swpoint = large_number_one
1789  dd = large_number_one
1790  lnum = -99
1791  lvars = large_number_one
1792  startrec = -99
1793  totalrec = -99
1794  knum = -99
1795  yrev = 'off'
1796  fendian = 'big'
1797  missval = large_number_one
1798 
1799  ! read namelist
1800  if ( io_fid_grads_nml > 0 ) then
1801  read(io_fid_grads_nml, nml=grdvar, iostat=ierr)
1802  if( ierr /= 0 ) exit
1803  endif
1804 
1805  if(item == item_list(ielem))then
1806  grads_item(ielem) = item
1807  grads_fname(ielem) = fname
1808  grads_dtype(ielem) = dtype
1809  grads_swpoint(ielem) = swpoint
1810  grads_dd(ielem) = dd
1811  grads_lnum(ielem) = lnum
1812  do k = 1, lvars_limit
1813  grads_lvars(k,ielem) = lvars(k)
1814  enddo
1815  grads_startrec(ielem) = startrec
1816  grads_totalrec(ielem) = totalrec
1817  grads_knum(ielem) = knum
1818  grads_yrev(ielem) = yrev
1819  grads_fendian(ielem) = fendian
1820  grads_missval(ielem) = missval
1821  data_available(ielem) = .true.
1822 
1823  exit
1824  endif
1825  enddo ! n
1826  if( io_l ) write(io_fid_log,*) 'GrADS data availability ',trim(item_list(ielem)),data_available(ielem)
1827  enddo ! ielem
1828 
1829  end subroutine read_namelist
1830 
1831  !-----------------------------------------------------------------------------
1832  subroutine open_grads_file(io_fid,filename,irecl)
1833  implicit none
1834  integer,intent(in) :: io_fid
1835  character(*),intent(in) :: filename
1836  integer,intent(in) :: irecl
1837  integer :: ierr
1838 
1839  open(io_fid, &
1840  file = trim(filename), &
1841  form = 'unformatted', &
1842  access = 'direct', &
1843  recl = irecl, &
1844  status = 'old', &
1845  iostat = ierr )
1846  if ( ierr /= 0 ) then
1847  write(*,*) 'xxx grads file does not found! ', trim(filename)
1848  call prc_mpistop
1849  endif
1850  return
1851  end subroutine open_grads_file
1852 
1853  !-----------------------------------------------------------------------------
1854  subroutine read_grads_file_2d( &
1855  io_fid, &
1856  gfile, &
1857  nx,ny,nz,it, &
1858  item, &
1859  startrec, &
1860  totalrec, &
1861  yrev, &
1862  gdata )
1863  implicit none
1864  integer, intent(in) :: io_fid
1865  character(*),intent(in) :: gfile
1866  integer, intent(in) :: nx,ny,nz,it
1867  character(*),intent(in) :: item
1868  integer, intent(in) :: startrec
1869  integer, intent(in) :: totalrec
1870  character(*),intent(in) :: yrev
1871  real(SP), intent(out) :: gdata(nx,ny)
1872 
1873  real(SP) :: work(nx,ny)
1874 
1875  integer :: ierr
1876  integer :: irec, irecl
1877  integer :: i,j,k
1878 
1879  irecl=nx*ny*4
1880  call open_grads_file(io_fid, gfile, irecl)
1881  irec = totalrec * (it-1) + startrec
1882  read(io_fid, rec=irec, iostat=ierr) gdata(:,:)
1883  if ( ierr /= 0 ) then
1884  write(*,*) 'xxx grads data is not found! ',trim(item),it
1885  write(*,*) 'xxx namelist or grads data might be wrong.'
1886  call prc_mpistop
1887  endif
1888 
1889  if( trim(yrev) == "on" )then
1890  work(:,:)=gdata(:,:)
1891  do j=1,ny
1892  do i=1,nx
1893  gdata(i,j)=work(i,ny-j+1)
1894  enddo
1895  enddo
1896  endif
1897 
1898  call close_grads_file(io_fid,gfile)
1899 
1900  return
1901  end subroutine read_grads_file_2d
1902 
1903  !-----------------------------------------------------------------------------
1904  subroutine read_grads_file_3d( &
1905  io_fid, &
1906  gfile, &
1907  nx,ny,nz,it, &
1908  item, &
1909  startrec, &
1910  totalrec, &
1911  yrev, &
1912  gdata )
1913  implicit none
1914  integer,intent(in) :: io_fid
1915  character(*),intent(in) :: gfile
1916  integer, intent(in) :: nx,ny,nz,it
1917  character(*),intent(in) :: item
1918  integer, intent(in) :: startrec
1919  integer, intent(in) :: totalrec
1920  character(*),intent(in) :: yrev
1921  real(SP),intent(out) :: gdata(nx,ny,nz)
1922 
1923  real(SP) :: work(nx,ny,nz)
1924 
1925  integer :: ierr
1926  integer :: irec,irecl
1927  integer :: i,j,k
1928 
1929  irecl=nx*ny*4
1930  call open_grads_file(io_fid, gfile, irecl)
1931  do k = 1, nz
1932  irec = totalrec * (it-1) + startrec + (k-1)
1933  read(io_fid, rec=irec, iostat=ierr) gdata(:,:,k)
1934  if ( ierr /= 0 ) then
1935  write(*,*) 'xxx grads data does not found! ',trim(item),', k=',k,', it=',it,' in ', trim(gfile)
1936  call prc_mpistop
1937  endif
1938  enddo
1939 
1940  if( trim(yrev) == "on" )then
1941  work(:,:,:)=gdata(:,:,:)
1942  do k=1,nz
1943  do j=1,ny
1944  do i=1,nx
1945  gdata(i,j,k)=work(i,ny-j+1,k)
1946  enddo
1947  enddo
1948  enddo
1949  endif
1950 
1951  call close_grads_file(io_fid,gfile)
1952 
1953  return
1954  end subroutine read_grads_file_3d
1955 
1956  !-----------------------------------------------------------------------------
1957  subroutine close_grads_file(io_fid,filename)
1958  implicit none
1959  integer,intent(in) :: io_fid
1960  character(*),intent(in) :: filename
1961  integer :: ierr
1962 
1963  close(io_fid, iostat=ierr)
1964  if ( ierr /= 0 ) then
1965  write(*,*) 'xxx grads file was not closed peacefully! ',trim(filename)
1966  call prc_mpistop
1967  endif
1968 
1969  return
1970  end subroutine close_grads_file
1971 
1972 end module mod_realinput_grads
subroutine, public parentatomsetupgrads(dims, basename)
Atmos Setup.
subroutine, public parentoceaninputgrads(tw_org, sst_org, omask_org, olon_org, olat_org, basename_num, odims, nt)
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
module ATMOSPHERE / Saturation adjustment
subroutine, public parentoceanopengrads
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module STDIO
Definition: scale_stdio.F90:12
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:95
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:60
subroutine, public parentlandsetupgrads(ldims, use_waterratio, use_file_landwater, basename)
Land Setup.
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
real(rp), public const_undef
Definition: scale_const.F90:43
module TRACER
subroutine read_namelist(grads_item, grads_fname, grads_dtype, grads_swpoint, grads_dd, grads_lnum, grads_lvars, grads_startrec, grads_totalrec, grads_knum, grads_yrev, grads_fendian, grads_missval, data_available, item_list, num_item_list, basename, io_fid_grads_nml)
integer function, public io_get_available_fid()
search & get available file ID
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:93
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
module PROCESS
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:71
subroutine, public parentlandinputgrads(tg_org, strg_org, smds_org, lst_org, llon_org, llat_org, lz_org, topo_org, lmask_org, basename_num, ldims, use_file_landwater, nt)
subroutine, public log(type, message)
Definition: dc_log.f90:133
subroutine, public parentoceansetupgrads(odims, timelen, basename)
Ocean Setup.
module CONSTANT
Definition: scale_const.F90:14
integer, public prc_myrank
process num in local communicator
real(rp), public const_eps
small number
Definition: scale_const.F90:36
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
module PRECISION
module REAL input GrADS
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public rp
subroutine, public parentatominputgrads(velz_org, velx_org, vely_org, pres_org, temp_org, qtrc_org, lon_org, lat_org, cz_org, basename_num, dims, nt)
subroutine, public parentatomopengrads