SCALE-RM
Functions/Subroutines
mod_realinput_grads Module Reference

module REAL input GrADS More...

Functions/Subroutines

subroutine, public parentatomsetupgrads (dims, basename)
 Atmos Setup. More...
 
subroutine, public parentatomopengrads
 
subroutine, public parentatominputgrads (velz_org, velx_org, vely_org, pres_org, temp_org, qtrc_org, lon_org, lat_org, cz_org, basename_num, dims, nt)
 
subroutine, public parentlandsetupgrads (ldims, use_waterratio, use_file_landwater, basename)
 Land Setup. More...
 
subroutine, public parentlandinputgrads (tg_org, strg_org, smds_org, lst_org, llon_org, llat_org, lz_org, topo_org, lmask_org, basename_num, ldims, use_file_landwater, nt)
 
subroutine, public parentoceansetupgrads (odims, timelen, basename)
 Ocean Setup. More...
 
subroutine, public parentoceanopengrads
 
subroutine, public parentoceaninputgrads (tw_org, sst_org, omask_org, olon_org, olat_org, basename_num, odims, nt)
 
subroutine read_namelist (grads_item, grads_fname, grads_dtype, grads_swpoint, grads_dd, grads_lnum, grads_lvars, grads_startrec, grads_totalrec, grads_knum, grads_yrev, grads_fendian, grads_missval, data_available, item_list, num_item_list, basename, io_fid_grads_nml)
 

Detailed Description

module REAL input GrADS

Description
read data from GrADS file for real atmospheric simulations
Author
Team SCALE
History
NAMELIST
  • nml_grads_grid
    nametypedefault valuecomment
    OUTER_NX integer -1
    OUTER_NY integer -1
    OUTER_NZ integer -1 number of atmos layers
    OUTER_NL integer -1 number of land layers
    OUTER_NX_SFC integer -1
    OUTER_NY_SFC integer -1
    OUTER_NX_SST integer -1
    OUTER_NY_SST integer -1

  • PARAM_MKINIT_REAL_GrADS
    nametypedefault valuecomment
    UPPER_QV_TYPE character(len=H_SHORT) "ZERO" how qv is given at higher level than outer model

History Output
No history output

Function/Subroutine Documentation

◆ parentatomsetupgrads()

subroutine, public mod_realinput_grads::parentatomsetupgrads ( integer, dimension(6), intent(out)  dims,
character(len=h_long), intent(in)  basename 
)

Atmos Setup.

Definition at line 183 of file mod_realinput_grads.f90.

References scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_stdio::io_lnml, scale_process::prc_mpistop(), and read_namelist().

Referenced by mod_realinput::parentatomsetup().

183  implicit none
184 
185  integer, intent(out) :: dims(6)
186  character(len=H_LONG), intent(in) :: basename
187 
188 
189  namelist / param_mkinit_real_grads / &
190  upper_qv_type
191 
192  integer :: ielem
193  integer :: k, n
194 
195  integer :: ierr
196  !---------------------------------------------------------------------------
197 
198  if( io_l ) write(io_fid_log,*) '+++ Real Case/Atom Input File Type: GrADS format'
199 
200  !--- read namelist
201  rewind(io_fid_conf)
202  read(io_fid_conf,nml=param_mkinit_real_grads,iostat=ierr)
203 
204  if( ierr > 0 ) then
205  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_GrADS. Check!'
206  call prc_mpistop
207  endif
208  if( io_lnml ) write(io_fid_log,nml=param_mkinit_real_grads)
209 
210 
211  if ( len_trim(basename) == 0 ) then
212  if( io_l ) write(io_fid_log,*) &
213  'xxx "BASENAME_ORG" is not specified in "PARAM_MKINIT_REAL_ATMOS"!',trim(basename)
214  call prc_mpistop
215  endif
216 
217  !--- read namelist
218  io_fid_grads_nml = io_get_available_fid()
219  open( io_fid_grads_nml, &
220  file = trim(basename), &
221  form = 'formatted', &
222  status = 'old', &
223  action = 'read', &
224  iostat = ierr )
225  if ( ierr /= 0 ) then
226  if( io_l ) write(io_fid_log,*) 'xxx Input file is not found! ', trim(basename)
227  call prc_mpistop
228  endif
229 
230  read(io_fid_grads_nml,nml=nml_grads_grid,iostat=ierr)
231  if( ierr /= 0 ) then !--- missing or fatal error
232  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in nml_grads_grid in ', trim(basename),'. Check!'
233  call prc_mpistop
234  endif
235  if( io_lnml ) write(io_fid_log,nml=nml_grads_grid)
236 
237  ! full level
238  dims(1) = outer_nz ! bottom_top
239  dims(2) = outer_nx ! west_east
240  dims(3) = outer_ny ! south_north
241  ! half level
242  dims(4) = outer_nz ! bottom_top_stag
243  dims(5) = outer_nx ! west_east for 2dim data
244  dims(6) = outer_ny ! south_north for 2dim data
245 
246  allocate( gdata2d( dims(2), dims(3) ) )
247  allocate( gdata3d( dims(2), dims(3), dims(1) ) )
248 
249  call read_namelist( &
250  grads_item(:,1), & ! (out)
251  grads_fname(:,1), & ! (out)
252  grads_dtype(:,1), & ! (out)
253  grads_swpoint(:,1), & ! (out)
254  grads_dd(:,1), & ! (out)
255  grads_lnum(:,1), & ! (out)
256  grads_lvars(:,:,1), & ! (out)
257  grads_startrec(:,1), & ! (out)
258  grads_totalrec(:,1), & ! (out)
259  grads_knum(:,1), & ! (out)
260  grads_yrev(:,1), & ! (out)
261  grads_fendian(:,1), & ! (out)
262  grads_missval(:,1), & ! (out)
263  data_available(:,1), & ! (out)
264  item_list_atom, & ! (in)
265  num_item_list_atom, & ! (in)
266  basename, & ! (in)
267  io_fid_grads_nml ) ! (in)
268 
269  do ielem = 1, num_item_list_atom
270  item = item_list_atom(ielem)
271  !--- check data
272  select case(trim(item))
273  case('QV')
274  if (.not. data_available(ia_qv,1)) then
275  if (.not.data_available(ia_rh,1)) then
276  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : QV and RH'
277  call prc_mpistop
278  else ! will read RH
279  cycle
280  endif
281  endif
282  case('RH')
283  if (.not. data_available(ia_qv,1))then
284  if(data_available(ia_rh,1)) then
285  if ((.not. data_available(ia_t,1)).or.(.not. data_available(ia_p,1))) then
286  if( io_l ) write(io_fid_log,*) 'xxx Temperature and pressure are required to convert from RH to QV ! '
287  call prc_mpistop
288  else
289  cycle ! read RH and estimate QV
290  endif
291  else
292  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : QV and RH'
293  call prc_mpistop
294  endif
295  endif
296  case('MSLP','PSFC','U10','V10','T2')
297  if (.not. data_available(ielem,1)) then
298  if (io_l) write(io_fid_log,*) 'warning: ',trim(item),' is not found & will be estimated.'
299  cycle
300  endif
301  case('Q2')
302  if ( .not. data_available(ia_q2,1) ) then
303  if (io_l) write(io_fid_log,*) 'warning: Q2 is not found & will be estimated.'
304  cycle
305  endif
306  case('RH2')
307  if ( data_available(ia_q2,1) ) then
308  cycle
309  else
310  if ( data_available(ia_rh2,1) ) then
311  if ((.not. data_available(ia_t2,1)).or.(.not. data_available(ia_ps,1))) then
312  if (io_l) write(io_fid_log,*) 'warning: T2 and PSFC are required to convert from RH2 to Q2 !'
313  if (io_l) write(io_fid_log,*) ' Q2 will be copied from data at above level.'
314  data_available(ia_rh2,1) = .false.
315  cycle
316  endif
317  else
318  if (io_l) write(io_fid_log,*) 'warning: Q2 and RH2 are not found, Q2 will be estimated.'
319  cycle
320  endif
321  endif
322  case('TOPO')
323  if ( .not. data_available(ielem,1) ) then
324  if (io_l) write(io_fid_log,*) 'warning: ',trim(item),' is not found & not used.'
325  cycle
326  endif
327  case default ! lon, lat, plev, U, V, T, HGT
328  if ( .not. data_available(ielem,1) ) then
329  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : ',trim(item_list_atom(ielem))
330  call prc_mpistop
331  endif
332  end select
333 
334  end do
335 
336  return
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer function, public io_get_available_fid()
search & get available file ID
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentatomopengrads()

subroutine, public mod_realinput_grads::parentatomopengrads ( )

Definition at line 341 of file mod_realinput_grads.f90.

References scale_stdio::io_fid_log, and scale_stdio::io_l.

Referenced by mod_realinput::parentatomsetup().

341  implicit none
342 
343  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[AtomOpenGrADS]'
344 
345  return
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the caller graph for this function:

◆ parentatominputgrads()

