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