SCALE-RM
Functions/Subroutines
mod_realinput_wrfarw Module Reference

module REAL input WRF-ARW More...

Functions/Subroutines

subroutine, public parentatmossetupwrfarw (dims, timelen, basename_org)
 Atmos Setup. More...
 
subroutine, public parentatmosopenwrfarw
 
subroutine, public parentatmosinputwrfarw (velz_org, llvelx_org, llvely_org, pres_org, temp_org, qv_org, qhyd_org, qnum_org, lon_org, lat_org, cz_org, basename, dims, it)
 
subroutine, public parentlandsetupwrfarw (ldims, basename_land)
 Land Setup. More...
 
subroutine, public parentlandinputwrfarw (tg_org, sh2o_org, lst_org, ust_org, albg_org, topo_org, lmask_org, llon_org, llat_org, lz_org, basename, ldims, use_file_landwater, it)
 
subroutine, public parentoceansetupwrfarw (odims, timelen, basename_org)
 Ocean Setup. More...
 
subroutine, public parentoceanopenwrfarw
 
subroutine, public parentoceaninputwrfarw (tw_org, sst_org, albw_org, z0w_org, omask_org, olon_org, olat_org, basename, odims, it)
 
subroutine wrf_arwpost_calc_uvmet (u_latlon, v_latlon, u_on_map, v_on_map, xlon, xlat, basename, K1, I1, J1)
 convert vector varibles from map-projected grid on wrf model to lat-lon grid More...
 

Detailed Description

module REAL input WRF-ARW

Description
read data from WRF-ARW file for real atmospheric simulations
Author
Team SCALE
NAMELIST
  • PARAM_MKINIT_REAL_WRFARW
    nametypedefault valuecomment
    WRF_FILE_TYPE logical .false. wrf filetype: T=wrfout, F=wrfrst

History Output
No history output

Function/Subroutine Documentation

◆ parentatmossetupwrfarw()

subroutine, public mod_realinput_wrfarw::parentatmossetupwrfarw ( integer, dimension(6), intent(out)  dims,
integer, intent(out)  timelen,
character(len=*), intent(in)  basename_org 
)

Atmos Setup.

Definition at line 83 of file mod_realinput_wrfarw.F90.

References scale_file::file_get_dimlength(), scale_file::file_open(), scale_io::io_fid_conf, and scale_prc::prc_abort().

Referenced by mod_realinput::realinput_surface().