subroutine, public mod_realinput_grads::parentatominputgrads ( real(rp), dimension(:,:,:), intent(out)  velz_org,
real(rp), dimension(:,:,:), intent(out)  velx_org,
real(rp), dimension(:,:,:), intent(out)  vely_org,
real(rp), dimension(:,:,:), intent(out)  pres_org,
real(rp), dimension(:,:,:), intent(out)  temp_org,
real(rp), dimension(:,:,:,:), intent(out)  qtrc_org,
real(rp), dimension(:,:), intent(out)  lon_org,
real(rp), dimension(:,:), intent(out)  lat_org,
real(rp), dimension(:,:,:), intent(out)  cz_org,
character(len=*), intent(in)  basename_num,
integer, dimension(6), intent(in)  dims,
integer, intent(in)  nt 
)

Definition at line 362 of file mod_realinput_grads.f90.

References scale_const::const_cpdry, scale_const::const_d2r, scale_const::const_eps, scale_const::const_epsvap, scale_const::const_grav, scale_const::const_laps, scale_const::const_pre00, scale_const::const_rdry, scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, dc_log::log(), and scale_process::prc_mpistop().

Referenced by mod_realinput::parentatomsetup().

362  use scale_const, only: &
363  d2r => const_d2r, &
364  eps => const_eps, &
365  epsvap => const_epsvap, &
366  grav => const_grav, &
367  laps => const_laps, &
368  p00 => const_pre00, &
369  rdry => const_rdry, &
370  cpdry => const_cpdry
371  use scale_atmos_saturation, only: &
372  psat => atmos_saturation_psat_liq
373  implicit none
374 
375 
376  real(RP), intent(out) :: velz_org(:,:,:)
377  real(RP), intent(out) :: velx_org(:,:,:)
378  real(RP), intent(out) :: vely_org(:,:,:)
379  real(RP), intent(out) :: pres_org(:,:,:)
380  real(RP), intent(out) :: temp_org(:,:,:)
381  real(RP), intent(out) :: qtrc_org(:,:,:,:)
382  real(RP), intent(out) :: lon_org(:,:)
383  real(RP), intent(out) :: lat_org(:,:)
384  real(RP), intent(out) :: cz_org(:,:,:)
385  character(len=*), intent(in) :: basename_num
386  integer, intent(in) :: dims(6)
387  integer, intent(in) :: nt
388 
389  real(RP) :: rhprs_org(dims(1)+2,dims(2),dims(3))
390  real(RP) :: pott(dims(2),dims(3))
391 
392  real(RP) :: rovcp
393  real(RP) :: cpovr
394 
395  ! data
396  character(len=H_LONG) :: gfile
397 
398  integer :: qa_outer = 1
399  real(RP) :: p_sat, qm, rhsfc
400 
401  integer :: i, j, k, ielem
402 
403  !---------------------------------------------------------------------------
404 
405  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[AtomInputGrADS]'
406 
407  qtrc_org = 0.0_rp
408 
409  !--- read grads data
410  loop_inputatomgrads : do ielem = 1, num_item_list_atom
411 
412  if ( .not. data_available(ielem,1) ) cycle
413 
414  item = grads_item(ielem,1)
415  dtype = grads_dtype(ielem,1)
416  fname = grads_fname(ielem,1)
417  lnum = grads_lnum(ielem,1)
418 
419  if ( dims(1) < grads_knum(ielem,1) ) then
420  write(*,*) 'xxx "knum" must be less than or equal to outer_nz. knum:',knum,'> outer_nz:',dims(1),trim(item)
421  call prc_mpistop
422  else if ( grads_knum(ielem,1) > 0 )then
423  knum = grads_knum(ielem,1) ! not missing
424  else
425  knum = dims(1)
426  endif
427 
428  select case (trim(dtype))
429  case("linear")
430  swpoint = grads_swpoint(ielem,1)
431  dd = grads_dd(ielem,1)
432  if( (abs(swpoint-large_number_one)<eps).or.(abs(dd-large_number_one)<eps) )then
433  write(*,*) 'xxx "swpoint" is required in grads namelist! ',swpoint
434  write(*,*) 'xxx "dd" is required in grads namelist! ',dd
435  call prc_mpistop
436  endif
437  case("levels")
438  if ( lnum < 0 )then
439  write(*,*) 'xxx "lnum" is required in grads namelist for levels data! '
440  call prc_mpistop
441  endif
442  do k=1, lnum
443  lvars(k)=grads_lvars(k,ielem,1)
444  enddo
445  if(abs(lvars(1)-large_number_one)<eps)then
446  write(*,*) 'xxx "lvars" must be specified in grads namelist for levels data! '
447  call prc_mpistop
448  endif
449  case("map")
450  startrec = grads_startrec(ielem,1)
451  totalrec = grads_totalrec(ielem,1)
452  fendian = grads_fendian(ielem,1)
453  yrev = grads_yrev(ielem,1)
454  if( (startrec<0).or.(totalrec<0) )then
455  write(*,*) 'xxx "startrec" is required in grads namelist! ',startrec
456  write(*,*) 'xxx "totalrec" is required in grads namelist! ',totalrec
457  call prc_mpistop
458  endif
459  ! get file_id
460  if(io_fid_grads_data < 0)then
461  io_fid_grads_data = io_get_available_fid()
462  endif
463  gfile=trim(fname)//trim(basename_num)//'.grd'
464  if( len_trim(fname)==0 )then
465  write(*,*) 'xxx "fname" is required in grads namelist for map data! ',trim(fname)
466  call prc_mpistop
467  endif
468  end select
469 
470  ! read data
471  select case (trim(item))
472  case("lon")
473  if ( trim(dtype) == "linear" ) then
474  do j = 1, dims(3)
475  do i = 1, dims(2)
476  lon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * d2r
477  enddo
478  enddo
479  else if ( trim(dtype) == "map" ) then
480  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,1,item,startrec,totalrec,yrev,gdata2d)
481  lon_org(:,:) = real(gdata2D(:,:), kind=RP) * d2r
482  endif
483  case("lat")
484  if ( trim(dtype) == "linear" ) then
485  do j = 1, dims(3)
486  do i = 1, dims(2)
487  lat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * d2r
488  enddo
489  enddo
490  else if ( trim(dtype) == "map" ) then
491  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,1,item,startrec,totalrec,yrev,gdata2d)
492  lat_org(:,:) = real(gdata2D(:,:), kind=RP) * d2r
493  endif
494  case("plev")
495  if(dims(1)/=knum)then
496  write(*,*) 'xxx "knum" must be equal to outer_nz for plev. knum:',knum,'> outer_nz:',dims(1)
497  call prc_mpistop
498  endif
499  if ( trim(dtype) == "levels" ) then
500  if(dims(1)/=lnum)then
501  write(*,*) 'xxx lnum must be same as the outer_nz for plev! ',dims(1),lnum
502  call prc_mpistop
503  endif
504  do j = 1, dims(3)
505  do i = 1, dims(2)
506  do k = 1, dims(1)
507  pres_org(k+2,i,j) = real(lvars(k), kind=rp)
508  enddo
509  enddo
510  enddo
511  else if ( trim(dtype) == "map" ) then
512  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),dims(1),nt,item,startrec,totalrec,yrev,gdata3d)
513  do j = 1, dims(3)
514  do i = 1, dims(2)
515  do k = 1, dims(1)
516  pres_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
517  enddo
518  enddo
519  enddo
520  endif
521  case('U')
522  if ( trim(dtype) == "map" ) then
523  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
524  do j = 1, dims(3)
525  do i = 1, dims(2)
526  velx_org(1:2,i,j) = 0.0_rp
527  do k = 1, knum
528  velx_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
529  enddo
530  if(dims(1)>knum)then
531  do k = knum+1, dims(1)
532  velx_org(k+2,i,j) = velx_org(knum+2,i,j)
533  enddo
534  endif
535  enddo
536  enddo
537  endif
538  case('V')
539  if ( trim(dtype) == "map" ) then
540  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
541  do j = 1, dims(3)
542  do i = 1, dims(2)
543  vely_org(1:2,i,j) = 0.0_rp
544  do k = 1, knum
545  vely_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
546  enddo
547  if(dims(1)>knum)then
548  do k = knum+1, dims(1)
549  vely_org(k+2,i,j) = vely_org(knum+2,i,j)
550  enddo
551  endif
552  enddo
553  enddo
554  endif
555  case('T')
556  if ( trim(dtype) == "map" ) then
557  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
558  do j = 1, dims(3)
559  do i = 1, dims(2)
560  do k = 1, knum
561  temp_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
562  enddo
563  if(dims(1)>knum)then
564  do k = knum+1, dims(1)
565  temp_org(k+2,i,j) = temp_org(knum+2,i,j)
566  enddo
567  endif
568  enddo
569  enddo
570  endif
571  case('HGT')
572  if(dims(1)/=knum)then
573  write(*,*) 'xxx The number of levels for HGT must be same as plevs! knum:',knum,'> outer_nz:',dims(1)
574  call prc_mpistop
575  endif
576  if ( trim(dtype) == "map" ) then
577  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),dims(1),nt,item,startrec,totalrec,yrev,gdata3d)
578  do j = 1, dims(3)
579  do i = 1, dims(2)
580  do k = 1, dims(1)
581  cz_org(k+2,i,j) = real(gdata3D(i,j,k), kind=rp)
582  enddo
583  cz_org(1,i,j) = 0.0_rp
584  enddo
585  enddo
586  endif
587  case('QV')
588  if ( trim(dtype) == "map" ) then
589  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
590  do j = 1, dims(3)
591  do i = 1, dims(2)
592  do k = 1, knum
593  qtrc_org(k+2,i,j,qa_outer) = real(gdata3D(i,j,k), kind=rp)
594  enddo
595  qtrc_org(1:2,i,j,qa_outer) = qtrc_org(3,i,j,qa_outer)
596  enddo
597  enddo
598  if( dims(1)>knum ) then
599  select case ( upper_qv_type )
600  case ("COPY")
601  do j = 1, dims(3)
602  do i = 1, dims(2)
603  do k = knum+1, dims(1)
604  qtrc_org(k+2,i,j,qa_outer) = qtrc_org(knum+2,i,j,qa_outer)
605  enddo
606  enddo
607  enddo
608  case ("ZERO")
609  ! do nothing
610  case default
611  write(*,*) 'xxx upper_qv_type in PARAM_MKINIT_REAL_GrADS is invalid! ', upper_qv_type
612  call prc_mpistop
613  end select
614  endif
615  endif
616  case('RH')
617  if (data_available(ia_qv,1)) cycle ! use QV
618  if ( trim(dtype) == "map" ) then
619  call read_grads_file_3d(io_fid_grads_data,gfile,dims(2),dims(3),knum,nt,item,startrec,totalrec,yrev,gdata3d)
620  do j = 1, dims(3)
621  do i = 1, dims(2)
622  do k = 1, knum
623  rhprs_org(k+2,i,j) = real(gdata3D(i,j,k), kind=RP) / 100.0_rp ! relative humidity
624  call psat( p_sat, temp_org(k+2,i,j) ) ! satulation pressure
625  qm = epsvap * rhprs_org(k+2,i,j) * p_sat &
626  / ( pres_org(k+2,i,j) - rhprs_org(k+2,i,j) * p_sat ) ! mixing ratio
627  qtrc_org(k+2,i,j,qa_outer) = qm / ( 1.0_rp + qm ) ! specific humidity
628  enddo
629  qtrc_org(1:2,i,j,qa_outer) = qtrc_org(3,i,j,qa_outer)
630  enddo
631  enddo
632  if( dims(3)>knum ) then
633  select case ( upper_qv_type )
634  case ("COPY")
635  do j = 1, dims(3)
636  do i = 1, dims(2)
637  do k = knum+1, dims(1)
638  rhprs_org(k+2,i,j) = rhprs_org(knum+2,i,j) ! relative humidity
639  call psat( p_sat, temp_org(k+2,i,j) ) ! satulated specific humidity
640  qm = epsvap * rhprs_org(k+2,i,j) * p_sat &
641  / ( pres_org(k+2,i,j) - rhprs_org(k+2,i,j) * p_sat ) ! mixing ratio
642  qtrc_org(k+2,i,j,qa_outer) = qm / ( 1.0_rp + qm ) ! specific humidity
643  qtrc_org(k+2,i,j,qa_outer) = min(qtrc_org(k+2,i,j,qa_outer),qtrc_org(k+1,i,j,qa_outer))
644  enddo
645  enddo
646  enddo
647  case ("ZERO")
648  ! do nothing
649  case default
650  write(*,*) 'xxx upper_qv_type in PARAM_MKINIT_REAL_GrADS is invalid! ', upper_qv_type
651  call prc_mpistop
652  end select
653  endif
654  endif
655  case('MSLP')
656  if ( trim(dtype) == "map" ) then
657  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
658  do j = 1, dims(3)
659  do i = 1, dims(2)
660  pres_org(1,i,j) = real(gdata2D(i,j), kind=rp)
661  enddo
662  enddo
663  endif
664  case('PSFC')
665  if ( trim(dtype) == "map" ) then
666  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
667  do j = 1, dims(3)
668  do i = 1, dims(2)
669  pres_org(2,i,j) = real(gdata2D(i,j), kind=rp)
670  enddo
671  enddo
672  endif
673  case('U10')
674  if ( trim(dtype) == "map" ) then
675  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
676  do j = 1, dims(3)
677  do i = 1, dims(2)
678  velx_org(2,i,j) = real(gdata2D(i,j), kind=rp)
679  enddo
680  enddo
681  endif
682  case('V10')
683  if ( trim(dtype) == "map" ) then
684  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
685  do j = 1, dims(3)
686  do i = 1, dims(2)
687  vely_org(2,i,j) = real(gdata2D(i,j), kind=rp)
688  enddo
689  enddo
690  endif
691  case('T2')
692  if ( trim(dtype) == "map" ) then
693  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
694  do j = 1, dims(3)
695  do i = 1, dims(2)
696  temp_org(2,i,j) = real(gdata2D(i,j), kind=rp)
697  enddo
698  enddo
699  endif
700  case('Q2')
701  if ( trim(dtype) == "map" ) then
702  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
703  do j = 1, dims(3)
704  do i = 1, dims(2)
705  qtrc_org(2,i,j,qa_outer) = real(gdata2D(i,j), kind=rp)
706  enddo
707  enddo
708  endif
709  case('RH2')
710  if (data_available(ia_q2,1)) cycle ! use QV
711  if ( trim(dtype) == "map" ) then
712  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
713  do j = 1, dims(3)
714  do i = 1, dims(2)
715  rhsfc = real(gdata2D(i,j), kind=RP) / 100.0_rp ! relative humidity
716  call psat( p_sat, temp_org(2,i,j) ) ! satulation pressure
717  qm = epsvap * rhsfc * p_sat &
718  / ( pres_org(2,i,j) - rhsfc * p_sat ) ! mixing ratio
719  qtrc_org(2,i,j,qa_outer) = qm / ( 1.0_rp + qm ) ! specific humidity
720  enddo
721  enddo
722  endif
723  case('TOPO')
724  if ( trim(dtype) == "map" ) then
725  call read_grads_file_2d(io_fid_grads_data,gfile,dims(2),dims(3),1,nt,item,startrec,totalrec,yrev,gdata2d)
726  do j = 1, dims(3)
727  do i = 1, dims(2)
728  cz_org(2,i,j) = real(gdata2D(i,j), kind=rp)
729  enddo
730  enddo
731  endif
732  end select
733  enddo loop_inputatomgrads
734 
735 
736  rovcp = rdry / cpdry
737  cpovr = cpdry / rdry
738 
739  if ( data_available(ia_t2,1) .and. data_available(ia_ps,1) ) then
740  do j = 1, dims(3)
741  do i = 1, dims(2)
742  pott(i,j) = temp_org(2,i,j) * (p00/pres_org(2,i,j))**rovcp
743  end do
744  end do
745  else
746  do j = 1, dims(3)
747  do i = 1, dims(2)
748  pott(i,j) = temp_org(3,i,j) * (p00/pres_org(3,i,j))**rovcp
749  end do
750  end do
751  end if
752 
753  if ( .not. data_available(ia_t2,1) ) then
754  if ( data_available(ia_ps,1) ) then
755  do j = 1, dims(3)
756  do i = 1, dims(2)
757  temp_org(2,i,j) = pott(i,j) * (pres_org(2,i,j)/p00)**rovcp
758  end do
759  end do
760  else
761  if ( data_available(ia_topo,1) ) then
762  do j = 1, dims(3)
763  do i = 1, dims(2)
764  temp_org(2,i,j) = temp_org(3,i,j) &
765  + laps * (cz_org(3,i,j)-cz_org(2,i,j))
766  end do
767  end do
768  else
769  do j = 1, dims(3)
770  do i = 1, dims(2)
771  temp_org(2,i,j) = temp_org(3,i,j)
772  end do
773  end do
774  end if
775  end if
776  end if
777 
778  if ( .not. data_available(ia_ps,1) ) then
779  do j = 1, dims(3)
780  do i = 1, dims(2)
781  pres_org(2,i,j) = p00 * (temp_org(2,i,j)/pott(i,j))**cpovr
782  end do
783  end do
784  end if
785 
786  if ( data_available(ia_slp,1) ) then
787  do j = 1, dims(3)
788  do i = 1, dims(2)
789  temp_org(1,i,j) = pott(i,j) * (pres_org(1,i,j)/p00)**rovcp
790  end do
791  end do
792  else
793  if ( data_available(ia_t2,1) .and. data_available(ia_topo,1) ) then
794  do j = 1, dims(3)
795  do i = 1, dims(2)
796  temp_org(1,i,j) = temp_org(2,i,j) + laps * cz_org(2,i,j)
797  end do
798  end do
799  else
800  do j = 1, dims(3)
801  do i = 1, dims(2)
802  temp_org(1,i,j) = temp_org(3,i,j) + laps * cz_org(3,i,j)
803  end do
804  end do
805  end if
806  do j = 1, dims(3)
807  do i = 1, dims(2)
808  pres_org(1,i,j) = p00 * (temp_org(1,i,j)/pott(i,j))**cpovr
809  end do
810  end do
811  end if
812 
813  if ( .not. data_available(ia_topo,1) ) then
814  ! guess surface height (elevation)
815  do j = 1, dims(3)
816  do i = 1, dims(2)
817  cz_org(2,i,j) = max( 0.0_rp, &
818  cz_org(3,i,j) &
819  * ( log(pres_org(2,i,j)/pres_org(1,i,j)) ) &
820  / ( log(pres_org(3,i,j)/pres_org(1,i,j)) ) )
821  end do
822  end do
823  end if
824 
825  velz_org = 0.0_rp
826 
827 
828  ! check verticaly extrapolated data in outer model
829  do j = 1, dims(3)
830  do i = 1, dims(2)
831  do k = 3, dims(1)+2
832 
833  if( pres_org(k,i,j)>pres_org(2,i,j) )then ! if Pressure is larger than Surface pressure
834  !velz_org(k,i,j)=velz_org(2,i,j)
835  velx_org(k,i,j)=velx_org(2,i,j)
836  vely_org(k,i,j)=vely_org(2,i,j)
837  temp_org(k,i,j)=temp_org(2,i,j)
838  qtrc_org(k,i,j,:)=qtrc_org(2,i,j,:)
839  cz_org(k,i,j)=cz_org(2,i,j)
840  end if
841 
842  enddo
843  enddo
844  enddo
845 
846 
847  !do it = 1, nt
848  ! k=1 ; j=int(dims(3)/2) ; i=int(dims(2)/2) ; iq = 1
849  ! write(*,*) "read 3D grads data",i,j,k
850  ! write(*,*) "lon_org ",lon_org (i,j)/D2R
851  ! write(*,*) "lat_org ",lat_org (i,j)/D2R
852  ! write(*,*) "pres_org ",pres_org (k,i,j)
853  ! write(*,*) "usfc_org ",usfc_org (i,j)
854  ! write(*,*) "vsfc_org ",vsfc_org (i,j)
855  ! write(*,*) "tsfc_org ",tsfc_org (i,j)
856  ! write(*,*) "qsfc_org ",qsfc_org (i,j,iq)
857  ! write(*,*) "rhsfc_org ",rhsfc_org (i,j)
858  ! write(*,*) "velx_org ",velx_org (k,i,j)
859  ! write(*,*) "vely_org ",vely_org (k,i,j)
860  ! write(*,*) "temp_org ",temp_org (k,i,j)
861  ! write(*,*) "hgt_org ",hgt_org (k,i,j)
862  ! write(*,*) "qtrc_org ",qtrc_org (k,i,j,iq)
863  ! write(*,*) "rhprs_org ",rhprs_org (k,i,j)
864  !enddo
865 
866  return
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:60
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
integer function, public io_get_available_fid()
search & get available file ID
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:93
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:71
subroutine, public log(type, message)
Definition: dc_log.f90:133
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public rp
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,
logical, intent(out)  use_waterratio,
logical, intent(in)  use_file_landwater,
character(len=*), intent(in)  basename 
)

