SCALE-RM
Functions/Subroutines
mod_realinput_grads Module Reference

module REAL input GrADS More...

Functions/Subroutines

subroutine, public parentatmossetupgrads (dims, timelen, qtrc_flag, LON_all, LAT_all, basename_org, basename_num)
 Atmos Setup. More...
 
subroutine, public parentatmosinputgrads (KA_org, KS_org, KE_org, IA_org, IS_org, IE_org, JA_org, JS_org, JE_org, QA, w_org, u_org, v_org, pres_org, dens_org, pt_org, temp_org, qv_org, rh_org, qhyd_org, qtrc_org, cz_org, nopres, nodens, temp2pt, rh2qv, basename_num, sfc_diagnoses, nt)
 
subroutine, public parentlandsetupgrads (ldims, timelen, lon_all, lat_all, basename, basename_num)
 Land Setup. More...
 
subroutine, public parentlandinputgrads (KA_org, KS_org, KE_org, IA_org, IS_org, IE_org, JA_org, JS_org, JE_org, tg_org, strg_org, smds_org, lst_org, lz_org, topo_org, lmask_org, use_waterratio, ldims, basename_num, use_file_landwater, nt)
 
subroutine, public parentoceansetupgrads (odims, timelen, lon_all, lat_all, basename, basename_num)
 Ocean Setup. More...
 
subroutine, public parentoceaninputgrads (IA_org, IS_org, IE_org, JA_org, JS_org, JE_org, tw_org, sst_org, omask_org, basename_num, odims, nt)
 
subroutine read2d (start, count, data, name, fid, postfix, exist, oneD, step)
 
subroutine read3d (start, count, data, name, fid, postfix, exist, step)
 

Detailed Description

module REAL input GrADS

Description
read data from GrADS file for real atmospheric simulations
Author
Team SCALE

Function/Subroutine Documentation

◆ parentatmossetupgrads()

subroutine, public mod_realinput_grads::parentatmossetupgrads ( integer, dimension(6), intent(out)  dims,
integer, intent(out)  timelen,
logical, dimension(qa), intent(out)  qtrc_flag,
real(rp), dimension(:,:), intent(out), allocatable  LON_all,
real(rp), dimension(:,:), intent(out), allocatable  LAT_all,
character(len=*), intent(in)  basename_org,
character(len=*), intent(in)  basename_num 
)

Atmos Setup.

Definition at line 62 of file mod_realinput_grads.F90.

