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