Land Setup.

Definition at line 876 of file mod_realinput_grads.f90.

References scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_stdio::io_lnml, scale_process::prc_mpistop(), and read_namelist().

Referenced by mod_realinput::parentatomsetup().

876  implicit none
877 
878  integer, intent(out) :: ldims(3)
879  logical, intent(out) :: use_waterratio
880  logical, intent(in) :: use_file_landwater ! use landwater data from files
881  character(len=*), intent(in) :: basename
882 
883  integer :: ielem
884  integer :: k, n
885 
886  integer :: ierr
887  !---------------------------------------------------------------------------
888 
889  if( io_l ) write(io_fid_log,*) '+++ Real Case/Land Input File Type: GrADS format'
890 
891  !--- initialization
892  use_waterratio = .false.
893 
894  if ( len_trim(basename) == 0 ) then
895  if( io_l ) write(io_fid_log,*) &
896  'xxx "BASEMAAME" is not specified in "PARAM_MKINIT_REAL_ATOMS"!',trim(basename)
897  call prc_mpistop
898  endif
899 
900  !--- read namelist
901  io_fid_grads_nml = io_get_available_fid()
902  open( io_fid_grads_nml, &
903  file = trim(basename), &
904  form = 'formatted', &
905  status = 'old', &
906  action = 'read', &
907  iostat = ierr )
908  if ( ierr /= 0 ) then
909  if( io_l ) write(io_fid_log,*) 'xxx Input file is not found! ', trim(basename)
910  call prc_mpistop
911  endif
912 
913  read(io_fid_grads_nml,nml=nml_grads_grid,iostat=ierr)
914  if( ierr /= 0 ) then !--- missing or fatal error
915  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in nml_grads_grid in ', trim(basename),'. Check!'
916  call prc_mpistop
917  endif
918  if( io_lnml ) write(io_fid_log,nml=nml_grads_grid)
919 
920  ! land
921  ldims(1) = outer_nl ! soil_layers_stag
922  if(outer_nx_sfc > 0)then
923  ldims(2) = outer_nx_sfc
924  else
925  ldims(2) = outer_nx
926  outer_nx_sfc = outer_nx
927  endif
928  if(outer_ny_sfc > 0)then
929  ldims(3) = outer_ny_sfc
930  else
931  ldims(3) = outer_ny
932  outer_ny_sfc = outer_ny
933  endif
934 
935  allocate( gland2d( ldims(2), ldims(3) ) )
936  allocate( gland3d( ldims(2), ldims(3), ldims(1) ) )
937 
938  call read_namelist( &
939  grads_item(:,2), & ! (out)
940  grads_fname(:,2), & ! (out)
941  grads_dtype(:,2), & ! (out)
942  grads_swpoint(:,2), & ! (out)
943  grads_dd(:,2), & ! (out)
944  grads_lnum(:,2), & ! (out)
945  grads_lvars(:,:,2), & ! (out)
946  grads_startrec(:,2), & ! (out)
947  grads_totalrec(:,2), & ! (out)
948  grads_knum(:,2), & ! (out)
949  grads_yrev(:,2), & ! (out)
950  grads_fendian(:,2), & ! (out)
951  grads_missval(:,2), & ! (out)
952  data_available(:,2), & ! (out)
953  item_list_land, & ! (in)
954  num_item_list_land, & ! (in)
955  basename, & ! (in)
956  io_fid_grads_nml ) ! (in)
957 
958  do ielem = 1, num_item_list_land
959  item = item_list_land(ielem)
960  !--- check data
961  select case(trim(item))
962  case('TOPO','lsmask')
963  if ( .not. data_available(ielem,2) ) then
964  if (io_l) write(io_fid_log,*) 'warning: ',trim(item),' is not found & not used.'
965  cycle
966  endif
967  case('lon', 'lat', 'lon_sfc', 'lat_sfc')
968  cycle
969  case('SMOISVC', 'SMOISDS')
970  if ( use_file_landwater ) then
971  if (.not. data_available(il_smoisvc,2) .and. .not. data_available(il_smoisds,2)) then
972  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : ',trim(item_list_land(ielem))
973  call prc_mpistop
974  end if
975  use_waterratio = data_available(il_smoisds,2)
976  else
977  cycle
978  end if
979  case default ! llev, SKINT, STEMP
980  if ( .not. data_available(ielem,2) ) then
981  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : ',trim(item_list_land(ielem))
982  call prc_mpistop
983  endif
984  end select
985 
986  end do
987 
988  return
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer function, public io_get_available_fid()
search & get available file ID
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentlandinputgrads()

subroutine, public mod_realinput_grads::parentlandinputgrads ( real(rp), dimension (:,:,:), intent(out)  tg_org,
real(rp), dimension (:,:,:), intent(out)  strg_org,
real(rp), dimension (:,:,:), intent(out)  smds_org,
real(rp), dimension (:,:), intent(out)  lst_org,
real(rp), dimension (:,:), intent(out)  llon_org,
real(rp), dimension (:,:), intent(out)  llat_org,
real(rp), dimension (:), intent(out)  lz_org,
real(rp), dimension(:,:), intent(out)  topo_org,
real(rp), dimension(:,:), intent(out)  lmask_org,
character(len=*), intent(in)  basename_num,
integer, dimension(3), intent(in)  ldims,
logical, intent(in)  use_file_landwater,
integer, intent(in)  nt 
)

grads data

Definition at line 1005 of file mod_realinput_grads.f90.

References scale_const::const_d2r, scale_const::const_eps, scale_const::const_tem00, scale_const::const_undef, scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, and scale_process::prc_mpistop().

Referenced by mod_realinput::parentatomsetup().

1005  use scale_const, only: &
1006  undef => const_undef, &
1007  d2r => const_d2r, &
1008  tem00 => const_tem00, &
1009  eps => const_eps
1010  implicit none
1011 
1012  real(RP), intent(out) :: tg_org (:,:,:)
1013  real(RP), intent(out) :: strg_org (:,:,:)
1014  real(RP), intent(out) :: smds_org (:,:,:)
1015  real(RP), intent(out) :: lst_org (:,:)
1016  real(RP), intent(out) :: llon_org (:,:)
1017  real(RP), intent(out) :: llat_org (:,:)
1018  real(RP), intent(out) :: lz_org (:)
1019  real(RP), intent(out) :: topo_org(:,:)
1020  real(RP), intent(out) :: lmask_org(:,:)
1021 
1022  character(len=*), intent(in) :: basename_num
1023  integer, intent(in) :: ldims(3)
1024  logical, intent(in) :: use_file_landwater ! use land water data from files
1025  integer, intent(in) :: nt
1026  ! ----------------------------------------------------------------
1027 
1029  character(len=H_LONG) :: gfile
1030 
1031  real(RP) :: qvsat, qm
1032 
1033  integer :: i, j, k, ielem, n
1034  integer :: ierr
1035 
1036  !---------------------------------------------------------------------------
1037 
1038  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[LandInputGrADS]'
1039 
1040 
1041  loop_inputlandgrads : do ielem = 1, num_item_list_land
1042 
1043  item = item_list_land(ielem)
1044 
1045  dtype = grads_dtype(ielem,2)
1046  fname = grads_fname(ielem,2)
1047  lnum = grads_lnum(ielem,2)
1048  missval = grads_missval(ielem,2)
1049 
1050  if ( grads_knum(ielem,2) > 0 )then
1051  knum = grads_knum(ielem,2)
1052  else
1053  knum = ldims(1)
1054  endif
1055 
1056  select case (trim(dtype))
1057  case("linear")
1058  swpoint = grads_swpoint(ielem,2)
1059  dd = grads_dd(ielem,2)
1060  if( (abs(swpoint-large_number_one)<eps).or.(abs(dd-large_number_one)<eps) )then
1061  write(*,*) 'xxx "swpoint" is required in grads namelist! ',swpoint
1062  write(*,*) 'xxx "dd" is required in grads namelist! ',dd
1063  call prc_mpistop
1064  endif
1065  case("levels")
1066  if ( lnum < 0 )then
1067  write(*,*) 'xxx "lnum" in grads namelist is required for levels data! '
1068  call prc_mpistop
1069  endif
1070  do k=1, lnum
1071  lvars(k)=grads_lvars(k,ielem,2)
1072  enddo
1073  if(abs(lvars(1)-large_number_one)<eps)then
1074  write(*,*) 'xxx "lvars" must be specified in grads namelist for levels data!',(lvars(k),k=1,lnum)
1075  call prc_mpistop
1076  endif
1077  case("map")
1078  startrec = grads_startrec(ielem,2)
1079  totalrec = grads_totalrec(ielem,2)
1080  fendian = grads_fendian(ielem,2)
1081  yrev = grads_yrev(ielem,2)
1082  if( (startrec<0).or.(totalrec<0) )then
1083  write(*,*) 'xxx "startrec" is required in grads namelist! ',startrec
1084  write(*,*) 'xxx "totalrec" is required in grads namelist! ',totalrec
1085  call prc_mpistop
1086  endif
1087  ! get file_io
1088  if(io_fid_grads_data < 0)then
1089  io_fid_grads_data = io_get_available_fid()
1090  endif
1091  gfile=trim(fname)//trim(basename_num)//'.grd'
1092  if( len_trim(fname)==0 )then
1093  write(*,*) 'xxx "fname" is required in grads namelist for map data! ',trim(fname)
1094  call prc_mpistop
1095  endif
1096  end select
1097 
1098  ! read data
1099  select case (trim(item))
1100  case("lsmask")
1101  if ( data_available(il_lsmask,2) ) then
1102  if ( trim(dtype) == "map" ) then
1103  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1104  lmask_org(:,:) = real(gland2D(:,:), kind=rp)
1105  endif
1106  else
1107  lmask_org = undef
1108  end if
1109  case("lon")
1110  if ( .not. data_available(il_lon_sfc,2) ) then
1111  if ( ldims(2).ne.outer_nx .or. ldims(3).ne.outer_ny ) then
1112  write(*,*) 'xxx namelist of "lon_sfc" is not found in grads namelist!'
1113  write(*,*) 'xxx dimension is different: outer_nx and outer_nx_sfc! ', outer_nx, ldims(2)
1114  write(*,*) ' : outer_ny and outer_ny_sfc! ', outer_ny, ldims(3)
1115  call prc_mpistop
1116  end if
1117  if ( trim(dtype) == "linear" ) then
1118  do j = 1, ldims(3)
1119  do i = 1, ldims(2)
1120  llon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * d2r
1121  enddo
1122  enddo
1123  else if ( trim(dtype) == "map" ) then
1124  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1125  llon_org(:,:) = real(gland2D(:,:), kind=RP) * d2r
1126  endif
1127  end if
1128  case("lon_sfc")
1129  if ( .not. data_available(il_lon_sfc,2) ) cycle
1130  if ( trim(dtype) == "linear" ) then
1131  do j = 1, ldims(3)
1132  do i = 1, ldims(2)
1133  llon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * d2r
1134  enddo
1135  enddo
1136  else if ( trim(dtype) == "map" ) then
1137  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1138  llon_org(:,:) = real(gland2D(:,:), kind=RP) * d2r
1139  endif
1140  case("lat")
1141  if ( .not. data_available(il_lat_sfc,2) ) then
1142  if ( ldims(2).ne.outer_nx .or. ldims(3).ne.outer_ny ) then
1143  write(*,*) 'xxx namelist of "lat_sfc" is not found in grads namelist!'
1144  write(*,*) 'xxx dimension is different: outer_nx and outer_nx_sfc! ', outer_nx, ldims(2)
1145  write(*,*) ' : outer_ny and outer_ny_sfc! ', outer_nx, ldims(3)
1146  call prc_mpistop
1147  end if
1148  if ( trim(dtype) == "linear" ) then
1149  do j = 1, ldims(3)
1150  do i = 1, ldims(2)
1151  llat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * d2r
1152  enddo
1153  enddo
1154  else if ( trim(dtype) == "map" ) then
1155  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1156  llat_org(:,:) = real(gland2D(:,:), kind=RP) * d2r
1157  endif
1158  end if
1159  case("lat_sfc")
1160  if ( .not. data_available(il_lat_sfc,2) ) cycle
1161  if ( trim(dtype) == "linear" ) then
1162  do j = 1, ldims(3)
1163  do i = 1, ldims(2)
1164  llat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * d2r
1165  enddo
1166  enddo
1167  else if ( trim(dtype) == "map" ) then
1168  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,1,item,startrec,totalrec,yrev,gland2d)
1169  llat_org(:,:) = real(gland2D(:,:), kind=RP) * d2r
1170  endif
1171  case("llev")
1172  if(ldims(1)/=knum)then
1173  write(*,*) 'xxx "knum" must be equal to outer_nl for llev. knum:',knum,'> outer_nl:',ldims(1)
1174  call prc_mpistop
1175  endif
1176  if ( trim(dtype) == "levels" ) then
1177  if(ldims(1)/=lnum)then
1178  write(*,*) 'xxx lnum must be same as the outer_nl for llev! ',ldims(1),lnum
1179  call prc_mpistop
1180  endif
1181  do k = 1, ldims(1)
1182  lz_org(k) = real(lvars(k), kind=rp)
1183  enddo
1184 ! else if ( trim(dtype) == "map" ) then
1185 ! call read_grads_file_3d(io_fid_grads_data,gfile,ldims(2),ldims(3),ldims(1),nt,item,startrec,totalrec,yrev,gland)
1186 ! do j = 1, ldims(3)
1187 ! do i = 1, ldims(2)
1188 ! do k = 1, ldims(1)
1189 ! lz_org(k,i,j) = real(gland(i,j,k), kind=RP)
1190 ! enddo
1191 ! enddo
1192 ! enddo
1193  endif
1194  case('STEMP')
1195  if(ldims(1)/=knum)then
1196  write(*,*) 'xxx The number of levels for STEMP must be same as llevs! ',ldims(1),knum
1197  call prc_mpistop
1198  endif
1199  if ( trim(dtype) == "map" ) then
1200  call read_grads_file_3d(io_fid_grads_data,gfile,ldims(2),ldims(3),ldims(1),nt,item,startrec,totalrec,yrev,gland3d)
1201  do j = 1, ldims(3)
1202  do i = 1, ldims(2)
1203  do k = 1, ldims(1)
1204  if ( abs(gland3d(i,j,k)-missval) < eps ) then
1205  tg_org(k,i,j) = undef
1206  else
1207  tg_org(k,i,j) = real(gland3D(i,j,k), kind=rp)
1208  end if
1209  enddo
1210  enddo
1211  enddo
1212  endif
1213  case('SMOISVC')
1214  if ( use_file_landwater ) then
1215  if(ldims(1)/=knum)then
1216  write(*,*) 'xxx The number of levels for SMOISVC must be same as llevs! ',ldims(1),knum
1217  call prc_mpistop
1218  endif
1219  if ( trim(dtype) == "map" ) then
1220  call read_grads_file_3d(io_fid_grads_data,gfile,ldims(2),ldims(3),ldims(1),nt,item,startrec,totalrec,yrev,gland3d)
1221  do j = 1, ldims(3)
1222  do i = 1, ldims(2)
1223  do k = 1, ldims(1)
1224  if ( abs(gland3d(i,j,k)-missval) < eps ) then
1225  strg_org(k,i,j) = undef
1226  else
1227  strg_org(k,i,j) = real(gland3D(i,j,k), kind=rp)
1228  end if
1229  enddo
1230  enddo
1231  enddo
1232  endif
1233  endif
1234  case('SMOISDS')
1235  if ( use_file_landwater ) then
1236  if(ldims(1)/=knum)then
1237  write(*,*) 'xxx The number of levels for SMOISDS must be same as llevs! ',ldims(1),knum
1238  call prc_mpistop
1239  endif
1240  if ( trim(dtype) == "map" ) then
1241  call read_grads_file_3d(io_fid_grads_data,gfile,ldims(2),ldims(3),ldims(1),nt,item,startrec,totalrec,yrev,gland3d)
1242  do j = 1, ldims(3)
1243  do i = 1, ldims(2)
1244  do k = 1, ldims(1)
1245  if ( abs(gland3d(i,j,k)-missval) < eps ) then
1246  smds_org(k,i,j) = undef
1247  else
1248  smds_org(k,i,j) = real(gland3D(i,j,k), kind=rp)
1249  end if
1250  enddo
1251  enddo
1252  enddo
1253  endif
1254  endif
1255  case('SKINT')
1256  if ( trim(dtype) == "map" ) then
1257  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,nt,item,startrec,totalrec,yrev,gland2d)
1258  do j = 1, ldims(3)
1259  do i = 1, ldims(2)
1260  if ( abs(gland2d(i,j)-missval) < eps ) then
1261  lst_org(i,j) = undef
1262  else
1263  lst_org(i,j) = real(gland2D(i,j), kind=rp)
1264  end if
1265  enddo
1266  enddo
1267  endif
1268  case('TOPO')
1269  if ( data_available(il_topo,2) ) then
1270  if ( trim(dtype) == "map" ) then
1271  call read_grads_file_2d(io_fid_grads_data,gfile,ldims(2),ldims(3),1,nt,item,startrec,totalrec,yrev,gland2d)
1272  do j = 1, ldims(3)
1273  do i = 1, ldims(2)
1274  if ( abs(gland2d(i,j)-missval) < eps ) then
1275  topo_org(i,j) = undef
1276  else
1277  topo_org(i,j) = real(gland2D(i,j), kind=rp)
1278  end if
1279  enddo
1280  enddo
1281  endif
1282  else
1283  topo_org = undef
1284  endif
1285  end select
1286  enddo loop_inputlandgrads
1287 
1288  !do it = 1, nt
1289  ! i=int(ldims(2)/2) ; j=int(ldims(3)/2)
1290  ! write(*,*) "read 2D grads data",ldims(2),ldims(3),i,j,it
1291  ! write(*,*) "lon_org ",lon_org (i,j)
1292  ! write(*,*) "lat_org ",lat_org (i,j)
1293  ! write(*,*) "lst_org ",lst_org(i,j)
1294  ! do k=1,dims(7)
1295  ! write(*,*) "tg_org ",tg_org (k,i,j)," k= ",k
1296  ! write(*,*) "strg_org ",strg_org (k,i,j)," k= ",k
1297  ! enddo
1298  !enddo
1299 
1300  return
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:95
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
real(rp), public const_undef
Definition: scale_const.F90:43
integer function, public io_get_available_fid()
search & get available file ID
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public rp
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,
character(len=*), intent(in)  basename 
)

Ocean Setup.

Definition at line 1309 of file mod_realinput_grads.f90.

References scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_stdio::io_lnml, scale_process::prc_mpistop(), and read_namelist().

Referenced by mod_realinput::parentatomsetup().

1309  implicit none
1310 
1311  integer, intent(out) :: odims(2)
1312  integer, intent(out) :: timelen
1313  character(len=*), intent(in) :: basename
1314 
1315  character(len=H_LONG) :: grads_ctl
1316  integer :: ielem
1317  integer :: n
1318 
1319  integer :: ierr
1320  !---------------------------------------------------------------------------
1321 
1322  if( io_l ) write(io_fid_log,*) '+++ Real Case/Ocean Input File Type: GrADS format'
1323 
1324  !--- read namelist
1325 
1326  if ( len_trim(basename) == 0 ) then
1327  grads_ctl = "namelist.grads_boundary"
1328  else
1329  grads_ctl = basename
1330  endif
1331 
1332  !--- read namelist
1333  io_fid_grads_nml = io_get_available_fid()
1334  open( io_fid_grads_nml, &
1335  file = trim(grads_ctl), &
1336  form = 'formatted', &
1337  status = 'old', &
1338  action = 'read', &
1339  iostat = ierr )
1340  if ( ierr /= 0 ) then
1341  if( io_l ) write(io_fid_log,*) 'xxx Input file is not found! ', trim(grads_ctl)
1342  call prc_mpistop
1343  endif
1344 
1345  read(io_fid_grads_nml,nml=nml_grads_grid,iostat=ierr)
1346  if( ierr /= 0 ) then !--- missing or fatal error
1347  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in nml_grads_grid in ', trim(grads_ctl),'. Check!'
1348  call prc_mpistop
1349  endif
1350  if( io_lnml ) write(io_fid_log,nml=nml_grads_grid)
1351 
1352  timelen = 0 ! will be replaced later
1353 
1354  ! sst
1355  if(outer_nx_sst > 0)then
1356  odims(1) = outer_nx_sst
1357  else if (outer_nx_sfc > 0) then
1358  odims(1) = outer_nx_sfc
1359  outer_nx_sst = outer_nx_sfc
1360  else
1361  odims(1) = outer_nx
1362  outer_nx_sst = outer_nx
1363  endif
1364  if(outer_ny_sst > 0)then
1365  odims(2) = outer_ny_sst
1366  else if(outer_ny_sfc > 0)then
1367  odims(2) = outer_ny_sfc
1368  outer_ny_sst = outer_ny_sfc
1369  else
1370  odims(2) = outer_ny
1371  outer_ny_sst = outer_ny
1372  endif
1373 
1374  allocate( gsst2d( odims(1), odims(2) ) )
1375 
1376 
1377  call read_namelist( &
1378  grads_item(:,3), & ! (out)
1379  grads_fname(:,3), & ! (out)
1380  grads_dtype(:,3), & ! (out)
1381  grads_swpoint(:,3), & ! (out)
1382  grads_dd(:,3), & ! (out)
1383  grads_lnum(:,3), & ! (out)
1384  grads_lvars(:,:,3), & ! (out)
1385  grads_startrec(:,3), & ! (out)
1386  grads_totalrec(:,3), & ! (out)
1387  grads_knum(:,3), & ! (out)
1388  grads_yrev(:,3), & ! (out)
1389  grads_fendian(:,3), & ! (out)
1390  grads_missval(:,3), & ! (out)
1391  data_available(:,3), & ! (out)
1392  item_list_ocean, & ! (in)
1393  num_item_list_ocean, & ! (in)
1394  grads_ctl, & ! (in)
1395  io_fid_grads_nml ) ! (in)
1396 
1397  do ielem = 1, num_item_list_ocean
1398  item = item_list_ocean(ielem)
1399  !--- check data
1400  select case(trim(item))
1401  case('lsmask')
1402  if ( .not. data_available(ielem,3) ) then
1403  if (io_l) write(io_fid_log,*) 'warning: ',trim(item),' is not found & not used.'
1404  cycle
1405  endif
1406  case('lon', 'lat', 'lon_sfc', 'lat_sfc', 'lon_sst', 'lat_sst')
1407  cycle
1408  case('SST')
1409  if (.not. data_available(io_sst,3) .and. .not. data_available(io_skint,3) ) then
1410  if (io_l) write(io_fid_log,*) 'xxx SST and SKINT are found in grads namelist!'
1411  call prc_mpistop
1412  endif
1413  if (.not. data_available(io_sst,3)) then
1414  if (io_l) write(io_fid_log,*) 'warning: SST is found in grads namelist. SKINT is used in place of SST.'
1415  cycle
1416  endif
1417  case('SKINT')
1418  cycle
1419  case default !
1420  if ( .not. data_available(ielem,3) ) then
1421  if( io_l ) write(io_fid_log,*) 'xxx Not found in grads namelist! : ',trim(item_list_ocean(ielem))
1422  call prc_mpistop
1423  endif
1424  end select
1425 
1426  end do
1427 
1428  return
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer function, public io_get_available_fid()
search & get available file ID
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentoceanopengrads()

subroutine, public mod_realinput_grads::parentoceanopengrads ( )

Definition at line 1433 of file mod_realinput_grads.f90.

References scale_stdio::io_fid_log, and scale_stdio::io_l.

Referenced by mod_realinput::parentatomsetup().

1433  implicit none
1434 
1435  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanOpenGrADS]'
1436 
1437  return
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the caller graph for this function:

◆ parentoceaninputgrads()

subroutine, public mod_realinput_grads::parentoceaninputgrads ( real(rp), dimension (:,:), intent(out)  tw_org,
real(rp), dimension (:,:), intent(out)  sst_org,
real(rp), dimension(:,:), intent(out)  omask_org,
real(rp), dimension (:,:), intent(out)  olon_org,
real(rp), dimension (:,:), intent(out)  olat_org,
character(len=*), intent(in)  basename_num,
integer, dimension(2), intent(in)  odims,
integer, intent(in)  nt 
)

grads data

Definition at line 1450 of file mod_realinput_grads.f90.

References scale_const::const_d2r, scale_const::const_eps, scale_const::const_tem00, scale_const::const_undef, scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, and scale_process::prc_mpistop().

Referenced by mod_realinput::parentatomsetup().

1450  use scale_const, only: &
1451  undef => const_undef, &
1452  d2r => const_d2r, &
1453  tem00 => const_tem00, &
1454  eps => const_eps
1455  implicit none
1456 
1457  real(RP), intent(out) :: tw_org (:,:)
1458  real(RP), intent(out) :: sst_org (:,:)
1459  real(RP), intent(out) :: omask_org(:,:)
1460  real(RP), intent(out) :: olon_org (:,:)
1461  real(RP), intent(out) :: olat_org (:,:)
1462 
1463  character(len=*), intent(in) :: basename_num
1464  integer, intent(in) :: odims(2)
1465  integer, intent(in) :: nt
1466  ! ----------------------------------------------------------------
1467 
1469  character(len=H_LONG) :: gfile
1470 
1471  real(RP) :: qvsat, qm
1472 
1473  integer :: i, j, ielem, n
1474  integer :: ierr
1475 
1476  !---------------------------------------------------------------------------
1477 
1478  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanInputGrADS]'
1479 
1480 
1481  loop_inputoceangrads : do ielem = 1, num_item_list_ocean
1482 
1483  item = item_list_ocean(ielem)
1484 
1485  dtype = grads_dtype(ielem,3)
1486  fname = grads_fname(ielem,3)
1487  lnum = grads_lnum(ielem,3)
1488  missval = grads_missval(ielem,3)
1489 
1490  select case (trim(dtype))
1491  case("linear")
1492  swpoint = grads_swpoint(ielem,3)
1493  dd = grads_dd(ielem,3)
1494  if( (abs(swpoint-large_number_one)<eps).or.(abs(dd-large_number_one)<eps) )then
1495  write(*,*) 'xxx "swpoint" is required in grads namelist! ',swpoint
1496  write(*,*) 'xxx "dd" is required in grads namelist! ',dd
1497  call prc_mpistop
1498  endif
1499  case("levels")
1500  write(*,*) 'xxx "lnum" in grads namelist is invalid for ocean data'
1501  call prc_mpistop
1502  case("map")
1503  startrec = grads_startrec(ielem,3)
1504  totalrec = grads_totalrec(ielem,3)
1505  fendian = grads_fendian(ielem,3)
1506  yrev = grads_yrev(ielem,3)
1507  if( (startrec<0).or.(totalrec<0) )then
1508  write(*,*) 'xxx "startrec" is required in grads namelist! ',startrec
1509  write(*,*) 'xxx "totalrec" is required in grads namelist! ',totalrec
1510  call prc_mpistop
1511  endif
1512  ! get file_io
1513  if(io_fid_grads_data < 0)then
1514  io_fid_grads_data = io_get_available_fid()
1515  endif
1516  gfile=trim(fname)//trim(basename_num)//'.grd'
1517  if( len_trim(fname)==0 )then
1518  write(*,*) 'xxx "fname" is required in grads namelist for map data! ',trim(fname)
1519  call prc_mpistop
1520  endif
1521  end select
1522 
1523  ! read data
1524  select case (trim(item))
1525  case("lsmask")
1526  if ( data_available(io_lsmask,3) ) then
1527  if ( trim(dtype) == "map" ) then
1528  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1529  omask_org(:,:) = real(gsst2D(:,:), kind=rp)
1530  endif
1531  else
1532  omask_org = undef
1533  end if
1534  case("lon")
1535  if ( .not. data_available(io_lon_sst,3) .and. .not. data_available(io_lon_sfc,3) ) then
1536  if ( odims(1).ne.outer_nx .or. odims(2).ne.outer_ny ) then
1537  write(*,*) 'xxx namelist of "lon_sst" is not found in grads namelist!'
1538  write(*,*) 'xxx dimension is different: outer_nx and outer_nx_sst! ', outer_nx, odims(1)
1539  write(*,*) ' : outer_ny and outer_ny_sst! ', outer_ny, odims(2)
1540  call prc_mpistop
1541  end if
1542  if ( trim(dtype) == "linear" ) then
1543  do j = 1, odims(2)
1544  do i = 1, odims(1)
1545  olon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * d2r
1546  enddo
1547  enddo
1548  else if ( trim(dtype) == "map" ) then
1549  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1550  olon_org(:,:) = real(gsst2D(:,:), kind=RP) * d2r
1551  endif
1552  end if
1553  case("lon_sfc")
1554  if ( .not. data_available(io_lon_sst,3) .and. data_available(io_lon_sfc,3) ) then
1555  if ( odims(1).ne.outer_nx_sfc .or. odims(2).ne.outer_ny_sfc ) then
1556  write(*,*) 'xxx namelist of "lon_sst" is not found in grads namelist!'
1557  write(*,*) 'xxx dimension is different: outer_nx_sfc and outer_nx_sst! ', outer_nx_sfc, odims(1)
1558  write(*,*) ' : outer_ny_sfc and outer_ny_sst! ', outer_ny_sfc, odims(2)
1559  call prc_mpistop
1560  end if
1561  if ( trim(dtype) == "linear" ) then
1562  do j = 1, odims(2)
1563  do i = 1, odims(1)
1564  olon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * d2r
1565  enddo
1566  enddo
1567  else if ( trim(dtype) == "map" ) then
1568  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1569  olon_org(:,:) = real(gsst2D(:,:), kind=RP) * d2r
1570  endif
1571  end if
1572  case("lon_sst")
1573  if ( .not. data_available(io_lon_sst,3) ) cycle
1574  if ( trim(dtype) == "linear" ) then
1575  do j = 1, odims(2)
1576  do i = 1, odims(1)
1577  olon_org(i,j) = real(swpoint+real(i-1)*dd, kind=RP) * d2r
1578  enddo
1579  enddo
1580  else if ( trim(dtype) == "map" ) then
1581  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1582  olon_org(:,:) = real(gsst2D(:,:), kind=RP) * d2r
1583  endif
1584  case("lat")
1585  if ( .not. data_available(io_lat_sfc,3) .and. .not. data_available(io_lat_sst,3) ) then
1586  if ( odims(1).ne.outer_nx .or. odims(2).ne.outer_ny ) then
1587  write(*,*) 'xxx namelist of "lat_sst" is not found in grads namelist!'
1588  write(*,*) 'xxx dimension is different: outer_nx and outer_nx_sst! ', outer_nx, odims(1)
1589  write(*,*) ' : outer_ny and outer_ny_sst! ', outer_ny, odims(2)
1590  call prc_mpistop
1591  end if
1592  if ( trim(dtype) == "linear" ) then
1593  do j = 1, odims(2)
1594  do i = 1, odims(1)
1595  olat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * d2r
1596  enddo
1597  enddo
1598  else if ( trim(dtype) == "map" ) then
1599  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1600  olat_org(:,:) = real(gsst2D(:,:), kind=RP) * d2r
1601  endif
1602  end if
1603  case("lat_sfc")
1604  if ( .not. data_available(io_lat_sst,3) .and. data_available(io_lat_sfc,3) ) then
1605  if ( odims(1).ne.outer_nx .or. odims(1).ne.outer_ny ) then
1606  write(*,*) 'xxx namelist of "lat_sst" is not found in grads namelist!'
1607  write(*,*) 'xxx dimension is different: outer_nx_sfc and outer_nx_sst! ', outer_nx_sfc, odims(1)
1608  write(*,*) ' : outer_ny_sfc and outer_ny_sst! ', outer_ny_sfc, odims(2)
1609  call prc_mpistop
1610  end if
1611  if ( trim(dtype) == "linear" ) then
1612  do j = 1, odims(2)
1613  do i = 1, odims(1)
1614  olat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * d2r
1615  enddo
1616  enddo
1617  else if ( trim(dtype) == "map" ) then
1618  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1619  olat_org(:,:) = real(gsst2D(:,:), kind=RP) * d2r
1620  endif
1621  end if
1622  case("lat_sst")
1623  if ( .not. data_available(io_lat_sst,3) ) cycle
1624  if ( trim(dtype) == "linear" ) then
1625  do j = 1, odims(2)
1626  do i = 1, odims(1)
1627  olat_org(i,j) = real(swpoint+real(j-1)*dd, kind=RP) * d2r
1628  enddo
1629  enddo
1630  else if ( trim(dtype) == "map" ) then
1631  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,1,item,startrec,totalrec,yrev,gsst2d)
1632  olat_org(:,:) = real(gsst2D(:,:), kind=RP) * d2r
1633  endif
1634  case('SKINT')
1635  if ( .not. data_available(io_sst,3) ) then
1636  if ( odims(1).ne.outer_nx_sfc .or. odims(2).ne.outer_ny_sfc ) then
1637  write(*,*) 'xxx dimsntion is different: outer_nx_sst/outer_nx_sfc and outer_nx_sst! ', odims(1), outer_nx_sfc
1638  write(*,*) ' : outer_ny_sst/outer_ny_sfc and outer_ny_sst! ', odims(2), outer_ny_sfc
1639  call prc_mpistop
1640  end if
1641  if ( trim(dtype) == "map" ) then
1642  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,nt,item,startrec,totalrec,yrev,gsst2d)
1643  do j = 1, odims(2)
1644  do i = 1, odims(1)
1645  if ( abs(gsst2d(i,j)-missval) < eps ) then
1646  sst_org(i,j) = undef
1647  else
1648  sst_org(i,j) = real(gsst2D(i,j), kind=rp)
1649  end if
1650  enddo
1651  enddo
1652  end if
1653  endif
1654  case('SST')
1655  if ( .not. data_available(io_sst,3) ) cycle
1656  if ( trim(dtype) == "map" ) then
1657  call read_grads_file_2d(io_fid_grads_data,gfile,odims(1),odims(2),1,nt,item,startrec,totalrec,yrev,gsst2d)
1658  do j = 1, odims(2)
1659  do i = 1, odims(1)
1660  if ( abs(gsst2d(i,j)-missval) < eps ) then
1661  sst_org(i,j) = undef
1662  else
1663  sst_org(i,j) = real(gsst2D(i,j), kind=rp)
1664  end if
1665  enddo
1666  enddo
1667  end if
1668  end select
1669  enddo loop_inputoceangrads
1670 
1671  tw_org = sst_org
1672 
1673  !do it = 1, nt
1674  ! i=int(dims(8)/2) ; j=int(dims(9)/2)
1675  ! write(*,*) "read 2D grads data",dims(8),dims(9),i,j,it
1676  ! write(*,*) "lon_org ",lon_org (i,j)
1677  ! write(*,*) "lat_org ",lat_org (i,j)
1678  ! write(*,*) "sst_org ",sst_org (i,j)
1679  ! write(*,*) "lst_org ",lst_org(i,j)
1680  ! do k=1,dims(7)
1681  ! write(*,*) "tg_org ",tg_org (k,i,j)," k= ",k
1682  ! write(*,*) "strg_org ",strg_org (k,i,j)," k= ",k
1683  ! enddo
1684  !enddo
1685 
1686  return
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:95
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
real(rp), public const_undef
Definition: scale_const.F90:43
integer function, public io_get_available_fid()
search & get available file ID
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public rp
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_namelist()

subroutine mod_realinput_grads::read_namelist ( character(len=h_short), dimension (:), intent(out)  grads_item,
character(len=h_long), dimension (:), intent(out)  grads_fname,
character(len=h_long), dimension (:), intent(out)  grads_dtype,
real(rp), dimension (:), intent(out)  grads_swpoint,
real(rp), dimension (:), intent(out)  grads_dd,
integer, dimension (:), intent(out)  grads_lnum,
real(rp), dimension (:,:), intent(out)  grads_lvars,
integer, dimension(:), intent(out)  grads_startrec,
integer, dimension(:), intent(out)  grads_totalrec,
integer, dimension (:), intent(out)  grads_knum,
character(len=h_short), dimension (:), intent(out)  grads_yrev,
character(len=h_short), dimension (:), intent(out)  grads_fendian,
real(sp), dimension (:), intent(out)  grads_missval,
logical, dimension(:), intent(out)  data_available,
character(len=h_short), dimension (:), intent(in)  item_list,
integer, intent(in)  num_item_list,
character(len=*), intent(in)  basename,
integer, intent(in)  io_fid_grads_nml 
)

Definition at line 1708 of file mod_realinput_grads.f90.

References scale_stdio::io_fid_log, scale_stdio::io_l, and scale_process::prc_mpistop().