62  use scale_const, only: &
63  d2r => const_d2r
64  use scale_file_grads, only: &
67  file_grads_get_shape
68  use mod_atmos_phy_mp_vars, only: &
69  qs_mp, &
70  qe_mp
71  implicit none
72  integer, intent(out) :: dims(6)
73  integer, intent(out) :: timelen
74  logical, intent(out) :: qtrc_flag(QA)
75  real(RP), allocatable, intent(out) :: LON_all(:,:)
76  real(RP), allocatable, intent(out) :: LAT_all(:,:)
77 
78  character(len=*), intent(in) :: basename_org
79  character(len=*), intent(in) :: basename_num
80 
81 ! namelist / PARAM_MKINIT_REAL_GrADS /
82 
83  integer :: var_id
84  real(RP), allocatable :: lon1d(:), lat1d(:)
85 
86 
87  integer :: i, j, iq
88  integer :: ierr
89  !---------------------------------------------------------------------------
90 
91  log_newline
92  log_info("ParentAtmosSetupGrADS",*) 'Setup'
93 
94 !!$ !--- read namelist
95 !!$ rewind(IO_FID_CONF)
96 !!$ read(IO_FID_CONF,nml=PARAM_MKINIT_REAL_GrADS,iostat=ierr)
97 !!$
98 !!$ if( ierr > 0 ) then
99 !!$ LOG_ERROR("ParentAtmosSetupGrADS",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_GrADS. Check!'
100 !!$ call PRC_abort
101 !!$ endif
102 !!$ LOG_NML(PARAM_MKINIT_REAL_GrADS)
103 
104 
105  if ( basename_org == "" ) then
106  log_error("ParentAtmosSetupGrADS",*) '"BASENAME_ORG" is not specified in "PARAM_MKINIT_REAL_ATMOS"!', trim(basename_org)
107  call prc_abort
108  endif
109 
110  call file_grads_open( basename_org, & ! (in)
111  file_id_atm ) ! (out)
112 
113 #ifdef QUICKDEBUG
114  qtrc_flag(:) = .false.
115 #endif
116  do iq = 1, qa
117  if ( iq >= qs_mp .and. iq <= qe_mp ) cycle
118  call file_grads_varid( file_id_atm, tracer_name(iq), var_id )
119  qtrc_flag(iq) = var_id > 0
120  end do
121 
122  call file_grads_get_shape( file_id_atm, "U", & ! (in)
123  dims(:3) ) ! (out)
124 
125  ! half level
126  dims(4) = dims(1)
127  dims(5) = dims(2)
128  dims(6) = dims(3)
129 
130 
131  allocate( lon_all(dims(2), dims(3)) )
132  allocate( lat_all(dims(2), dims(3)) )
133 
134 
135  call read2d( (/1,1/), dims(2:3), lon_all(:,:), "lon", file_id_atm, basename_num, oned=1 )
136  lon_all(:,:) = lon_all(:,:) * d2r
137 
138  call read2d( (/1,1/), dims(2:3), lat_all(:,:), "lat", file_id_atm, basename_num, oned=2 )
139  lat_all(:,:) = lat_all(:,:) * d2r
140 
141  ! tentative
142  timelen = 1
143 
144  return

References scale_const::const_d2r, scale_file_grads::file_grads_open(), scale_file_grads::file_grads_varid(), scale_prc::prc_abort(), scale_tracer::qa, mod_atmos_phy_mp_vars::qe_mp, mod_atmos_phy_mp_vars::qs_mp, read2d(), and scale_tracer::tracer_name.

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentatmosinputgrads()

subroutine, public mod_realinput_grads::parentatmosinputgrads ( integer, intent(in)  KA_org,
integer, intent(in)  KS_org,
integer, intent(in)  KE_org,
integer, intent(in)  IA_org,
integer, intent(in)  IS_org,
integer, intent(in)  IE_org,
integer, intent(in)  JA_org,
integer, intent(in)  JS_org,
integer, intent(in)  JE_org,
integer, intent(in)  QA,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  w_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  u_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  v_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  pres_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  dens_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  pt_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  temp_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  qv_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  rh_org,
real(rp), dimension(ka_org,ia_org,ja_org,n_hyd), intent(out)  qhyd_org,
real(rp), dimension(ka_org,ia_org,ja_org,qa), intent(out)  qtrc_org,
real(rp), dimension(ka_org,ia_org,ja_org), intent(out)  cz_org,
logical, intent(out)  nopres,
logical, intent(out)  nodens,
logical, intent(out)  temp2pt,
logical, intent(out)  rh2qv,
character(len=*), intent(in)  basename_num,
logical, intent(in)  sfc_diagnoses,
integer, intent(in)  nt 
)

Definition at line 170 of file mod_realinput_grads.F90.

170  use scale_const, only: &
171  undef => const_undef, &
172  d2r => const_d2r, &
173  eps => const_eps, &
174  epsvap => const_epsvap, &
175  epstvap => const_epstvap, &
176  grav => const_grav, &
177  laps => const_laps, &
178  p00 => const_pre00, &
179  rdry => const_rdry, &
180  cpdry => const_cpdry
181  use scale_atmos_hydrometeor, only: &
182  hyd_name, &
183  n_hyd
184  use mod_atmos_phy_mp_vars, only: &
185  qs_mp, &
186  qe_mp
187  implicit none
188  integer, intent(in) :: KA_org
189  integer, intent(in) :: KS_org
190  integer, intent(in) :: KE_org
191  integer, intent(in) :: IA_org
192  integer, intent(in) :: IS_org
193  integer, intent(in) :: IE_org
194  integer, intent(in) :: JA_org
195  integer, intent(in) :: JS_org
196  integer, intent(in) :: JE_org
197  integer, intent(in) :: QA
198  real(RP), intent(out) :: w_org(KA_org,IA_org,JA_org)
199  real(RP), intent(out) :: u_org(KA_org,IA_org,JA_org)
200  real(RP), intent(out) :: v_org(KA_org,IA_org,JA_org)
201  real(RP), intent(out) :: pres_org(KA_org,IA_org,JA_org)
202  real(RP), intent(out) :: dens_org(KA_org,IA_org,JA_org)
203  real(RP), intent(out) :: pt_org(KA_org,IA_org,JA_org)
204  real(RP), intent(out) :: temp_org(KA_org,IA_org,JA_org)
205  real(RP), intent(out) :: qv_org(KA_org,IA_org,JA_org)
206  real(RP), intent(out) :: rh_org(KA_org,IA_org,JA_org)
207  real(RP), intent(out) :: qhyd_org(KA_org,IA_org,JA_org,N_HYD)
208  real(RP), intent(out) :: qtrc_org(KA_org,IA_org,JA_org,QA)
209  real(RP), intent(out) :: cz_org(KA_org,IA_org,JA_org)
210  logical, intent(out) :: nopres
211  logical, intent(out) :: nodens
212  logical, intent(out) :: temp2pt
213  logical, intent(out) :: rh2qv
214  character(len=*), intent(in) :: basename_num
215  logical, intent(in) :: sfc_diagnoses
216  integer, intent(in) :: nt
217 
218  character(len=H_SHORT) :: item
219 
220  integer :: lm_layer(IA_org,JA_org)
221 
222  real(RP) :: work(KA_org-2,IA_org,JA_org)
223  real(RP) :: work2d(IA_org,JA_org)
224  real(RP) :: coef
225  integer :: start(3), count(3)
226 
227  logical :: exist
228 
229  integer :: i, j, k, iq
230  !---------------------------------------------------------------------------
231 
232  start(:) = (/1,is_org,js_org/)
233  count(:) = (/ka_org-2,ia_org,ja_org/)
234 
235  ! pressure
236  nopres = .false.
237  call read3d( start(:), count(:), work(:,:,:), "pressure", file_id_atm, basename_num, exist=exist, step=nt )
238  if ( .not. exist ) then
239  call read3d( start(:), count(:), work(:,:,:), "plev", file_id_atm, basename_num, exist=exist, step=nt )
240  if ( .not. exist ) nopres = .true.
241  end if
242  if ( .not. nopres ) then
243  !$omp parallel do collapse(2)
244  do j = 1, ja_org
245  do i = 1, ia_org
246  do k = 1, ka_org-2
247  pres_org(k+2,i,j) = work(k,i,j)
248  enddo
249  enddo
250  enddo
251  end if
252 
253 
254  ! dens
255  call read3d( start(:), count(:), work(:,:,:), "DENS", file_id_atm, basename_num, exist=exist, step=nt )
256  if ( exist ) then
257  !$omp parallel do collapse(2)
258  do j = 1, ja_org
259  do i = 1, ia_org
260  do k = 1, ka_org-2
261  dens_org(k+2,i,j) = work(k,i,j)
262  enddo
263  enddo
264  enddo
265  nodens = .false.
266  else
267  nodens = .true.
268  end if
269 
270  ! W
271  call read3d( start(:), count(:), work(:,:,:), "W", file_id_atm, basename_num, exist=exist, step=nt )
272  if ( exist ) then
273  !$omp parallel do collapse(2)
274  do j = 1, ja_org
275  do i = 1, ia_org
276  do k = 1, ka_org-2
277  w_org(k+2,i,j) = work(k,i,j)
278  enddo
279  enddo
280  enddo
281  else
282  !$omp parallel do collapse(2)
283  do j = 1, ja_org
284  do i = 1, ia_org
285  do k = 1, ka_org-2
286  w_org(k+2,i,j) = 0.0_rp
287  enddo
288  enddo
289  enddo
290  end if
291 
292  ! U
293  call read3d( start(:), count(:), work(:,:,:), "U", file_id_atm, basename_num, exist=exist, step=nt )
294  if ( exist ) then
295  !$omp parallel do collapse(2)
296  do j = 1, ja_org
297  do i = 1, ia_org
298  do k = 1, ka_org-2
299  u_org(k+2,i,j) = work(k,i,j)
300  enddo
301  enddo
302  enddo
303  else
304  log_error("ParentAtmosInputGrADS",*) '"U" is requierd'
305  call prc_abort
306  end if
307 
308  ! V
309  call read3d( start(:), count(:), work(:,:,:), "V", file_id_atm, basename_num, exist=exist, step=nt )
310  if ( exist ) then
311  !$omp parallel do collapse(2)
312  do j = 1, ja_org
313  do i = 1, ia_org
314  do k = 1, ka_org-2
315  v_org(k+2,i,j) = work(k,i,j)
316  enddo
317  enddo
318  enddo
319  else
320  log_error("ParentAtmosInputGrADS",*) '"V" is requierd'
321  call prc_abort
322  end if
323 
324  ! T
325  call read3d( start(:), count(:), work(:,:,:), "PT", file_id_atm, basename_num, exist=exist, step=nt )
326  if ( exist ) then
327  !$omp parallel do collapse(2)
328  do j = 1, ja_org
329  do i = 1, ia_org
330  do k = 1, ka_org-2
331  pt_org(k+2,i,j) = work(k,i,j)
332  enddo
333  enddo
334  enddo
335  temp2pt = .false.
336  else
337  call read3d( start(:), count(:), work(:,:,:), "T", file_id_atm, basename_num, exist=exist, step=nt )
338  if ( exist ) then
339  !$omp parallel do collapse(2)
340  do j = 1, ja_org
341  do i = 1, ia_org
342  do k = 1, ka_org-2
343  temp_org(k+2,i,j) = work(k,i,j)
344  enddo
345  enddo
346  enddo
347  temp2pt = .true.
348  else
349  log_error("ParentAtmosInputGrADS",*) '"PT" or "T" is requierd'
350  call prc_abort
351  end if
352  end if
353 
354  ! height
355  coef = 1.0_rp
356  call read3d( start(:), count(:), work(:,:,:), "height", file_id_atm, basename_num, exist=exist, step=nt )
357  if ( .not. exist ) then
358  call read3d( start(:), count(:), work(:,:,:), "HGT", file_id_atm, basename_num, exist=exist, step=nt )
359  if ( .not. exist ) then
360  ! geopotential to height
361  coef = 1.0_rp / grav
362  call read3d( start(:), count(:), work(:,:,:), "GP", file_id_atm, basename_num, exist=exist, step=nt )
363  if ( .not. exist ) then
364  log_error("ParentAtmosInputGrADS",*) '"height", "HGT", or "GP" is requierd'
365  call prc_abort
366  end if
367  end if
368  end if
369  !$omp parallel do collapse(2)
370  do j = 1, ja_org
371  do i = 1, ia_org
372  do k = 1, ka_org-2
373  cz_org(k+2,i,j) = work(k,i,j) * coef
374  enddo
375  enddo
376  enddo
377 
378  ! QV
379  call read3d( start(:), count(:), work(:,:,:), "QV", file_id_atm, basename_num, exist=exist, step=nt )
380  if ( exist ) then
381  !$omp parallel do collapse(2)
382  do j = 1, ja_org
383  do i = 1, ia_org
384  do k = 1, ka_org-2
385  qv_org(k+2,i,j) = work(k,i,j)
386  enddo
387  enddo
388  enddo
389  rh2qv = .false.
390  else
391  call read3d( start(:), count(:), work(:,:,:), "RH", file_id_atm, basename_num, exist=exist, step=nt )
392  if ( exist ) then
393  !$omp parallel do collapse(2)
394  do j = 1, ja_org
395  do i = 1, ia_org
396  do k = 1, ka_org-2
397  rh_org(k+2,i,j) = work(k,i,j)
398  enddo
399  enddo
400  enddo
401  rh2qv = .true.
402  else
403  !$omp parallel do collapse(2)
404  do j = 1, ja_org
405  do i = 1, ia_org
406  do k = 1, ka_org-2
407  qv_org(k+2,i,j) = 0.0_rp
408  enddo
409  enddo
410  enddo
411  rh2qv = .false.
412  end if
413  end if
414 
415  ! QC, QR, QI, QS, QG, QH
416  do iq = 1, n_hyd
417  call read3d( start(:), count(:), work(:,:,:), hyd_name(iq), file_id_atm, basename_num, exist=exist, step=nt )
418  if ( exist ) then
419  !$omp parallel do collapse(2)
420  do j = 1, ja_org
421  do i = 1, ia_org
422  do k = 1, ka_org-2
423  qhyd_org(k+2,i,j,iq) = work(k,i,j)
424  enddo
425  enddo
426  enddo
427  else
428  !$omp parallel do collapse(2)
429  do j = 1, ja_org
430  do i = 1, ia_org
431  do k = 1, ka_org-2
432  qhyd_org(k+2,i,j,iq) = 0.0_rp
433  enddo
434  enddo
435  enddo
436  end if
437 #ifdef QUICKDEBUG
438  !$omp parallel do collapse(2)
439  do j = 1, ja_org
440  do i = 1, ia_org
441  do k = 1, 2
442  qhyd_org(k,i,j,iq) = undef
443  enddo
444  enddo
445  enddo
446 #endif
447  end do
448 
449  ! QTRC
450  do iq = 1, qa
451  if ( iq >= qs_mp .and. iq <= qe_mp ) cycle
452  call read3d( start(:), count(:), work(:,:,:), tracer_name(iq), file_id_atm, basename_num, exist=exist, step=nt )
453  if ( exist ) then
454  !$omp parallel do collapse(2)
455  do j = 1, ja_org
456  do i = 1, ia_org
457  do k = 1, ka_org-2
458  qtrc_org(k+2,i,j,iq) = work(k,i,j)
459  enddo
460  enddo
461  enddo
462  else
463  !$omp parallel do collapse(2)
464  do j = 1, ja_org
465  do i = 1, ia_org
466  do k = 1, ka_org-2
467  qtrc_org(k+2,i,j,iq) = undef
468  enddo
469  enddo
470  enddo
471  end if
472  end do
473 
474 
475  if ( sfc_diagnoses ) then
476 
477  ! MSLP
478  call read2d( start(2:), count(2:), work2d(:,:), "MSLP", file_id_atm, basename_num, exist=exist, step=nt )
479  if ( exist ) then
480  !$omp parallel do
481  do j = 1, ja_org
482  do i = 1, ia_org
483  pres_org(1,i,j) = work2d(i,j)
484  enddo
485  enddo
486  else
487  !$omp parallel do
488  do j = 1, ja_org
489  do i = 1, ia_org
490  pres_org(1,i,j) = undef
491  enddo
492  enddo
493  end if
494 
495  ! SFC_PRES, PSFC
496  call read2d( start(2:), count(2:), work2d(:,:), "SFC_PRES", file_id_atm, basename_num, exist=exist, step=nt )
497  if ( .not. exist ) then
498  call read2d( start(2:), count(2:), work2d(:,:), "PSFC", file_id_atm, basename_num, exist=exist, step=nt )
499  end if
500  if ( exist ) then
501  !$omp parallel do
502  do j = 1, ja_org
503  do i = 1, ia_org
504  pres_org(2,i,j) = work2d(i,j)
505  enddo
506  enddo
507  else
508  !$omp parallel do
509  do j = 1, ja_org
510  do i = 1, ia_org
511  pres_org(2,i,j) = undef
512  enddo
513  enddo
514  end if
515 
516  ! U10
517  call read2d( start(2:), count(2:), work2d(:,:), "U10", file_id_atm, basename_num, exist=exist, step=nt )
518  if ( exist ) then
519  !$omp parallel do
520  do j = 1, ja_org
521  do i = 1, ia_org
522  u_org(2,i,j) = work2d(i,j)
523  enddo
524  enddo
525  else
526  !$omp parallel do
527  do j = 1, ja_org
528  do i = 1, ia_org
529  u_org(2,i,j) = undef
530  enddo
531  enddo
532  end if
533 
534  ! V10
535  call read2d( start(2:), count(2:), work2d(:,:), "V10", file_id_atm, basename_num, exist=exist, step=nt )
536  if ( exist ) then
537  !$omp parallel do
538  do j = 1, ja_org
539  do i = 1, ia_org
540  v_org(2,i,j) = work2d(i,j)
541  enddo
542  enddo
543  else
544  !$omp parallel do
545  do j = 1, ja_org
546  do i = 1, ia_org
547  v_org(2,i,j) = undef
548  enddo
549  enddo
550  end if
551 
552  ! T2
553  call read2d( start(2:), count(2:), work2d(:,:), "T2", file_id_atm, basename_num, exist=exist, step=nt )
554  if ( exist ) then
555  !$omp parallel do
556  do j = 1, ja_org
557  do i = 1, ia_org
558  temp_org(2,i,j) = work2d(i,j)
559  enddo
560  enddo
561  else
562  !$omp parallel do
563  do j = 1, ja_org
564  do i = 1, ia_org
565  temp_org(2,i,j) = undef
566  enddo
567  enddo
568  end if
569 
570  ! Q2, RH2
571  if ( rh2qv ) then
572  call read2d( start(2:), count(2:), work2d(:,:), "RH2", file_id_atm, basename_num, exist=exist, step=nt )
573  if ( exist ) then
574  !$omp parallel do
575  do j = 1, ja_org
576  do i = 1, ia_org
577  rh_org(2,i,j) = work2d(i,j)
578  enddo
579  enddo
580  else
581  !$omp parallel do
582  do j = 1, ja_org
583  do i = 1, ia_org
584  rh_org(2,i,j) = undef
585  enddo
586  enddo
587  end if
588  else
589  call read2d( start(2:), count(2:), work2d(:,:), "Q2", file_id_atm, basename_num, exist=exist, step=nt )
590  if ( exist ) then
591  !$omp parallel do
592  do j = 1, ja_org
593  do i = 1, ia_org
594  qv_org(2,i,j) = work2d(i,j)
595  enddo
596  enddo
597  else
598  !$omp parallel do
599  do j = 1, ja_org
600  do i = 1, ia_org
601  qv_org(2,i,j) = undef
602  enddo
603  enddo
604  end if
605 
606  end if
607 
608  ! topo
609  coef = 1.0_rp
610  call read2d( start(2:), count(2:), work2d(:,:), "topo", file_id_atm, basename_num, exist=exist, step=nt )
611  if ( .not. exist ) then
612  ! surface geopotential
613  coef = 1.0_rp / grav
614  call read2d( start(2:), count(2:), work2d(:,:), "SGP", file_id_atm, basename_num, exist=exist, step=nt )
615  end if
616  if ( exist ) then
617  !$omp parallel do
618  do j = 1, ja_org
619  do i = 1, ia_org
620  cz_org(2,i,j) = work2d(i,j) * coef
621  enddo
622  enddo
623  else
624  !$omp parallel do
625  do j = 1, ja_org
626  do i = 1, ia_org
627  cz_org(2,i,j) = undef
628  enddo
629  enddo
630  end if
631 
632  end if
633 
634  return

References scale_const::const_cpdry, scale_const::const_d2r, scale_const::const_eps, scale_const::const_epstvap, scale_const::const_epsvap, scale_const::const_grav, scale_const::const_laps, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_undef, scale_atmos_hydrometeor::hyd_name, scale_tracer::k, scale_atmos_hydrometeor::n_hyd, scale_prc::prc_abort(), scale_tracer::qa, mod_atmos_phy_mp_vars::qe_mp, mod_atmos_phy_mp_vars::qs_mp, read2d(), read3d(), and scale_tracer::tracer_name.

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentlandsetupgrads()

subroutine, public mod_realinput_grads::parentlandsetupgrads ( integer, dimension(3), intent(out)  ldims,
integer, intent(out)  timelen,
real(rp), dimension(:,:), intent(out), allocatable  lon_all,
real(rp), dimension(:,:), intent(out), allocatable  lat_all,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  basename_num 
)

Land Setup.

Definition at line 645 of file mod_realinput_grads.F90.

645  use scale_const, only: &
646  d2r => const_d2r
647  use scale_file_grads, only: &
648  file_grads_open, &
649  file_grads_get_shape, &
652  implicit none
653  integer, intent(out) :: ldims(3)
654  integer, intent(out) :: timelen
655  real(RP), allocatable, intent(out) :: lon_all(:,:)
656  real(RP), allocatable, intent(out) :: lat_all(:,:)
657 
658  character(len=*), intent(in) :: basename
659  character(len=*), intent(in) :: basename_num
660 
661  character(len=7) :: vname
662  integer :: vid
663  integer :: shape(2)
664  logical :: exist
665  !---------------------------------------------------------------------------
666 
667  log_info("ParentLandSetupGrADS",*) 'Real Case/Land Input File Type: GrADS format'
668 
669 
670  if ( basename == "" ) then
671  log_error("ParentLandSetupGrADS",*) '"BASENAME_ORG" is not specified in "PARAM_MKINIT_REAL_LAND"!', trim(basename)
672  call prc_abort
673  endif
674 
675  ! open
676  call file_grads_open( basename, & ! (in)
677  file_id_lnd ) ! (out)
678 
679  ! get shape
680  call file_grads_varid( file_id_lnd, "LAND_SFC_TEMP", & ! (in)
681  vid ) ! (out)
682  if ( vid < 0 ) then
683  call file_grads_varid( file_id_lnd, "LAND_TEMP", & ! (in)
684  vid ) ! (out)
685  end if
686  if ( vid < 0 ) then
687  call file_grads_varid( file_id_lnd, "SFC_TEMP", & ! (in)
688  vid ) ! (out)
689  end if
690  if ( vid < 0 ) then
691  call file_grads_varid( file_id_lnd, "STEMP", & ! (in)
692  vid ) ! (out)
693  end if
694  if ( vid < 0 ) then
695  log_error("ParentLandSetupGrADS",*) '"LAND_SFC_TEMP", "LAND_TEMP", "SFC_TEMP" or "STEMP" is necessary'
696  call prc_abort
697  end if
698  call file_grads_get_shape( file_id_lnd, vid, & ! (in)
699  ldims(:) ) ! (out)
700 
701 
702  ! get lon, lat
703  allocate( lon_all(ldims(2), ldims(3)) )
704  allocate( lat_all(ldims(2), ldims(3)) )
705 
706 
707  call file_grads_varid( file_id_lnd, "lon_sfc", & ! (in)
708  vid ) ! (out)
709  if ( vid > 0 ) then
710  vname = "lon_sfc"
711  else
712  call file_grads_varid( file_id_lnd, "lon", & ! (in)
713  vid ) ! (out)
714  if ( vid < 0 ) then
715  log_error("ParentLandSetupGrADS",*) '"lon_sfc" or "lon" is necessary'
716  call prc_abort
717  end if
718  call file_grads_get_shape( file_id_lnd, vid, & ! (in)
719  shape(:) ) ! (out)
720  if ( file_grads_isoned( file_id_lnd, vid ) ) then
721  if ( ldims(2) .ne. shape(1) .and. shape(1) .ne. -1 ) then
722  log_error("ParentLandSetupGrADS",*) 'dimension of "lon" is different! ', ldims(2), shape(1)
723  call prc_abort
724  end if
725  else
726  if ( ldims(2) .ne. shape(1) .or. ldims(3) .ne. shape(2) ) then
727  log_error("ParentLandSetupGrADS",*) 'dimension of "lon" is different! ', ldims(2), shape(1), ldims(3), shape(2)
728  call prc_abort
729  end if
730  end if
731  vname = "lon"
732  end if
733  call read2d( (/1,1/), ldims(2:3), lon_all(:,:), vname, file_id_lnd, basename_num, oned=1 )
734  lon_all(:,:) = lon_all(:,:) * d2r
735 
736 
737  call file_grads_varid( file_id_lnd, "lat_sfc", & ! (in)
738  vid ) ! (out)
739  if ( vid > 0 ) then
740  vname = "lat_sfc"
741  else
742  call file_grads_varid( file_id_lnd, "lat", & ! (in)
743  vid ) ! (out)
744  if ( vid < 0 ) then
745  log_error("ParentLandSetupGrADS",*) '"lat_sfc" or "lat" is necessary'
746  call prc_abort
747  end if
748  call file_grads_get_shape( file_id_lnd, vid, & ! (in)
749  shape(:) ) ! (out)
750  if ( file_grads_isoned( file_id_lnd, vid ) ) then
751  if ( ldims(3) .ne. shape(1) .and. shape(1) .ne. -1 ) then
752  log_error("ParentLandSetupGrADS",*) 'dimension of "lat" is different! ', ldims(2), shape(1)
753  call prc_abort
754  end if
755  else
756  if ( ldims(2) .ne. shape(1) .or. ldims(3) .ne. shape(2) ) then
757  log_error("ParentLandSetupGrADS",*) 'dimension of "lat" is different! ', ldims(2), shape(1), ldims(3), shape(2)
758  call prc_abort
759  end if
760  end if
761  vname = "lat"
762  end if
763  call read2d( (/1,1/), ldims(2:3), lat_all(:,:), vname, file_id_atm, basename_num, oned=2 )
764  lat_all(:,:) = lat_all(:,:) * d2r
765 
766 
767  ! tentative
768  timelen = 1
769 
770 
771  return

References scale_const::const_d2r, scale_file_grads::file_grads_isoned(), scale_file_grads::file_grads_open(), scale_file_grads::file_grads_varid(), scale_prc::prc_abort(), and read2d().

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentlandinputgrads()

subroutine, public mod_realinput_grads::parentlandinputgrads ( integer, intent(in)  KA_org,
integer, intent(in)  KS_org,
integer, intent(in)  KE_org,
integer, intent(in)  IA_org,
integer, intent(in)  IS_org,
integer, intent(in)  IE_org,
integer, intent(in)  JA_org,
integer, intent(in)  JS_org,
integer, intent(in)  JE_org,
real(rp), dimension (ka_org,ia_org,ja_org), intent(out)  tg_org,
real(rp), dimension (ka_org,ia_org,ja_org), intent(out)  strg_org,
real(rp), dimension (ka_org,ia_org,ja_org), intent(out)  smds_org,
real(rp), dimension (ia_org,ja_org), intent(out)  lst_org,
real(rp), dimension (ka_org), intent(out)  lz_org,
real(rp), dimension (ia_org,ja_org), intent(out)  topo_org,
real(rp), dimension(ia_org,ja_org), intent(out)  lmask_org,
logical, intent(out)  use_waterratio,
integer, dimension(3), intent(in)  ldims,
character(len=*), intent(in)  basename_num,
logical, intent(in)  use_file_landwater,
integer, intent(in)  nt 
)

Definition at line 790 of file mod_realinput_grads.F90.

790  use scale_prc, only: &
791  prc_abort
792  use scale_const, only: &
793  undef => const_undef, &
794  d2r => const_d2r, &
795  grav => const_grav, &
796  tem00 => const_tem00, &
797  eps => const_eps
798  use scale_file_grads, only: &
800  file_grads_get_shape, &
801  file_grads_read
802  implicit none
803  integer, intent(in) :: KA_org, KS_org, KE_org
804  integer, intent(in) :: IA_org, IS_org, IE_org
805  integer, intent(in) :: JA_org, JS_org, JE_org
806 
807  real(RP), intent(out) :: tg_org (KA_org,IA_org,JA_org)
808  real(RP), intent(out) :: strg_org (KA_org,IA_org,JA_org)
809  real(RP), intent(out) :: smds_org (KA_org,IA_org,JA_org)
810  real(RP), intent(out) :: lst_org (IA_org,JA_org)
811  real(RP), intent(out) :: lz_org (KA_org)
812  real(RP), intent(out) :: topo_org (IA_org,JA_org)
813  real(RP), intent(out) :: lmask_org(IA_org,JA_org)
814  logical, intent(out) :: use_waterratio
815 
816  integer, intent(in) :: ldims(3)
817  character(len=*), intent(in) :: basename_num
818  logical, intent(in) :: use_file_landwater ! use land water data from files
819  integer, intent(in) :: nt
820 
821  integer :: start(3), count(3), shape(2)
822  logical :: exist
823 
824  integer :: vid
825  integer :: i, j, k
826  !---------------------------------------------------------------------------
827 
828  start(:) = (/ks_org,is_org,js_org/)
829  count(:) = (/ka_org,ia_org,ja_org/)
830 
831  ! lsmask
832  call read2d( start(2:), count(2:), lmask_org(:,:), "lsmask", file_id_lnd, basename_num, exist=exist, step=nt )
833  if ( .not. exist ) then
834  !$omp parallel do
835  do j = 1, ja_org
836  do i = 1, ia_org
837  lmask_org(i,j) = undef
838  end do
839  end do
840  end if
841 
842  ! llev
843  call file_grads_get_shape( file_id_lnd, "llev", & ! (in)
844  shape(:) ) ! (out)
845  if ( ldims(1) .ne. shape(1) )then
846  log_error("ParentLandInputGrADS",*) '"nz" must be equal to nz of "STEMP" for llev. :', shape(1), ldims(1)
847  call prc_abort
848  endif
849  call file_grads_read( file_id_lnd, "llev", & ! (in)
850  lz_org(:) ) ! (out)
851 
852  ! LAND_TEMP, STEMP
853  call read3d( start(:), count(:), tg_org(:,:,:), "LAND_TEMP", file_id_lnd, basename_num, exist=exist, step=nt )
854  if ( .not. exist ) then
855  call read3d( start(:), count(:), tg_org(:,:,:), "STEMP", file_id_lnd, basename_num, exist=exist )
856  end if
857  if ( .not. exist ) then
858  log_error("ParentAtmosInputGrADS",*) '"LAND_TEMP" or "STEMP" is necessary'
859  call prc_abort
860  endif
861 
862  if ( use_file_landwater ) then
863  ! LAND_WATER, SMOISVC, SMOISDS
864  call read3d( start(:), count(:), strg_org(:,:,:), "LAND_WATER", file_id_lnd, basename_num, exist=exist, step=nt )
865  if ( .not. exist ) then
866  call read3d( start(:), count(:), strg_org(:,:,:), "SMOISVC", file_id_lnd, basename_num, exist=exist, step=nt )
867  end if
868  if ( exist ) then
869  !$omp parallel do collapse(2)
870  do j = 1, ja_org
871  do i = 1, ia_org
872  do k = 2, ka_org
873  if ( strg_org(k,i,j) == undef ) strg_org(k,i,j) = strg_org(k-1,i,j)
874  end do
875  end do
876  end do
877  use_waterratio = .false.
878  else
879  call read3d( start(:), count(:), smds_org(:,:,:), "SMOISDS", file_id_lnd, basename_num, exist=exist, step=nt )
880  if ( exist ) then
881  !$omp parallel do collapse(2)
882  do j = 1, ja_org
883  do i = 1, ia_org
884  do k = 2, ka_org
885  if ( smds_org(k,i,j) == undef ) smds_org(k,i,j) = smds_org(k-1,i,j)
886  end do
887  end do
888  end do
889  use_waterratio = .true.
890  else
891  log_error("ParentAtmosInputGrADS",*) '"LAND_WATER", "SMOISVC", or "SMOISDS" is necessary'
892  call prc_abort
893  end if
894  endif
895 
896  end if ! use_file_landwater
897 
898 
899  ! LAND_SFC_TEMP, SFC_TEMP, SKINT
900  call read2d( start(2:), count(2:), lst_org(:,:), "LAND_SFC_TEMP", file_id_lnd, basename_num, exist=exist, step=nt )
901  if ( .not. exist ) then
902  call read2d( start(2:), count(2:), lst_org(:,:), "SFC_TEMP", file_id_lnd, basename_num, exist=exist )
903  end if
904  if ( .not. exist ) then
905  call read2d( start(2:), count(2:), lst_org(:,:), "SKINT", file_id_lnd, basename_num, exist=exist )
906  end if
907  if ( .not. exist ) then
908  log_error("ParentAtmosInputGrADS",*) '"LAND_SFC_TEMP", "SFC_TEMP", or "SKINT" is necessary'
909  call prc_abort
910  endif
911 
912 
913  ! topo_sfc, topo, SGP
914 
915  call read2d( start(2:), count(2:), topo_org(:,:), "topo_sfc", file_id_lnd, basename_num, exist=exist, step=nt )
916  if ( .not. exist ) then
917  call file_grads_varcheck( file_id_lnd, "topo", & ! (in)
918  exist ) ! (out)
919  if ( exist ) then
920  call file_grads_get_shape( file_id_lnd, "topo", & ! (in)
921  shape(:) ) ! (out)
922  if ( ldims(2).ne.shape(1) .or. ldims(3).ne.shape(2) ) then
923  log_warn_cont(*) 'dimension of "topo" is different! ', ldims(2), shape(1), ldims(3), shape(2)
924  exist = .false.
925  else
926  call read2d( start(2:), count(2:), topo_org(:,:), "topo", file_id_lnd, basename_num, exist=exist )
927  end if
928  end if
929  end if
930  if ( .not. exist ) then
931  call file_grads_varcheck( file_id_lnd, "SGP", & ! (in)
932  exist ) ! (out)
933  if ( exist ) then
934  call file_grads_get_shape( file_id_lnd, "SGP", & ! (in)
935  shape(:) ) ! (out)
936  if ( ldims(2).ne.shape(1) .or. ldims(3).ne.shape(2) ) then
937  log_warn_cont(*) 'dimension of "SGP" is different! ', ldims(2), shape(1), ldims(3), shape(2)
938  exist = .false.
939  else
940  call read2d( start(2:), count(2:), topo_org(:,:), "SGP", file_id_lnd, basename_num, exist=exist )
941  !$omp parallel do
942  do j = 1, ja_org
943  do i = 1, ia_org
944  topo_org(i,j) = topo_org(i,j) / grav
945  end do
946  end do
947  end if
948  end if
949  end if
950  if ( .not. exist ) then
951  log_warn("ParentLandInputGrADS",*) '"topo_sfc", "topo", or "SGP" is not found in grads namelist'
952  end if
953 
954  return

References scale_const::const_d2r, scale_const::const_eps, scale_const::const_grav, scale_const::const_tem00, scale_const::const_undef, scale_file_grads::file_grads_varcheck(), scale_tracer::k, scale_prc::prc_abort(), read2d(), and read3d().

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentoceansetupgrads()

subroutine, public mod_realinput_grads::parentoceansetupgrads ( integer, dimension(2), intent(out)  odims,
integer, intent(out)  timelen,
real(rp), dimension(:,:), intent(out), allocatable  lon_all,
real(rp), dimension(:,:), intent(out), allocatable  lat_all,
character(len=*), intent(in)  basename,
character(len=*), intent(in)  basename_num 
)

Ocean Setup.

Definition at line 965 of file mod_realinput_grads.F90.

965  use scale_const, only: &
966  d2r => const_d2r
967  use scale_file_grads, only: &
968  file_grads_open, &
970  file_grads_get_shape, &
972  implicit none
973 
974  integer, intent(out) :: odims(2)
975  integer, intent(out) :: timelen
976  real(RP), allocatable, intent(out) :: lon_all(:,:)
977  real(RP), allocatable, intent(out) :: lat_all(:,:)
978 
979  character(len=*), intent(in) :: basename
980  character(len=*), intent(in) :: basename_num
981 
982  character(len=7) :: vname
983  integer :: vid
984  integer :: shape(2)
985  logical :: exist
986  !---------------------------------------------------------------------------
987 
988  log_info("ParentOceanSetupGrADS",*) 'Real Case/Ocean Input File Type: GrADS format'
989 
990  if ( basename == "" ) then
991  log_error("ParentOceanSetupGrADS",*) '"BASENAME_ORG" is not specified in "PARAM_MKINIT_REAL_OCEAN"!', trim(basename)
992  call prc_abort
993  endif
994  call file_grads_open( basename, & ! (in)
995  file_id_ocn ) ! (out)
996 
997  call file_grads_varid( file_id_ocn, "OCEAN_SFC_TEMP", & ! (in)
998  vid ) ! (out)
999  if ( vid < 0 ) then
1000  call file_grads_varid( file_id_ocn, "SST", & ! (in)
1001  vid ) ! (out)
1002  end if
1003  if ( vid < 0 ) then
1004  call file_grads_varid( file_id_ocn, "SFC_TEMP", & ! (in)
1005  vid ) ! (out)
1006  end if
1007  if ( vid < 0 ) then
1008  call file_grads_varid( file_id_ocn, "SKINT", & ! (in)
1009  vid ) ! (out)
1010  end if
1011  if ( vid < 0 ) then
1012  log_error("ParentOceanSetupGrADS",*) '"OCEAN_SFC_TEMP", "SST", "SFC_TEMP", or "SKINT" is necessary'
1013  call prc_abort
1014  end if
1015 
1016  call file_grads_get_shape( file_id_ocn, vid, & ! (in)
1017  odims(:) ) ! (out)
1018 
1019 
1020  ! get lon, lat
1021 
1022  allocate( lon_all(odims(1), odims(2)) )
1023  allocate( lat_all(odims(1), odims(2)) )
1024 
1025 
1026  call file_grads_varid( file_id_ocn, "lon_sst", & ! (in)
1027  vid ) ! (out)
1028  if ( vid > 0 ) then
1029  vname = "lon_sst"
1030  else
1031  call file_grads_varid( file_id_ocn, "lon_sfc", & ! (in)
1032  vid ) ! (out)
1033  if ( vid > 0 ) then
1034  call file_grads_get_shape( file_id_ocn, vid, & ! (in)
1035  shape(:) ) ! (out)
1036  if ( file_grads_isoned( file_id_lnd, vid ) ) then
1037  if ( odims(1) .eq. shape(1) .or. shape(1) .eq. -1 ) then
1038  vname = "lon_sfc"
1039  else
1040  vid = -1
1041  end if
1042  else
1043  if ( odims(1) .eq. shape(1) .and. odims(2) .eq. shape(2) ) then
1044  vname = "lon_sfc"
1045  else
1046  vid = -1
1047  end if
1048  end if
1049  end if
1050  end if
1051  if ( vid < 0 ) then
1052  call file_grads_varid( file_id_ocn, "lon", & ! (in)
1053  vid ) ! (out)
1054  if ( vid < 0 ) then
1055  log_error("ParentLandSetupGrADS",*) '"lon_sst", "lon_sfc", or "lon" is necessary'
1056  call prc_abort
1057  end if
1058  call file_grads_get_shape( file_id_ocn, vid, & ! (in)
1059  shape(:) ) ! (out)
1060  if ( file_grads_isoned( file_id_ocn, vid ) ) then
1061  if ( odims(1) .eq. shape(1) .or. shape(1) .eq. -1 ) then
1062  vname = "lon"
1063  else
1064  vid = -1
1065  end if
1066  else
1067  if ( odims(1) .eq. shape(1) .and. odims(2) .eq. shape(2) ) then
1068  vname = "lon"
1069  else
1070  vid = -1
1071  end if
1072  end if
1073  if ( vid < 0 ) then
1074  log_error("ParentOceanSetupGrADS",*) 'dimension of "lon_sfc" and "lon" is different! ', odims(:), shape(:)
1075  call prc_abort
1076  end if
1077  end if
1078  call read2d( (/1,1/), odims(:), lon_all(:,:), vname, file_id_ocn, basename_num, oned=1 )
1079  lon_all(:,:) = lon_all(:,:) * d2r
1080 
1081 
1082  call file_grads_varid( file_id_ocn, "lat_sst", & ! (in)
1083  vid ) ! (out)
1084  if ( vid > 0 ) then
1085  vname = "lat_sst"
1086  else
1087  call file_grads_varid( file_id_ocn, "lat_sfc", & ! (in)
1088  vid ) ! (out)
1089  if ( vid > 0 ) then
1090  call file_grads_get_shape( file_id_ocn, vid, & ! (in)
1091  shape(:) ) ! (out)
1092  if ( file_grads_isoned( file_id_lnd, vid ) ) then
1093  if ( odims(2) .eq. shape(1) .or. shape(1) .eq. -1 ) then
1094  vname = "lat_sfc"
1095  else
1096  vid = -1
1097  end if
1098  else
1099  if ( odims(1) .eq. shape(1) .and. odims(2) .eq. shape(2) ) then
1100  vname = "lat_sfc"
1101  else
1102  vid = -1
1103  end if
1104  end if
1105  end if
1106  end if
1107  if ( vid < 0 ) then
1108  call file_grads_varid( file_id_lnd, "lat", & ! (in)
1109  vid ) ! (out)
1110  if ( vid < 0 ) then
1111  log_error("ParentLandSetupGrADS",*) '"lat_sst", "lat_sfc", or "lat" is necessary'
1112  call prc_abort
1113  end if
1114  call file_grads_get_shape( file_id_ocn, vid, & ! (in)
1115  shape(:) ) ! (out)
1116  if ( file_grads_isoned( file_id_lnd, vid ) ) then
1117  if ( odims(2) .eq. shape(1) .or. shape(1) .eq. -1 ) then
1118  vname = "lat"
1119  else
1120  vid = -1
1121  end if
1122  else
1123  if ( odims(1) .eq. shape(1) .and. odims(2) .eq. shape(2) ) then
1124  vname = "lat"
1125  else
1126  vid = -1
1127  end if
1128  end if
1129  if ( vid < 0 ) then
1130  log_error("ParentOceanSetupGrADS",*) 'dimension of "lat_sfc" and "lat" is different! ', odims(:), shape(:)
1131  call prc_abort
1132  end if
1133  end if
1134  call read2d( (/1,1/), odims(:), lat_all(:,:), vname, file_id_ocn, basename_num, oned=2 )
1135  lat_all(:,:) = lat_all(:,:) * d2r
1136 
1137 
1138  ! tentative
1139  timelen = 1
1140 
1141  return

References scale_const::const_d2r, scale_file_grads::file_grads_isoned(), scale_file_grads::file_grads_open(), scale_file_grads::file_grads_varid(), scale_prc::prc_abort(), and read2d().

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentoceaninputgrads()

subroutine, public mod_realinput_grads::parentoceaninputgrads ( integer, intent(in)  IA_org,
integer, intent(in)  IS_org,
integer, intent(in)  IE_org,
integer, intent(in)  JA_org,
integer, intent(in)  JS_org,
integer, intent(in)  JE_org,
real(rp), dimension (ia_org,ja_org), intent(out)  tw_org,
real(rp), dimension (ia_org,ja_org), intent(out)  sst_org,
real(rp), dimension(ia_org,ja_org), intent(out)  omask_org,
character(len=*), intent(in)  basename_num,
integer, dimension(2), intent(in)  odims,
integer, intent(in)  nt 
)

Definition at line 1154 of file mod_realinput_grads.F90.

1154  use scale_const, only: &
1155  undef => const_undef, &
1156  d2r => const_d2r, &
1157  tem00 => const_tem00, &
1158  eps => const_eps
1159  use scale_file_grads, only: &
1160  file_grads_get_shape, &
1161  file_grads_varid, &
1162  file_grads_read
1163  implicit none
1164  integer, intent(in) :: IA_org, IS_org, IE_org
1165  integer, intent(in) :: JA_org, JS_org, JE_org
1166 
1167  real(RP), intent(out) :: tw_org (IA_org,JA_org)
1168  real(RP), intent(out) :: sst_org (IA_org,JA_org)
1169  real(RP), intent(out) :: omask_org(IA_org,JA_org)
1170 
1171  character(len=*), intent(in) :: basename_num
1172  integer, intent(in) :: odims(2)
1173  integer, intent(in) :: nt
1174 
1175  integer :: start(2), count(2), shape(2)
1176  integer :: vid
1177  logical :: exist
1178  integer :: i, j
1179  !---------------------------------------------------------------------------
1180 
1181 
1182  start(:) = (/is_org,js_org/)
1183  count(:) = (/ia_org,ja_org/)
1184 
1185  ! lsmask
1186  call read2d( start(:), count(:), omask_org(:,:), "lsmask_sst", file_id_ocn, basename_num, exist=exist, step=nt )
1187  if ( .not. exist ) then
1188  call file_grads_varid( file_id_ocn, "lsmask", & ! (in)
1189  vid ) ! (out)
1190  if ( vid > 0 ) then
1191  call file_grads_get_shape( file_id_ocn, vid, & ! (in)
1192  shape(:) ) ! (out)
1193  if ( odims(1) .ne. shape(1) .or. odims(2) .ne. shape(2) ) then
1194  log_warn("ParentOceanInputGrADS",*) 'dimension of lsmask is different. not use'
1195  else
1196  call read2d( start(:), count(:), omask_org(:,:), "lsmask", file_id_ocn, basename_num, exist=exist, step=nt )
1197  end if
1198  end if
1199  end if
1200  if ( .not. exist ) then
1201  !$omp parallel do
1202  do j = 1, ja_org
1203  do i = 1, ia_org
1204  omask_org(i,j) = undef
1205  end do
1206  end do
1207  end if
1208 
1209  ! OCEAN_SFC_TEMP, SST, SFC_TEMP, SKINT
1210  call read2d( start(:), count(:), sst_org(:,:), "OCEAN_SFC_TEMP", file_id_ocn, basename_num, exist=exist, step=nt )
1211  if ( .not. exist ) then
1212  call read2d( start(:), count(:), sst_org(:,:), "SST", file_id_ocn, basename_num, exist=exist, step=nt )
1213  end if
1214  if ( .not. exist ) then
1215  call file_grads_varid( file_id_ocn, "SFC_TEMP", & ! (in)
1216  vid ) ! (out)
1217  if ( vid > 0 ) then
1218  call file_grads_get_shape( file_id_ocn, vid, & ! (in)
1219  shape(:) ) ! (out)
1220  if ( odims(1).eq.shape(1) .and. odims(2).eq.shape(2) ) then
1221  call read2d( start(:), count(:), sst_org(:,:), "SFC_TEMP", file_id_ocn, basename_num, exist=exist, step=nt )
1222  else
1223  exist = .false.
1224  end if
1225  else
1226  exist = .false.
1227  end if
1228  end if
1229  if ( .not. exist ) then
1230  call file_grads_varid( file_id_ocn, "SKINT", & ! (in)
1231  vid ) ! (out)
1232  call file_grads_get_shape( file_id_ocn, vid, & ! (in)
1233  shape(:) ) ! (out)
1234  if ( odims(1).eq.shape(1) .and. odims(2).eq.shape(2) ) then
1235  call read2d( start(:), count(:), sst_org(:,:), "SKINT", file_id_ocn, basename_num, step=nt )
1236  else
1237  log_error("ParentOceanInputGrADS",*) 'dimension of "SFC_TEMP" and/or "SKINT" is different'
1238  call prc_abort
1239  end if
1240  end if
1241 
1242  tw_org = sst_org
1243 
1244  return

References scale_const::const_d2r, scale_const::const_eps, scale_const::const_tem00, scale_const::const_undef, scale_file_grads::file_grads_varid(), scale_prc::prc_abort(), and read2d().

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ read2d()

subroutine mod_realinput_grads::read2d ( integer, dimension(2), intent(in)  start,
integer, dimension(2), intent(in)  count,
real(rp), dimension(count(1),count(2)), intent(out)  data,
character(len=*), intent(in)  name,
integer, intent(in)  fid,
character(len=*), intent(in)  postfix,
logical, intent(out), optional  exist,
integer, intent(in), optional  oneD,
integer, intent(in), optional  step 
)

Definition at line 1250 of file mod_realinput_grads.F90.

1250  use scale_file_grads, only: &
1251  file_grads_varid, &
1252  file_grads_read, &
1254  integer, intent(in) :: start(2)
1255  integer, intent(in) :: count(2)
1256  real(RP), intent(out) :: data(count(1),count(2))
1257  character(len=*), intent(in) :: name
1258  integer, intent(in) :: fid
1259  character(len=*), intent(in) :: postfix
1260 
1261  logical, intent(out), optional :: exist
1262  integer, intent(in), optional :: oneD
1263  integer, intent(in), optional :: step
1264 
1265  integer :: vid
1266  real(RP), allocatable :: v1d(:)
1267  integer :: oneD_
1268  integer :: i, j
1269 
1270  call file_grads_varid( fid, name, & ! (in)
1271  vid ) ! (out)
1272  if ( vid < 0 ) then
1273  if ( present(exist) ) then
1274  exist = .false.
1275  return
1276  end if
1277  log_error("read2d",*) '"', trim(name), '" is required'
1278  call prc_abort
1279  end if
1280  if ( present(exist) ) exist = .true.
1281 
1282  if( file_grads_isoned( fid, vid ) ) then
1283  if ( present(oned) ) then
1284  oned_ = oned
1285  else
1286  oned_ = 1
1287  end if
1288  allocate( v1d(count(oned_)) )
1289  call file_grads_read( fid, vid, & ! (in)
1290  v1d(:), & ! (out)
1291  step=step, & ! (in)
1292  start=start(oned_:oned_), & ! (in)
1293  count=count(oned_:oned_) ) ! (in)
1294  if ( oned_ == 1 ) then
1295  !$omp parallel do
1296  do j = 1, count(2)
1297  do i = 1, count(1)
1298  data(i,j) = v1d(i)
1299  enddo
1300  enddo
1301  else
1302  !$omp parallel do
1303  do j = 1, count(2)
1304  do i = 1, count(1)
1305  data(i,j) = v1d(j)
1306  enddo
1307  enddo
1308  end if
1309  deallocate( v1d )
1310  else
1311  call file_grads_read( fid, vid, & ! (in)
1312  data(:,:), & ! (out)
1313  step=step, & ! (in)
1314  start=start(:), & ! (in)
1315  count=count(:), & ! (in)
1316  postfix=postfix ) ! (in)
1317  end if
1318 
1319  return

References scale_file_grads::file_grads_isoned(), scale_file_grads::file_grads_varid(), and scale_prc::prc_abort().

Referenced by parentatmosinputgrads(), parentatmossetupgrads(), parentlandinputgrads(), parentlandsetupgrads(), parentoceaninputgrads(), and parentoceansetupgrads().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ read3d()

subroutine mod_realinput_grads::read3d ( integer, dimension(3), intent(in)  start,
integer, dimension(3), intent(in)  count,
real(rp), dimension(count(1),count(2),count(3)), intent(out)  data,
character(len=*), intent(in)  name,
integer, intent(in)  fid,
character(len=*), intent(in)  postfix,
logical, intent(out), optional  exist,
integer, intent(in), optional  step 
)

Definition at line 1323 of file mod_realinput_grads.F90.

1323  use scale_file_grads, only: &
1324  file_grads_varid, &
1325  file_grads_read, &
1327  integer, intent(in) :: start(3)
1328  integer, intent(in) :: count(3)
1329  real(RP), intent(out) :: data(count(1),count(2),count(3))
1330  character(len=*), intent(in) :: name
1331  integer, intent(in) :: fid
1332  character(len=*), intent(in) :: postfix
1333 
1334  logical, intent(out), optional :: exist
1335  integer, intent(in), optional :: step
1336 
1337  integer :: vid
1338  real(RP), allocatable :: v1d(:)
1339  integer :: k, i, j
1340 
1341  call file_grads_varid( fid, name, & ! (in)
1342  vid ) ! (out)
1343  if ( vid < 0 ) then
1344  if ( present(exist) ) then
1345  exist = .false.
1346  return
1347  end if
1348  log_error("read3d",*) '"', trim(name), '" is required'
1349  call prc_abort
1350  end if
1351  if ( present(exist) ) exist = .true.
1352 
1353  if( file_grads_isoned( fid, vid ) ) then
1354  allocate( v1d(count(1)) )
1355  call file_grads_read( fid, vid, & ! (in)
1356  v1d(:), & ! (out)
1357  step=step, &
1358  start=start(1:1), & ! (in)
1359  count=count(1:1) ) ! (in)
1360  !$omp parallel do collapse(2)
1361  do j = 1, count(3)
1362  do i = 1, count(2)
1363  do k = 1, count(1)
1364  data(k,i,j) = v1d(k)
1365  enddo
1366  enddo
1367  enddo
1368  deallocate( v1d )
1369  else
1370  call file_grads_read( fid, vid, & ! (in)
1371  data(:,:,:), & ! (out)
1372  step=step, & ! (in)
1373  start=start(:), & ! (in)
1374  count=count(:), & ! (in)
1375  postfix=postfix ) ! (in)
1376  end if
1377 
1378  return

References scale_file_grads::file_grads_isoned(), scale_file_grads::file_grads_varid(), and scale_prc::prc_abort().

Referenced by parentatmosinputgrads(), and parentlandinputgrads().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
mod_atmos_phy_mp_vars
module Atmosphere / Physics Cloud Microphysics
Definition: mod_atmos_phy_mp_vars.F90:12
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_file_grads::file_grads_isoned
logical function, public file_grads_isoned(file_id, var_id)
Definition: scale_file_grads.F90:373
scale_const::const_epstvap
real(rp), public const_epstvap
1 / epsilon - 1
Definition: scale_const.F90:76
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
scale_file_grads::file_grads_open
subroutine, public file_grads_open(file_name, file_id)
Open.
Definition: scale_file_grads.F90:104
scale_atmos_hydrometeor::hyd_name
character(len=h_short), dimension(n_hyd), parameter, public hyd_name
Definition: scale_atmos_hydrometeor.F90:104
scale_const::const_epsvap
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:75
mod_atmos_phy_mp_vars::qs_mp
integer, public qs_mp
Definition: mod_atmos_phy_mp_vars.F90:79
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_file_grads::file_grads_varcheck
subroutine, public file_grads_varcheck(file_id, var_name, exist)
Definition: scale_file_grads.F90:344
scale_file_grads::file_grads_varid
subroutine, public file_grads_varid(file_id, var_name, var_id)
Definition: scale_file_grads.F90:312
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
scale_file_grads
module file_grads
Definition: scale_file_grads.F90:11
scale_const::const_tem00
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:99
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:59
mod_atmos_phy_mp_vars::qe_mp
integer, public qe_mp
Definition: mod_atmos_phy_mp_vars.F90:80
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
scale_const::const_laps
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:62
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:95