83  use scale_file, only: &
84  file_open, &
86  implicit none
87 
88  integer, intent(out) :: dims(6)
89  integer, intent(out) :: timelen
90  character(len=*), intent(in) :: basename_org
91 
92  logical :: wrf_file_type = .false. ! wrf filetype: T=wrfout, F=wrfrst
93 
94  namelist / param_mkinit_real_wrfarw / &
95  wrf_file_type
96 
97  integer :: fid
98  integer :: ierr
99  logical :: error
100  !---------------------------------------------------------------------------
101 
102  log_newline
103  log_info("ParentAtmosSetupWRFARW",*) 'Setup'
104 
105  !--- read namelist
106  rewind(io_fid_conf)
107  read(io_fid_conf,nml=param_mkinit_real_wrfarw,iostat=ierr)
108  if( ierr > 0 ) then
109  log_error("ParentAtmosSetupWRFARW",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_WRFARW. Check!'
110  call prc_abort
111  endif
112  log_nml(param_mkinit_real_wrfarw)
113 
114  call file_open( basename_org, fid, rankid=myrank, single=.true., postfix="" )
115 
116  call file_get_dimlength( fid, "bottom_top", dims(1) )
117  call file_get_dimlength( fid, "west_east", dims(2) )
118  call file_get_dimlength( fid, "south_north", dims(3) )
119  call file_get_dimlength( fid, "bottom_top_stag", dims(4) )
120  call file_get_dimlength( fid, "west_east_stag", dims(5) )
121  call file_get_dimlength( fid, "south_north_stag", dims(6) )
122 
123  call file_get_dimlength( fid, "Time", timelen, error=error )
124  if ( error ) call file_get_dimlength( fid, "time", timelen, error=error)
125  if ( error ) timelen = 0
126 
127  if ( wrf_file_type ) then
128  wrfout = .true.
129  log_info("ParentAtmosSetupWRFARW",*) 'WRF-ARW FILE-TYPE: WRF History Output'
130  else
131  wrfout = .false.
132  log_info("ParentAtmosSetupWRFARW",*) 'WRF-ARW FILE-TYPE: WRF Restart'
133  endif
134 
135 
136  allocate( read_xy(dims(2),dims(3)) )
137  allocate( read_xyz(dims(2),dims(3),dims(1)) )
138  allocate( read_xyw(dims(2),dims(3),dims(4)) )
139 
140  allocate( p_org(dims(1),dims(2),dims(3)) )
141  allocate( pb_org(dims(1),dims(2),dims(3)) )
142  allocate( ph_org(dims(4),dims(2),dims(3)) )
143  allocate( phb_org(dims(4),dims(2),dims(3)) )
144 
145  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
module file
Definition: scale_file.F90:15
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
subroutine, public file_get_dimlength(fid, dimname, len, error)
get length of dimension
Definition: scale_file.F90:525
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentatmosopenwrfarw()

subroutine, public mod_realinput_wrfarw::parentatmosopenwrfarw ( )

Definition at line 150 of file mod_realinput_wrfarw.F90.

Referenced by mod_realinput::realinput_surface().

150  implicit none
151 
152  return
Here is the caller graph for this function:

◆ parentatmosinputwrfarw()

subroutine, public mod_realinput_wrfarw::parentatmosinputwrfarw ( real(rp), dimension(:,:,:), intent(out)  velz_org,
real(rp), dimension(:,:,:), intent(out)  llvelx_org,
real(rp), dimension(:,:,:), intent(out)  llvely_org,
real(rp), dimension(:,:,:), intent(out)  pres_org,
real(rp), dimension(:,:,:), intent(out)  temp_org,
real(rp), dimension(:,:,:), intent(out)  qv_org,
real(rp), dimension(:,:,:,:), intent(out)  qhyd_org,
real(rp), dimension(:,:,:,:), intent(out)  qnum_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,
integer, dimension(6), intent(in)  dims,
integer, intent(in)  it 
)

Definition at line 171 of file mod_realinput_wrfarw.F90.

References scale_const::const_d2r, scale_const::const_grav, scale_const::const_laps, scale_const::const_rdry, scale_file::file_open(), scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hg, scale_atmos_hydrometeor::i_hi, scale_atmos_hydrometeor::i_hr, scale_atmos_hydrometeor::i_hs, scale_atmos_hydrometeor::n_hyd, and wrf_arwpost_calc_uvmet().

Referenced by mod_realinput::realinput_surface().

171  use scale_const, only: &
172  d2r => const_d2r, &
173  laps => const_laps, &
174  rdry => const_rdry, &
175  grav => const_grav
176  use scale_file, only: &
177  file_open, &
178  file_read
179  use scale_atmos_hydrometeor, only: &
180  n_hyd, &
181  i_hc, &
182  i_hr, &
183  i_hi, &
184  i_hs, &
185  i_hg
186  implicit none
187 
188  real(RP), intent(out) :: velz_org(:,:,:)
189  real(RP), intent(out) :: llvelx_org(:,:,:)
190  real(RP), intent(out) :: llvely_org(:,:,:)
191  real(RP), intent(out) :: pres_org(:,:,:)
192  real(RP), intent(out) :: temp_org(:,:,:)
193  real(RP), intent(out) :: qv_org(:,:,:)
194  real(RP), intent(out) :: qhyd_org(:,:,:,:)
195  real(RP), intent(out) :: qnum_org(:,:,:,:)
196  real(RP), intent(out) :: lon_org(:,:)
197  real(RP), intent(out) :: lat_org(:,:)
198  real(RP), intent(out) :: cz_org(:,:,:)
199  character(len=*), intent(in) :: basename
200  integer, intent(in) :: dims(6)
201  integer, intent(in) :: it
202 
203  ! k, i, j
204  real(RP) :: velx_org(dims(1)+2,dims(2),dims(3))
205  real(RP) :: vely_org(dims(1)+2,dims(2),dims(3))
206  real(RP) :: pott_org(dims(1)+2,dims(2),dims(3))
207  real(RP) :: topo_org( dims(2),dims(3))
208  real(RP) :: geof_org(dims(4) ,dims(2),dims(3))
209 
210 
211  ! i, j, k
212  real(RP) :: velzs_org(dims(2),dims(3),dims(4))
213  real(RP) :: velxs_org(dims(5),dims(3),dims(1))
214  real(RP) :: velys_org(dims(2),dims(6),dims(1))
215 
216  real(RP) :: dens
217  real(RP) :: qtot
218 
219  integer :: k, i, j, iq
220 
221  character(len=H_MID) :: varname_t
222  character(len=H_MID) :: varname_w
223  character(len=H_MID) :: varname_u
224  character(len=H_MID) :: varname_v
225 
226  integer :: fid
227  !---------------------------------------------------------------------------
228 
229  if ( wrfout ) then
230  varname_t = "T"
231  varname_w = "W"
232  varname_u = "U"
233  varname_v = "V"
234  else
235  varname_t = "T_1"
236  varname_w = "W_1"
237  varname_u = "U_1"
238  varname_v = "V_1"
239  endif
240 
241 
242  call file_open( basename, fid, rankid=myrank, single=.true., postfix="" )
243 
244  call file_read( fid, "XLAT", lat_org(:,:), step=it )
245  lat_org(:,:) = lat_org(:,:) * d2r
246 
247  call file_read( fid, "XLONG", lon_org(:,:), step=it )
248  lon_org(:,:) = lon_org(:,:) * d2r
249 
250  call file_read( fid, "HGT", topo_org(:,:), step=it )
251 
252  call file_read( fid, "PH", read_xyw(:,:,:), step=it )
253  do j = 1, dims(3)
254  do i = 1, dims(2)
255  do k = 1, dims(4)
256  ph_org(k,i,j) = read_xyw(i,j,k)
257  end do
258  end do
259  end do
260 
261  call file_read( fid, "PHB", read_xyw(:,:,:), step=it )
262  do j = 1, dims(3)
263  do i = 1, dims(2)
264  do k = 1, dims(4)
265  phb_org(k,i,j) = read_xyw(i,j,k)
266  end do
267  end do
268  end do
269 
270  call file_read( fid, "P", read_xyz(:,:,:), step=it )
271  do j = 1, dims(3)
272  do i = 1, dims(2)
273  do k = 1, dims(1)
274  p_org(k,i,j) = read_xyz(i,j,k)
275  end do
276  end do
277  end do
278 
279  call file_read( fid, "PB", read_xyz(:,:,:), step=it )
280  do j = 1, dims(3)
281  do i = 1, dims(2)
282  do k = 1, dims(1)
283  pb_org(k,i,j) = read_xyz(i,j,k)
284  end do
285  end do
286  end do
287 
288  call file_read( fid, varname_w, velzs_org(:,:,:), step=it )
289 
290  call file_read( fid, varname_u, velxs_org(:,:,:), step=it )
291 
292  call file_read( fid, varname_v, velys_org(:,:,:), step=it )
293 
294 
295  ! from half level to full level
296  do j = 1, dims(3)
297  do i = 1, dims(2)
298  do k = 1, dims(1)
299  velz_org(k+2,i,j) = ( velzs_org(i,j,k) + velzs_org(i,j,k+1) ) * 0.5_rp
300  velx_org(k+2,i,j) = ( velxs_org(i,j,k) + velxs_org(i+1,j,k) ) * 0.5_rp
301  vely_org(k+2,i,j) = ( velys_org(i,j,k) + velys_org(i,j+1,k) ) * 0.5_rp
302  end do
303  velz_org(1:2,i,j) = 0.0_rp
304  velx_org(1:2,i,j) = 0.0_rp
305  vely_org(1:2,i,j) = 0.0_rp
306  end do
307  end do
308 
309  call wrf_arwpost_calc_uvmet( llvelx_org, llvely_org, & ! (out)
310  velx_org, vely_org, & ! (in)
311  lon_org, lat_org, & ! (in)
312  basename, & ! (in)
313  dims(1)+2, dims(2), dims(3) ) ! (in)
314 
315  qhyd_org(:,:,:,:) = 0.0_rp
316 
317  call file_read( fid, "Q2", read_xy(:,:), step=it )
318  do j = 1, dims(3)
319  do i = 1, dims(2)
320  qv_org(1,i,j) = read_xy(i,j)
321  qv_org(2,i,j) = read_xy(i,j)
322  end do
323  end do
324 
325  call file_read( fid, "QVAPOR", read_xyz(:,:,:), step=it )
326  do j = 1, dims(3)
327  do i = 1, dims(2)
328  do k = 1, dims(1)
329  qv_org(k+2,i,j) = read_xyz(i,j,k)
330  end do
331  end do
332  end do
333 
334 
335  call file_read( fid, "QCLOUD", read_xyz(:,:,:), step=it, allow_missing=.true. )
336  do j = 1, dims(3)
337  do i = 1, dims(2)
338  do k = 1, dims(1)
339  qhyd_org(k+2,i,j,i_hc) = read_xyz(i,j,k)
340  end do
341  end do
342  end do
343 
344  call file_read( fid, "QRAIN", read_xyz(:,:,:), step=it, allow_missing=.true. )
345  do j = 1, dims(3)
346  do i = 1, dims(2)
347  do k = 1, dims(1)
348  qhyd_org(k+2,i,j,i_hc) = read_xyz(i,j,k)
349  end do
350  end do
351  end do
352 
353  call file_read( fid, "QICE", read_xyz(:,:,:), step=it, allow_missing=.true. )
354  do j = 1, dims(3)
355  do i = 1, dims(2)
356  do k = 1, dims(1)
357  qhyd_org(k+2,i,j,i_hi) = read_xyz(i,j,k)
358  end do
359  end do
360  end do
361 
362  call file_read( fid, "QSNOW", read_xyz(:,:,:), step=it, allow_missing=.true. )
363  do j = 1, dims(3)
364  do i = 1, dims(2)
365  do k = 1, dims(1)
366  qhyd_org(k+2,i,j,i_hs) = read_xyz(i,j,k)
367  end do
368  end do
369  end do
370 
371  call file_read( fid, "QGRAUP", read_xyz(:,:,:), step=it, allow_missing=.true. )
372  do j = 1, dims(3)
373  do i = 1, dims(2)
374  do k = 1, dims(1)
375  qhyd_org(k+2,i,j,i_hg) = read_xyz(i,j,k)
376  end do
377  end do
378  end do
379 
380 
381  ! convert mixing ratio to specific ratio
382  do j = 1, dims(3)
383  do i = 1, dims(2)
384  do k = 1, dims(1)+2
385  qtot = qv_org(k,i,j)
386  do iq = 1, n_hyd
387  qtot = qtot + qhyd_org(k,i,j,iq)
388  end do
389  qv_org(k,i,j) = qv_org(k,i,j) / ( 1.0_rp + qtot )
390  do iq = 1, n_hyd
391  qhyd_org(k,i,j,iq) = qhyd_org(k,i,j,iq) / ( 1.0_rp + qtot )
392  end do
393  end do
394  end do
395  end do
396 
397  call file_read( fid, "NC", read_xyz(:,:,:), step=it, allow_missing=.true. )
398  do j = 1, dims(3)
399  do i = 1, dims(2)
400  do k = 1, dims(1)
401  qnum_org(k+2,i,j,i_hc) = read_xyz(i,j,k)
402  end do
403  end do
404  end do
405 
406  call file_read( fid, "NR", read_xyz(:,:,:), step=it, allow_missing=.true. )
407  do j = 1, dims(3)
408  do i = 1, dims(2)
409  do k = 1, dims(1)
410  qnum_org(k+2,i,j,i_hr) = read_xyz(i,j,k)
411  end do
412  end do
413  end do
414 
415  call file_read( fid, "NI", read_xyz(:,:,:), step=it, allow_missing=.true. )
416  do j = 1, dims(3)
417  do i = 1, dims(2)
418  do k = 1, dims(1)
419  qnum_org(k+2,i,j,i_hi) = read_xyz(i,j,k)
420  end do
421  end do
422  end do
423 
424  call file_read( fid, "NS", read_xyz(:,:,:), step=it, allow_missing=.true. )
425  do j = 1, dims(3)
426  do i = 1, dims(2)
427  do k = 1, dims(1)
428  qnum_org(k+2,i,j,i_hs) = read_xyz(i,j,k)
429  end do
430  end do
431  end do
432 
433  call file_read( fid, "NG", read_xyz(:,:,:), step=it, allow_missing=.true. )
434  do j = 1, dims(3)
435  do i = 1, dims(2)
436  do k = 1, dims(1)
437  qnum_org(k+2,i,j,i_hg) = read_xyz(i,j,k)
438  end do
439  end do
440  end do
441 
442  do iq = 1, n_hyd
443  do j = 1, dims(3)
444  do i = 1, dims(2)
445  do k = 1, dims(1)+2
446  qhyd_org(k,i,j,iq) = max( qhyd_org(k,i,j,iq), 0.0_rp )
447  qnum_org(k,i,j,iq) = max( qnum_org(k,i,j,iq), 0.0_rp )
448  end do
449  end do
450  end do
451  end do
452 
453 
454  call file_read( fid, varname_t, read_xyz(:,:,:), step=it, allow_missing=.true. )
455  do j = 1, dims(3)
456  do i = 1, dims(2)
457  do k = 1, dims(1)
458  pott_org(k+2,i,j) = read_xyz(i,j,k) + t0
459  end do
460  end do
461  end do
462  call file_read( fid, "T2", read_xy(:,:), step=it, allow_missing=.true. )
463  do j = 1, dims(3)
464  do i = 1, dims(2)
465  temp_org(2,i,j) = read_xy(i,j)
466  end do
467  end do
468 
469  call file_read( fid, "PSFC", read_xy(:,:), step=it, allow_missing=.true. )
470  do j = 1, dims(3)
471  do i = 1, dims(2)
472  pres_org(2,i,j) = read_xy(i,j)
473  end do
474  end do
475 
476  do j = 1, dims(3)
477  do i = 1, dims(2)
478  do k = 3, dims(1)+2
479  pres_org(k,i,j) = p_org(k-2,i,j) + pb_org(k-2,i,j)
480  temp_org(k,i,j) = pott_org(k,i,j) * ( pres_org(k,i,j) / p0 )**rcp
481  end do
482  pott_org(2,i,j) = temp_org(2,i,j) * ( p0/pres_org(2,i,j) )**rcp
483  temp_org(1,i,j) = temp_org(2,i,j) + laps * topo_org(i,j)
484  dens = pres_org(2,i,j) / ( rdry * temp_org(2,i,j) )
485  pres_org(1,i,j) = ( pres_org(2,i,j) + grav * dens * cz_org(2,i,j) * 0.5_rp ) &
486  / ( rdry * temp_org(1,i,j) - grav * cz_org(2,i,j) * 0.5_rp ) &
487  * rdry * temp_org(1,i,j)
488  end do
489  end do
490 
491  do j = 1, dims(3)
492  do i = 1, dims(2)
493  ! convert to geopotential height to use as real height in WRF
494  do k = 1, dims(4)
495  geof_org(k,i,j) = ( ph_org(k,i,j) + phb_org(k,i,j) ) / grav
496  end do
497  ! make half level of geopotential height from face level
498  do k = 1, dims(1)
499  cz_org(k+2,i,j) = ( geof_org(k,i,j) + geof_org(k+1,i,j) ) * 0.5_rp
500  end do
501  cz_org(2,i,j) = topo_org(i,j)
502  cz_org(1,i,j) = 0.0_rp
503  end do
504  end do
505 
506 
507 #ifdef DEBUG
508  !k=1 ; i=int(dims(2)/2) ; j=int(dims(3)/2)
509  k=2 ; i=3 ; j=3
510  log_info("ParentAtmosInputWRFARW",*) "read 3D wrf data",i,j,k
511  log_info("ParentAtmosInputWRFARW",*) "lon_org ",lon_org(i,j)/d2r
512  log_info("ParentAtmosInputWRFARW",*) "lat_org ",lat_org(i,j)/d2r
513  log_info("ParentAtmosInputWRFARW",*) "cz_org ",cz_org(k,i,j)
514  log_info("ParentAtmosInputWRFARW",*) "pres_org ",pres_org(k,i,j)
515  log_info("ParentAtmosInputWRFARW",*) "velx_org ",llvelx_org(k,i,j)
516  log_info("ParentAtmosInputWRFARW",*) "vely_org ",llvely_org(k,i,j)
517  log_info("ParentAtmosInputWRFARW",*) "velz_org ",velz_org(k,i,j)
518  log_info("ParentAtmosInputWRFARW",*) "temp_org ",temp_org(k,i,j)
519  log_info("ParentAtmosInputWRFARW",*) "qv_org ",qv_org(k,i,j)
520  k=3 ; i=3 ; j=3 ; iq = 1
521  log_info("ParentAtmosInputWRFARW",*) "read 3D wrf data",i,j,k
522  log_info("ParentAtmosInputWRFARW",*) "lon_org ",lon_org(i,j)/d2r
523  log_info("ParentAtmosInputWRFARW",*) "lat_org ",lat_org(i,j)/d2r
524  log_info("ParentAtmosInputWRFARW",*) "cz_org ",cz_org(k,i,j)
525  log_info("ParentAtmosInputWRFARW",*) "pres_org ",pres_org(k,i,j)
526  log_info("ParentAtmosInputWRFARW",*) "velx_org ",llvelx_org(k,i,j)
527  log_info("ParentAtmosInputWRFARW",*) "vely_org ",llvely_org(k,i,j)
528  log_info("ParentAtmosInputWRFARW",*) "velz_org ",velz_org(k,i,j)
529  log_info("ParentAtmosInputWRFARW",*) "temp_org ",temp_org(k,i,j)
530  log_info("ParentAtmosInputWRFARW",*) "qv_org ",qv_org(k,i,j)
531 #endif
532 
533  return
integer, parameter, public i_hs
snow
integer, parameter, public i_hr
liquid water rain
integer, parameter, public i_hi
ice water cloud
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:58
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
module file
Definition: scale_file.F90:15
module atmosphere / hydrometeor
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
module CONSTANT
Definition: scale_const.F90:11
integer, parameter, public i_hc
liquid water cloud
integer, parameter, public n_hyd
integer, parameter, public i_hg
graupel
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentlandsetupwrfarw()

subroutine, public mod_realinput_wrfarw::parentlandsetupwrfarw ( integer, dimension(3), intent(out)  ldims,
character(len=*), intent(in)  basename_land 
)

Land Setup.

Definition at line 541 of file mod_realinput_wrfarw.F90.

References scale_file::file_get_dimlength(), scale_file::file_open(), scale_io::io_fid_conf, and scale_prc::prc_abort().

Referenced by mod_realinput::realinput_surface().

541  use scale_file, only: &
542  file_open, &
544  implicit none
545 
546  integer, intent(out) :: ldims(3)
547  character(len=*), intent(in) :: basename_land
548 
549  logical :: wrf_file_type = .false. ! wrf filetype: T=wrfout, F=wrfrst
550 
551  namelist / param_mkinit_real_wrfarw / &
552  wrf_file_type
553 
554  integer :: fid
555  integer :: ierr
556  !---------------------------------------------------------------------------
557 
558  log_info("ParentLandSetupWRFARW",*) 'Real Case/Atmos Input File Type: WRF-ARW'
559 
560  !--- read namelist
561  rewind(io_fid_conf)
562  read(io_fid_conf,nml=param_mkinit_real_wrfarw,iostat=ierr)
563  if( ierr > 0 ) then
564  log_error("ParentLandSetupWRFARW",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_WRFARW. Check!'
565  call prc_abort
566  endif
567  log_nml(param_mkinit_real_wrfarw)
568 
569 
570  call file_open( basename_land, fid, rankid=myrank, single=.true., postfix="" )
571 
572  call file_get_dimlength( fid, "soil_layers_stag", ldims(1) )
573  call file_get_dimlength( fid, "west_east", ldims(2) )
574  call file_get_dimlength( fid, "south_north", ldims(3) )
575 
576  if ( wrf_file_type ) then
577  wrfout = .true.
578  log_info("ParentLandSetupWRFARW",*) 'WRF-ARW FILE-TYPE: WRF History Output'
579  else
580  wrfout = .false.
581  log_info("ParentLandSetupWRFARW",*) 'WRF-ARW FILE-TYPE: WRF Restart'
582  endif
583 
584 
585  if ( .not. allocated(read_xy) ) then
586  allocate( read_xy(ldims(2),ldims(3)) )
587  end if
588 
589  allocate( read_xyl(ldims(2),ldims(3),ldims(1)) )
590 
591  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
module file
Definition: scale_file.F90:15
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
subroutine, public file_get_dimlength(fid, dimname, len, error)
get length of dimension
Definition: scale_file.F90:525
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentlandinputwrfarw()

subroutine, public mod_realinput_wrfarw::parentlandinputwrfarw ( real(rp), dimension(:,:,:), intent(out)  tg_org,
real(rp), dimension(:,:,:), intent(out)  sh2o_org,
real(rp), dimension(:,:), intent(out)  lst_org,
real(rp), dimension(:,:), intent(out)  ust_org,
real(rp), dimension(:,:,:,:), intent(out)  albg_org,
real(rp), dimension(:,:), intent(out)  topo_org,
real(rp), dimension(:,:), intent(out)  lmask_org,
real(rp), dimension(:,:), intent(out)  llon_org,
real(rp), dimension(:,:), intent(out)  llat_org,
real(rp), dimension(:), intent(out)  lz_org,
character(len=*), intent(in)  basename,
integer, dimension(3), intent(in)  ldims,
logical, intent(in)  use_file_landwater,
integer, intent(in)  it 
)

Definition at line 610 of file mod_realinput_wrfarw.F90.

References scale_const::const_d2r, scale_const::const_undef, scale_file::file_open(), scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, and scale_cpl_sfc_index::i_r_vis.

Referenced by mod_realinput::realinput_surface().

610  use scale_const, only: &
611  d2r => const_d2r, &
612  undef => const_undef
613  use scale_file, only: &
614  file_open, &
615  file_read
616  implicit none
617  real(RP), intent(out) :: tg_org(:,:,:)
618  real(RP), intent(out) :: sh2o_org(:,:,:)
619  real(RP), intent(out) :: lst_org(:,:)
620  real(RP), intent(out) :: ust_org(:,:)
621  real(RP), intent(out) :: albg_org(:,:,:,:)
622  real(RP), intent(out) :: topo_org(:,:)
623  real(RP), intent(out) :: lmask_org(:,:)
624  real(RP), intent(out) :: llon_org(:,:)
625  real(RP), intent(out) :: llat_org(:,:)
626  real(RP), intent(out) :: lz_org(:)
627  character(len=*), intent( in) :: basename
628  integer, intent( in) :: ldims(3)
629  logical, intent( in) :: use_file_landwater ! use land water data from files
630  integer, intent( in) :: it
631 
632  integer :: fid
633  integer :: k, i, j
634  !---------------------------------------------------------------------------
635 
636  call file_open( basename, fid, rankid=myrank, single=.true., postfix="" )
637 
638  call file_read( fid, "XLAT", llat_org(:,:), step=it )
639  llat_org(:,:) = llat_org(:,:) * d2r
640 
641  call file_read( fid, "XLONG", llon_org(:,:), step=it )
642  llon_org(:,:) = llon_org(:,:) * d2r
643 
644  call file_read( fid, "HGT", topo_org(:,:), step=it )
645 
646 
647  ! depth
648  call file_read( fid, "ZS", lz_org(:), step=it )
649 
650  ! land mask (1:land, 0:water)
651  call file_read( fid, "LANDMASK", lmask_org(:,:), step=it )
652 
653  ! soil temperature [K]
654  call file_read( fid, "TSLB", read_xyl(:,:,:), step=it )
655  do j = 1, ldims(3)
656  do i = 1, ldims(2)
657  do k = 1, ldims(1)
658  tg_org(k,i,j) = read_xyl(i,j,k)
659  end do
660  end do
661  end do
662 
663  ! soil liquid water [m3 m-3] (no wrfout-default)
664  if( use_file_landwater ) then
665  call file_read( fid, "SH2O", read_xyl(:,:,:), step=it, allow_missing=.true., missing_value=undef )
666  do j = 1, ldims(3)
667  do i = 1, ldims(2)
668  do k = 1, ldims(1)
669  sh2o_org(k,i,j) = read_xyl(i,j,k)
670  end do
671  end do
672  end do
673  endif
674 
675 ! ! surface runoff [mm]
676 ! call FILE_read( fid, "SFROFF", org_3D(:,:), step=it )
677 ! do j = 1, ldims(3)
678 ! do i = 1, ldims(2)
679 ! org_3D(k,i,j) = org_3D(i,j,k) * 1000.0_RP * dwatr
680 ! end do
681 ! end do
682 
683 
684  ! SURFACE SKIN TEMPERATURE [K]
685  call file_read( fid, "TSK", lst_org(:,:), step=it )
686 
687  ust_org(:,:) = lst_org(:,:)
688 
689  ! ALBEDO [-]
690  call file_read( fid, "ALBEDO", albg_org(:,:,i_r_direct ,i_r_vis), step=it )
691  albg_org(:,:,i_r_direct ,i_r_nir) = albg_org(:,:,i_r_direct ,i_r_vis)
692  albg_org(:,:,i_r_diffuse,i_r_nir) = albg_org(:,:,i_r_direct ,i_r_vis)
693  albg_org(:,:,i_r_diffuse,i_r_vis) = albg_org(:,:,i_r_direct ,i_r_vis)
694 
695  ! SURFACE EMISSIVITY [-]
696  call file_read( fid, "EMISS", read_xy(:,:), step=it )
697  do j = 1, ldims(3)
698  do i = 1, ldims(2)
699  albg_org(i,j,i_r_diffuse,i_r_ir) = 1.0_rp - read_xy(i,j)
700  end do
701  end do
702  albg_org(:,:,i_r_direct,i_r_ir) = albg_org(:,:,i_r_diffuse,i_r_ir)
703 
704 ! ! SNOW WATER EQUIVALENT [kg m-2] (no wrfout-default)
705 ! call FILE_read( fid, "SNOW", snowq_org(:,:), step=it, allow_missing=.true., missing_value=UNDEF )
706 
707 ! ! AVERAGE SNOW TEMPERATURE [C] (no wrfout-default)
708 ! call FILE_read( fid, "TSNAV", snowt_org(:,:), step=it, allow_missing=.true., missing_value=UNDEF )
709 ! do j = 1, ldims(3)
710 ! do i = 1, ldims(2)
711 ! if ( snowt_org(k,i,j) /= UNDEF ) snowt_org(k,i,j) = snowt_org(i,j,k) + TEM00
712 ! end do
713 ! end do
714 
715  return
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
real(rp), public const_undef
Definition: scale_const.F90:41
module file
Definition: scale_file.F90:15
module CONSTANT
Definition: scale_const.F90:11
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentoceansetupwrfarw()

subroutine, public mod_realinput_wrfarw::parentoceansetupwrfarw ( integer, dimension(2), intent(out)  odims,
integer, intent(out)  timelen,
character(len=*), intent(in)  basename_org 
)

Ocean Setup.

Definition at line 724 of file mod_realinput_wrfarw.F90.

References scale_file::file_get_dimlength(), scale_file::file_open(), scale_io::io_fid_conf, and scale_prc::prc_abort().

Referenced by mod_realinput::realinput_surface().

724  use scale_file, only: &
725  file_open, &
727  implicit none
728 
729  integer, intent(out) :: odims(2)
730  integer, intent(out) :: timelen
731  character(len=*), intent(in) :: basename_org
732 
733  logical :: wrf_file_type = .false. ! wrf filetype: T=wrfout, F=wrfrst
734 
735  namelist / param_mkinit_real_wrfarw / &
736  wrf_file_type
737 
738  integer :: fid
739  integer :: ierr
740  logical :: error
741  !---------------------------------------------------------------------------
742 
743  log_info("ParentOceanSetupWRFARW",*) 'Real Case/Ocean Input File Type: WRF-ARW'
744 
745  !--- read namelist
746  rewind(io_fid_conf)
747  read(io_fid_conf,nml=param_mkinit_real_wrfarw,iostat=ierr)
748  if( ierr > 0 ) then
749  log_error("ParentOceanSetupWRFARW",*) 'Not appropriate names in namelist PARAM_MKINIT_REAL_WRFARW. Check!'
750  call prc_abort
751  endif
752  log_nml(param_mkinit_real_wrfarw)
753 
754 
755  call file_open( basename_org, fid, rankid=myrank, single=.true., postfix="" )
756 
757  call file_get_dimlength(fid, "west_east", odims(1) )
758  call file_get_dimlength(fid, "south_north", odims(2) )
759 
760  call file_get_dimlength(fid, "Time", timelen, error=error )
761  if ( error ) call file_get_dimlength(fid, "time", timelen, error=error)
762  if ( error ) timelen = 0
763 
764  if ( wrf_file_type ) then
765  wrfout = .true.
766  log_info("ParentOceanSetupWRFARW",*) 'WRF-ARW FILE-TYPE: WRF History Output'
767  else
768  wrfout = .false.
769  log_info("ParentOceanSetupWRFARW",*) 'WRF-ARW FILE-TYPE: WRF Restart'
770  endif
771 
772 
773  if ( .not. allocated(read_xy) ) then
774  allocate( read_xy(odims(1),odims(2)) )
775  end if
776 
777  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
module file
Definition: scale_file.F90:15
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
subroutine, public file_get_dimlength(fid, dimname, len, error)
get length of dimension
Definition: scale_file.F90:525
Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentoceanopenwrfarw()

subroutine, public mod_realinput_wrfarw::parentoceanopenwrfarw ( )

Definition at line 782 of file mod_realinput_wrfarw.F90.

Referenced by mod_realinput::realinput_surface().

782  implicit none
783 
784  return
Here is the caller graph for this function:

◆ parentoceaninputwrfarw()

subroutine, public mod_realinput_wrfarw::parentoceaninputwrfarw ( real(rp), dimension(:,:), intent(out)  tw_org,
real(rp), dimension(:,:), intent(out)  sst_org,
real(rp), dimension(:,:,:,:), intent(out)  albw_org,
real(rp), dimension(:,:), intent(out)  z0w_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,
integer, dimension(2), intent(in)  odims,
integer, intent(in)  it 
)

Definition at line 799 of file mod_realinput_wrfarw.F90.

References scale_const::const_d2r, scale_const::const_undef, scale_file::file_open(), scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, and scale_cpl_sfc_index::i_r_vis.

Referenced by mod_realinput::realinput_surface().

799  use scale_const, only: &
800  d2r => const_d2r, &
801  undef => const_undef
802  use scale_file, only: &
803  file_open, &
804  file_read
805  implicit none
806  real(RP), intent(out) :: tw_org(:,:)
807  real(RP), intent(out) :: sst_org(:,:)
808  real(RP), intent(out) :: albw_org(:,:,:,:)
809  real(RP), intent(out) :: z0w_org(:,:)
810  real(RP), intent(out) :: omask_org(:,:)
811  real(RP), intent(out) :: olon_org(:,:)
812  real(RP), intent(out) :: olat_org(:,:)
813  character(len=*), intent( in) :: basename
814  integer, intent( in) :: odims(2)
815  integer, intent( in) :: it
816 
817  integer :: fid
818  integer :: i, j
819  !---------------------------------------------------------------------------
820 
821  call file_open( basename, fid, rankid=myrank, single=.true., postfix="" )
822 
823  call file_read( fid, "XLAT", olat_org(:,:), step=it )
824  olat_org(:,:) = olat_org(:,:) * d2r
825 
826  call file_read( fid, "XLONG", olon_org(:,:), step=it )
827  olon_org(:,:) = olon_org(:,:) * d2r
828 
829 
830  ! land mask (1:land, 0:water)
831  call file_read( fid, "LANDMASK", omask_org(:,:), step=it )
832 
833  ! SEA SURFACE TEMPERATURE [K]
834  call file_read( fid, "SST", sst_org(:,:), step=it )
835 
836  tw_org(:,:) = sst_org(:,:)
837 
838  ! ALBEDO [-]
839  call file_read( fid, "ALBEDO", albw_org(:,:,i_r_direct ,i_r_vis), step=it )
840  albw_org(:,:,i_r_direct ,i_r_nir) = albw_org(:,:,i_r_direct ,i_r_vis)
841  albw_org(:,:,i_r_diffuse,i_r_nir) = albw_org(:,:,i_r_direct ,i_r_vis)
842  albw_org(:,:,i_r_diffuse,i_r_vis) = albw_org(:,:,i_r_direct ,i_r_vis)
843 
844  ! SURFACE EMISSIVITY [-]
845  call file_read( fid, "EMISS", read_xy(:,:), step=it )
846  do j = 1, odims(2)
847  do i = 1, odims(1)
848  albw_org(i,j,i_r_diffuse,i_r_ir) = 1.0_rp - read_xy(i,j)
849  enddo
850  enddo
851  albw_org(:,:,i_r_direct,i_r_ir) = albw_org(:,:,i_r_diffuse,i_r_ir)
852 
853  ! TIME-VARYING ROUGHNESS LENGTH [m] (no wrfout-default)
854  call file_read( fid, "ZNT", z0w_org(:,:), step=it, allow_missing=.true., missing_value=undef )
855 
856 
857  return
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
real(rp), public const_undef
Definition: scale_const.F90:41
module file
Definition: scale_file.F90:15
module CONSTANT
Definition: scale_const.F90:11
Here is the call graph for this function:
Here is the caller graph for this function:

◆ wrf_arwpost_calc_uvmet()

subroutine mod_realinput_wrfarw::wrf_arwpost_calc_uvmet ( real(rp), dimension(:,:,:), intent(out)  u_latlon,
real(rp), dimension(:,:,:), intent(out)  v_latlon,
real(rp), dimension(:,:,:), intent(in)  u_on_map,
real(rp), dimension(:,:,:), intent(in)  v_on_map,
real(rp), dimension(:,:), intent(in)  xlon,
real(rp), dimension(:,:), intent(in)  xlat,
character(len=*), intent(in)  basename,
integer, intent(in)  K1,
integer, intent(in)  I1,
integer, intent(in)  J1 
)

convert vector varibles from map-projected grid on wrf model to lat-lon grid

Definition at line 872 of file mod_realinput_wrfarw.F90.

References scale_const::const_d2r, scale_const::const_pi, and scale_file::file_open().

Referenced by parentatmosinputwrfarw().

872  use scale_const, only: &
873  d2r => const_d2r, &
874  pi => const_pi
875  use scale_file, only: &
876  file_open, &
877  file_get_attribute
878  implicit none
879  real(RP), intent(out) :: u_latlon(:,:,:)
880  real(RP), intent(out) :: v_latlon(:,:,:)
881  real(RP), intent(in ) :: u_on_map(:,:,:)
882  real(RP), intent(in ) :: v_on_map(:,:,:)
883  real(RP), intent(in ) :: xlon(:,:)
884  real(RP), intent(in ) :: xlat(:,:)
885  integer , intent(in ) :: k1, i1, j1
886 
887  character(len=*), intent( in) :: basename
888 
889  integer :: fid
890 
891  real(RP) :: truelat1, truelat2
892  real(RP) :: stand_lon
893  real(RP) :: diff
894  real(RP) :: alpha
895  real(RP) :: sine(i1,j1)
896  real(RP) :: cose(i1,j1)
897  real(RP) :: cone
898  integer :: map_proj
899 
900  real(RP) :: dum_r(1)
901  integer :: dum_i(1)
902 
903 
904  integer :: k, i, j
905  !---------------------------------------------------------------------------
906 
907  call file_open( basename, fid, rankid=myrank, single=.true., postfix="" )
908 
909  call file_get_attribute( fid, "global", "MAP_PROJ", dum_i(:) )
910  map_proj = dum_i(1)
911 
912  call file_get_attribute( fid, "global", "TRUELAT1", dum_r(:) )
913  truelat1 = dum_r(1) * d2r
914  call file_get_attribute( fid, "global", "TRUELAT2", dum_r(:) )
915  truelat2 = dum_r(1) * d2r
916  call file_get_attribute( fid, "global", "STAND_LON", dum_r(:) )
917  stand_lon = dum_r(1) * d2r
918 
919  ! No need to rotate
920  if ( map_proj .ge. 3 ) then
921  u_latlon(:,:,:) = u_on_map(:,:,:)
922  v_latlon(:,:,:) = v_on_map(:,:,:)
923 
924  return
925  endif
926 
927  ! Lambert Conformal mapping
928  cone = 1.0_rp ! PS
929  if ( map_proj .eq. 1 ) then
930  if ( abs(truelat1-truelat2) .gt. 0.1_rp*d2r ) then
931  cone = ( log(cos(truelat1)) - &
932  log(cos(truelat2)) ) / &
933  ( log(tan((pi*0.5_rp-abs(truelat1))*0.5_rp )) - &
934  log(tan((pi*0.5_rp-abs(truelat2))*0.5_rp )) )
935  else
936  cone = sin( abs(truelat1) )
937  endif
938  endif
939 
940  do j = 1, j1
941  do i = 1, i1
942  diff = xlon(i,j) - stand_lon
943  if ( diff .gt. pi ) then
944  diff = diff - pi*2.0_rp
945  endif
946  if ( diff .lt. -pi ) then
947  diff = diff + pi*2.0_rp
948  endif
949  alpha = diff * cone * sign(1.0_rp, xlat(i,j))
950  sine(i,j) = sin( alpha )
951  cose(i,j) = cos( alpha )
952  enddo
953  enddo
954 
955  do j = 1, j1
956  do i = 1, i1
957  do k = 1, k1
958  u_latlon(k,i,j) = v_on_map(k,i,j)*sine(i,j) + u_on_map(k,i,j)*cose(i,j)
959  v_latlon(k,i,j) = v_on_map(k,i,j)*cose(i,j) - u_on_map(k,i,j)*sine(i,j)
960  enddo
961  enddo
962  enddo
963 
964  return
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
module file
Definition: scale_file.F90:15
module CONSTANT
Definition: scale_const.F90:11
real(rp), public const_pi
pi
Definition: scale_const.F90:31
Here is the call graph for this function:
Here is the caller graph for this function: