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