Referenced by parentatomsetupgrads(), parentlandsetupgrads(), and parentoceansetupgrads().

1708  implicit none
1709  character(len=H_SHORT), intent(out) :: grads_item (:)
1710  character(len=H_LONG), intent(out) :: grads_fname (:)
1711  character(len=H_LONG), intent(out) :: grads_dtype (:)
1712  real(RP), intent(out) :: grads_swpoint (:)
1713  real(RP), intent(out) :: grads_dd (:)
1714  integer, intent(out) :: grads_lnum (:)
1715  real(RP), intent(out) :: grads_lvars (:,:)
1716  integer, intent(out) :: grads_startrec(:)
1717  integer, intent(out) :: grads_totalrec(:)
1718  integer, intent(out) :: grads_knum (:)
1719  character(len=H_SHORT), intent(out) :: grads_yrev (:)
1720  character(len=H_SHORT), intent(out) :: grads_fendian (:)
1721  real(SP), intent(out) :: grads_missval (:)
1722  logical, intent(out) :: data_available(:)
1723  character(len=H_SHORT), intent(in) :: item_list (:)
1724  integer, intent(in) :: num_item_list
1725  character(len=*), intent(in) :: basename
1726  integer, intent(in) :: io_fid_grads_nml
1727 
1728  integer :: grads_vars_nmax
1729  integer :: k, n, ielem, ierr
1730 
1731  namelist /grdvar/ &
1732  item, & ! necessary
1733  dtype, & ! necessary
1734  fname, & ! necessary except for linear data
1735  swpoint, & ! for linear data
1736  dd, & ! for linear data
1737  lnum, & ! for levels data
1738  lvars, & ! for levels data
1739  startrec, & ! for map data
1740  totalrec, & ! for map data
1741  missval, & ! option
1742  knum, & ! option
1743  yrev, & ! option
1744  fendian ! option
1745 
1746  ! listup variables
1747  if ( io_fid_grads_nml > 0 ) then
1748  rewind( io_fid_grads_nml )
1749  grads_vars_nmax = 0
1750  do n = 1, grads_vars_limit
1751  read(io_fid_grads_nml, nml=grdvar, iostat=ierr)
1752  if( ierr > 0 )then
1753  if( io_l ) write(io_fid_log,*) 'xxx Not appropriate names in grdvar in ', trim(basename),'. Check!'
1754  call prc_mpistop
1755  else if( ierr < 0 )then
1756  exit
1757  endif
1758  grads_vars_nmax = grads_vars_nmax + 1
1759  enddo
1760  else
1761  if( io_l ) write(io_fid_log,*) 'xxx namelist file is not open! ', trim(basename)
1762  call prc_mpistop
1763  endif
1764 
1765  if ( grads_vars_nmax > grads_vars_limit ) then
1766  if( io_l ) write(io_fid_log,*) &
1767  'xxx The number of grads vars exceeds grads_vars_limit! ',grads_vars_nmax,' >', grads_vars_limit
1768  call prc_mpistop
1769  endif
1770 
1771 
1772  ! check data availability
1773  data_available(:) = .false.
1774  do ielem = 1, num_item_list
1775  if ( io_fid_grads_nml > 0 ) rewind( io_fid_grads_nml )
1776  do n = 1, grads_vars_nmax
1777 
1778  ! set default
1779  item = ''
1780  dtype = ''
1781  fname = ''
1782  swpoint = large_number_one
1783  dd = large_number_one
1784  lnum = -99
1785  lvars = large_number_one
1786  startrec = -99
1787  totalrec = -99
1788  knum = -99
1789  yrev = 'off'
1790  fendian = 'big'
1791  missval = large_number_one
1792 
1793  ! read namelist
1794  if ( io_fid_grads_nml > 0 ) then
1795  read(io_fid_grads_nml, nml=grdvar, iostat=ierr)
1796  if( ierr /= 0 ) exit
1797  endif
1798 
1799  if(item == item_list(ielem))then
1800  grads_item(ielem) = item
1801  grads_fname(ielem) = fname
1802  grads_dtype(ielem) = dtype
1803  grads_swpoint(ielem) = swpoint
1804  grads_dd(ielem) = dd
1805  grads_lnum(ielem) = lnum
1806  do k = 1, lvars_limit
1807  grads_lvars(k,ielem) = lvars(k)
1808  enddo
1809  grads_startrec(ielem) = startrec
1810  grads_totalrec(ielem) = totalrec
1811  grads_knum(ielem) = knum
1812  grads_yrev(ielem) = yrev
1813  grads_fendian(ielem) = fendian
1814  grads_missval(ielem) = missval
1815  data_available(ielem) = .true.
1816 
1817  exit
1818  endif
1819  enddo ! n
1820  if( io_l ) write(io_fid_log,*) 'GrADS data availability ',trim(item_list(ielem)),data_available(ielem)
1821  enddo ! ielem
1822 
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the call graph for this function:
Here is the caller graph for this function: