SCALE-RM
mod_realinput_wrfarw.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_stdio
20  use scale_tracer
21  use scale_process, only: &
22  myrank => prc_myrank, &
24  use scale_external_io, only: &
25  iwrfarw
26 
27  !-----------------------------------------------------------------------------
28  implicit none
29  private
30  !-----------------------------------------------------------------------------
31  !
32  !++ Public procedure
33  !
34  public :: parentatomsetupwrfarw
35  public :: parentatomopenwrfarw
36  public :: parentatominputwrfarw
37  public :: parentlandsetupwrfarw
38  public :: parentlandinputwrfarw
39  public :: parentoceansetupwrfarw
40  public :: parentoceanopenwrfarw
41  public :: parentoceaninputwrfarw
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Public parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private procedure
50  !
51  !-----------------------------------------------------------------------------
52  !
53  !++ Private parameters & variables
54  !
55  ! Defined parameters in WRF
56  integer, parameter :: mdlid = iwrfarw
57 
58  real(RP), parameter :: t0 = 300.0_rp
59  real(RP), parameter :: p0 = 1000.0e+2_rp
60  real(RP), parameter :: rd = 287.04_rp
61  real(RP), parameter :: cp = 7.0_rp * rd / 2.0_rp
62  real(RP), parameter :: rcp = rd / cp
63 
64  integer, parameter :: cosin = 1
65  integer, parameter :: sine = 2
66 
67  real(SP), allocatable :: read_xy (:,:,:)
68  real(SP), allocatable :: read_uy (:,:,:)
69  real(SP), allocatable :: read_xv (:,:,:)
70  real(SP), allocatable :: read_zxy(:,:,:,:)
71  real(SP), allocatable :: read_wxy(:,:,:,:)
72  real(SP), allocatable :: read_zuy(:,:,:,:)
73  real(SP), allocatable :: read_zxv(:,:,:,:)
74  real(SP), allocatable :: read_lzxy(:,:,:,:)
75  real(SP), allocatable :: read_lz(:,:)
76 
77  real(RP), allocatable :: p_org (:,:,:)
78  real(RP), allocatable :: pb_org (:,:,:)
79  real(RP), allocatable :: ph_org (:,:,:)
80  real(RP), allocatable :: phb_org (:,:,:)
81 
82  logical, private :: wrfout = .false. ! file type switch (wrfout or wrfrst)
83 
84  !-----------------------------------------------------------------------------
85 contains
86  !-----------------------------------------------------------------------------
88  subroutine parentatomsetupwrfarw( &
89  dims, &
90  timelen, &
91  basename_org )
92  use scale_external_io, only: &
93  iwrfarw, &
95  implicit none
96 
97  integer, intent(out) :: dims(6)
98  integer, intent(out) :: timelen
99  character(len=H_LONG), intent(in) :: basename_org
100 
101  logical :: WRF_FILE_TYPE = .false. ! wrf filetype: T=wrfout, F=wrfrst
102 
103  namelist / param_mkinit_real_wrfarw / &
104  wrf_file_type
105 
106  integer :: dims_wrf(7)
107  integer :: ierr
108  !---------------------------------------------------------------------------
109 
110  if( io_l ) write(io_fid_log,*) '+++ Real Case/Atom Input File Type: WRF-ARW'
111 
112  !--- read namelist
113  rewind(io_fid_conf)
114  read(io_fid_conf,nml=param_mkinit_real_wrfarw,iostat=ierr)
115  if( ierr > 0 ) then
116  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_WRFARW. Check!'
117  call prc_mpistop
118  endif
119  if( io_lnml ) write(io_fid_log,nml=param_mkinit_real_wrfarw)
120 
121  call externalfilegetshape( dims_wrf, timelen, mdlid, basename_org, myrank, single=.true. )
122  dims(1:6) = dims_wrf(1:6)
123 
124  if ( wrf_file_type ) then
125  wrfout = .true.
126  if( io_l ) write(io_fid_log,*) '+++ WRF-ARW FILE-TYPE: WRF History Output'
127  else
128  wrfout = .false.
129  if( io_l ) write(io_fid_log,*) '+++ WRF-ARW FILE-TYPE: WRF Restart'
130  endif
131 
132 
133  allocate( read_xy( dims(2),dims(3),1) )
134  allocate( read_uy( dims(5),dims(3),1) )
135  allocate( read_xv( dims(2),dims(6),1) )
136  allocate( read_zxy(dims(1),dims(2),dims(3),1) )
137  allocate( read_wxy(dims(4),dims(2),dims(3),1) )
138  allocate( read_zuy(dims(1),dims(5),dims(3),1) )
139  allocate( read_zxv(dims(1),dims(2),dims(6),1) )
140 
141  allocate( p_org(dims(1),dims(2),dims(3)) )
142  allocate( pb_org(dims(1),dims(2),dims(3)) )
143  allocate( ph_org(dims(4),dims(2),dims(3)) )
144  allocate( phb_org(dims(4),dims(2),dims(3)) )
145 
146  return
147  end subroutine parentatomsetupwrfarw
148 
149  !-----------------------------------------------------------------------------
150  subroutine parentatomopenwrfarw
151  implicit none
152 
153  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[AtomOpenWRFARW]'
154  return
155  end subroutine parentatomopenwrfarw
156 
157  !-----------------------------------------------------------------------------
158  subroutine parentatominputwrfarw( &
159  velz_org, &
160  llvelx_org, &
161  llvely_org, &
162  pres_org, &
163  temp_org, &
164  qtrc_org, &
165  lon_org, &
166  lat_org, &
167  cz_org, &
168  basename, &
169  mptype_parent, &
170  dims, &
171  it ) ! (in)
172  use scale_const, only: &
173  d2r => const_d2r, &
174  laps => const_laps, &
175  grav => const_grav
176  use scale_external_io, only: &
177  externalfileread
178  use scale_atmos_thermodyn, only: &
179  thermodyn_pott => atmos_thermodyn_pott
180  implicit none
181  real(RP), intent(out) :: velz_org(:,:,:)
182  real(RP), intent(out) :: llvelx_org(:,:,:)
183  real(RP), intent(out) :: llvely_org(:,:,:)
184  real(RP), intent(out) :: pres_org(:,:,:)
185  real(RP), intent(out) :: temp_org(:,:,:)
186  real(RP), intent(out) :: qtrc_org(:,:,:,:)
187  real(RP), intent(out) :: lon_org(:,:)
188  real(RP), intent(out) :: lat_org(:,:)
189  real(RP), intent(out) :: cz_org(:,:,:)
190  character(len=*), intent(in) :: basename
191  integer, intent(in) :: mptype_parent
192  integer, intent(in) :: dims(6)
193  integer, intent(in) :: it
194 
195  ! full level
196  real(RP) :: velx_org(dims(1)+2,dims(2),dims(3))
197  real(RP) :: vely_org(dims(1)+2,dims(2),dims(3))
198  real(RP) :: pott_org(dims(1)+2,dims(2),dims(3))
199  real(RP) :: topo_org( dims(2),dims(3))
200 
201  ! half level
202  real(RP) :: velzs_org(dims(4),dims(2),dims(3))
203  real(RP) :: velxs_org(dims(1),dims(5),dims(3))
204  real(RP) :: velys_org(dims(1),dims(2),dims(6))
205  real(RP) :: geof_org (dims(4),dims(2),dims(3))
206 
207  real(RP) :: qhyd
208 
209  integer :: k, i, j, iq
210 
211  character(len=H_MID) :: varname_T
212  character(len=H_MID) :: varname_W
213  character(len=H_MID) :: varname_U
214  character(len=H_MID) :: varname_V
215 
216  logical :: lack_of_val
217 
218  integer :: ierr
219  !---------------------------------------------------------------------------
220 
221  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[InputWRF]'
222 
223  if ( wrfout ) then
224  varname_t = "T"
225  varname_w = "W"
226  varname_u = "U"
227  varname_v = "V"
228  else
229  varname_t = "T_1"
230  varname_w = "W_1"
231  varname_u = "U_1"
232  varname_v = "V_1"
233  endif
234 
235 
236  call externalfileread( read_xy(:,:,:), basename, "XLAT", it, it, myrank, mdlid, single=.true. )
237  lat_org(:,:) = real( read_xy(:,:,1), kind=RP ) * D2R
238 
239  call externalfileread( read_xy(:,:,:), basename, "XLONG", it, it, myrank, mdlid, single=.true. )
240  lon_org(:,:) = real( read_xy(:,:,1), kind=RP ) * D2R
241 
242  call externalfileread( read_xy(:,:,:), basename, "HGT", it, it, myrank, mdlid, single=.true. )
243  topo_org(:,:) = real( read_xy(:,:,1), kind=rp )
244 
245  call externalfileread( read_wxy(:,:,:,:), basename, "PH", it, it, myrank, mdlid, single=.true., zstag=.true. )
246  ph_org(:,:,:) = real( read_wxy(:,:,:,1), kind=rp )
247 
248  call externalfileread( read_wxy(:,:,:,:), basename, "PHB", it, it, myrank, mdlid, single=.true., zstag=.true. )
249  phb_org(:,:,:) = real( read_wxy(:,:,:,1), kind=rp )
250 
251  call externalfileread( read_zxy(:,:,:,:), basename, "P", it, it, myrank, mdlid, single=.true. )
252  p_org(:,:,:) = real( read_zxy(:,:,:,1), kind=rp )
253 
254  call externalfileread( read_zxy(:,:,:,:), basename, "PB", it, it, myrank, mdlid, single=.true. )
255  pb_org(:,:,:) = real( read_zxy(:,:,:,1), kind=rp )
256 
257  call externalfileread( read_wxy(:,:,:,:), basename, varname_w, it, it, myrank, mdlid, single=.true., zstag=.true. )
258  velzs_org(:,:,:) = real( read_wxy(:,:,:,1), kind=rp )
259 
260  call externalfileread( read_zuy(:,:,:,:), basename, varname_u, it, it, myrank, mdlid, single=.true., xstag=.true. )
261 
262  velxs_org(:,:,:) = real( read_zuy(:,:,:,1), kind=rp )
263 
264  call externalfileread( read_zxv(:,:,:,:), basename, varname_v, it, it, myrank, mdlid, single=.true., ystag=.true. )
265  velys_org(:,:,:) = real( read_zxv(:,:,:,1), kind=rp )
266 
267  ! from half level to full level
268  do j = 1, dims(3)
269  do i = 1, dims(2)
270  do k = 1, dims(1)
271  velz_org(k+2,i,j) = ( velzs_org(k,i,j) + velzs_org(k+1,i,j) ) * 0.5_rp
272  velx_org(k+2,i,j) = ( velxs_org(k,i,j) + velxs_org(k,i+1,j) ) * 0.5_rp
273  vely_org(k+2,i,j) = ( velys_org(k,i,j) + velys_org(k,i,j+1) ) * 0.5_rp
274  end do
275  velz_org(1:2,i,j) = 0.0_rp
276  velx_org(1:2,i,j) = 0.0_rp
277  vely_org(1:2,i,j) = 0.0_rp
278  end do
279  end do
280 
281  call wrf_arwpost_calc_uvmet( llvelx_org, llvely_org, & ! (out)
282  velx_org, vely_org, & ! (in)
283  lon_org, lat_org, & ! (in)
284  basename, & ! (in)
285  dims(1)+2, dims(2), dims(3) ) ! (in)
286 
287  qtrc_org(:,:,:,:) = 0.0_rp
288  call externalfileread( read_xy(:,:,:), basename, "Q2", it, it, myrank, mdlid, single=.true. )
289  qtrc_org(2,:,:,i_qv) = real(read_xy(:,:,1),kind=rp)
290  qtrc_org(1,:,:,i_qv) = qtrc_org(2,:,:,i_qv)
291 
292  call externalfileread( read_zxy(:,:,:,:), basename, "QVAPOR", it, it, myrank, mdlid, single=.true. )
293  qtrc_org(3:,:,:,i_qv) = real(read_zxy(:,:,:,1),kind=rp)
294 
295 #ifndef DRY
296  if( mptype_parent > 0 ) then
297  call externalfileread( read_zxy(:,:,:,:), basename, "QCLOUD", it, it, myrank, mdlid, single=.true. )
298  qtrc_org(3:,:,:,i_qc) = real( read_zxy(:,:,:,1), kind=rp)
299 
300  call externalfileread( read_zxy(:,:,:,:), basename, "QRAIN", it, it, myrank, mdlid, single=.true. )
301  qtrc_org(3:,:,:,i_qr) = real( read_zxy(:,:,:,1), kind=rp)
302  endif
303 
304  if( mptype_parent > 3 ) then
305  call externalfileread( read_zxy(:,:,:,:), basename, "QICE", it, it, myrank, mdlid, single=.true. )
306  qtrc_org(3:,:,:,i_qi) = real( read_zxy(:,:,:,1), kind=rp)
307 
308  call externalfileread( read_zxy(:,:,:,:), basename, "QSNOW", it, it, myrank, mdlid, single=.true. )
309  qtrc_org(3:,:,:,i_qs) = real( read_zxy(:,:,:,1), kind=rp)
310  endif
311 
312  if( mptype_parent > 5 ) then
313  call externalfileread( read_zxy(:,:,:,:), basename, "QGRAUP", it, it, myrank, mdlid, single=.true. )
314  qtrc_org(3:,:,:,i_qg) = real( read_zxy(:,:,:,1), kind=rp)
315  endif
316 
317  ! convert mixing ratio to specific ratio
318  do j = 1, dims(3)
319  do i = 1, dims(2)
320  do k = 1, dims(1)+2
321  qhyd = 0.0_rp
322  do iq = 1, min( mptype_parent, 6)
323  qhyd = qhyd + qtrc_org(k,i,j,iq)
324  end do
325  do iq = 1, min( mptype_parent, 6)
326  qtrc_org(k,i,j,iq) = qtrc_org(k,i,j,iq) / ( 1.0_rp + qhyd )
327  end do
328  end do
329  end do
330  end do
331 
332  if( mptype_parent > 6 ) then
333  call externalfileread( read_zxy(:,:,:,:), basename, "NC", it, it, myrank, mdlid, single=.true. )
334  qtrc_org(3:,:,:,i_nc) = real( read_zxy(:,:,:,1), kind=rp )
335 
336  call externalfileread( read_zxy(:,:,:,:), basename, "NR", it, it, myrank, mdlid, single=.true. )
337  qtrc_org(3:,:,:,i_nr) = real( read_zxy(:,:,:,1), kind=rp )
338 
339  call externalfileread( read_zxy(:,:,:,:), basename, "NI", it, it, myrank, mdlid, single=.true. )
340  qtrc_org(3:,:,:,i_ni) = real( read_zxy(:,:,:,1), kind=rp )
341 
342  call externalfileread( read_zxy(:,:,:,:), basename, "NS", it, it, myrank, mdlid, single=.true. )
343  qtrc_org(3:,:,:,i_ns) = real( read_zxy(:,:,:,1), kind=rp )
344 
345  call externalfileread( read_zxy(:,:,:,:), basename, "NG", it, it, myrank, mdlid, single=.true. )
346  qtrc_org(3:,:,:,i_ng) = real( read_zxy(:,:,:,1), kind=rp )
347  endif
348 #endif
349 
350  do iq = 1, qa
351  do j = 1, dims(3)
352  do i = 1, dims(2)
353  do k = 1, dims(1)+2
354  qtrc_org(k,i,j,iq) = max( qtrc_org(k,i,j,iq), 0.0_rp )
355  end do
356  end do
357  end do
358  end do
359 
360 
361  call externalfileread( read_zxy(:,:,:,:), basename, varname_t, it, it, myrank, mdlid, single=.true. )
362  pott_org(3:,:,:) = real( read_zxy(:,:,:,1), kind=RP ) + t0
363  call externalfileread( read_xy(:,:,:), basename, "T2", it, it, myrank, mdlid, single=.true. )
364  temp_org(2,:,:) = real( read_xy(:,:,1), kind=rp )
365 
366  call externalfileread( read_xy(:,:,:), basename, "PSFC", it, it, myrank, mdlid, single=.true. )
367  pres_org(2,:,:) = real( read_xy(:,:,1), kind=rp )
368 
369  do j = 1, dims(3)
370  do i = 1, dims(2)
371  do k = 3, dims(1)+2
372  pres_org(k,i,j) = p_org(k-2,i,j) + pb_org(k-2,i,j)
373  temp_org(k,i,j) = pott_org(k,i,j) * ( pres_org(k,i,j) / p0 )**rcp
374  end do
375  pott_org(2,i,j) = temp_org(2,i,j) * ( p0/pres_org(2,i,j) )**rcp
376  temp_org(1,i,j) = temp_org(2,i,j) + laps * topo_org(i,j)
377  pres_org(1,i,j) = p0 * ( temp_org(1,i,j) / pott_org(2,i,j) )**(1.0_rp/rcp) ! pott_org(1,i,j) == pott_org(2,i,j)
378  end do
379  end do
380 
381  do j = 1, dims(3)
382  do i = 1, dims(2)
383  ! convert to geopotential height to use as real height in WRF
384  do k = 1, dims(4)
385  geof_org(k,i,j) = ( ph_org(k,i,j) + phb_org(k,i,j) ) / grav
386  end do
387  ! make half level of geopotential height from face level
388  do k = 1, dims(1)
389  cz_org(k+2,i,j) = ( geof_org(k,i,j) + geof_org(k+1,i,j) ) * 0.5_rp
390  end do
391  cz_org(2,i,j) = topo_org(i,j)
392  cz_org(1,i,j) = 0.0_rp
393  end do
394  end do
395 
396 
397 #ifdef DEBUG
398  !k=1 ; i=int(dims(2)/2) ; j=int(dims(3)/2) ; iq = 1
399  k=2 ; i=3 ; j=3 ; iq = 1
400  write(*,*) "read 3D wrf data",i,j,k
401  write(*,*) "lon_org ",lon_org(i,j)/d2r
402  write(*,*) "lat_org ",lat_org(i,j)/d2r
403  write(*,*) "cz_org ",cz_org(k,i,j)
404  write(*,*) "pres_org ",pres_org(k,i,j)
405  write(*,*) "velx_org ",llvelx_org(k,i,j)
406  write(*,*) "vely_org ",llvely_org(k,i,j)
407  write(*,*) "velz_org ",velz_org(k,i,j)
408  write(*,*) "temp_org ",temp_org(k,i,j)
409  write(*,*) "qtrc_org ",qtrc_org(k,i,j,iq)
410  k=3 ; i=3 ; j=3 ; iq = 1
411  write(*,*) "read 3D wrf data",i,j,k
412  write(*,*) "lon_org ",lon_org(i,j)/d2r
413  write(*,*) "lat_org ",lat_org(i,j)/d2r
414  write(*,*) "cz_org ",cz_org(k,i,j)
415  write(*,*) "pres_org ",pres_org(k,i,j)
416  write(*,*) "velx_org ",llvelx_org(k,i,j)
417  write(*,*) "vely_org ",llvely_org(k,i,j)
418  write(*,*) "velz_org ",velz_org(k,i,j)
419  write(*,*) "temp_org ",temp_org(k,i,j)
420  write(*,*) "qtrc_org ",qtrc_org(k,i,j,iq)
421 #endif
422 
423  return
424  end subroutine parentatominputwrfarw
425 
426  !-----------------------------------------------------------------------------
428  subroutine parentlandsetupwrfarw( &
429  ldims, &
430  basename_land )
431  use scale_external_io, only: &
432  iwrfarw, &
434  implicit none
435 
436  integer, intent(out) :: ldims(3)
437  character(len=H_LONG), intent(in) :: basename_land
438 
439  logical :: WRF_FILE_TYPE = .false. ! wrf filetype: T=wrfout, F=wrfrst
440 
441  namelist / param_mkinit_real_wrfarw / &
442  wrf_file_type
443 
444  integer :: dims_wrf(7)
445  integer :: timelen
446  integer :: ierr
447  !---------------------------------------------------------------------------
448 
449  if( io_l ) write(io_fid_log,*) '+++ Real Case/Atom Input File Type: WRF-ARW'
450 
451  !--- read namelist
452  rewind(io_fid_conf)
453  read(io_fid_conf,nml=param_mkinit_real_wrfarw,iostat=ierr)
454  if( ierr > 0 ) then
455  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_WRFARW. Check!'
456  call prc_mpistop
457  endif
458  if( io_lnml ) write(io_fid_log,nml=param_mkinit_real_wrfarw)
459 
460 
461  call externalfilegetshape( dims_wrf, timelen, mdlid, basename_land, myrank, single=.true. )
462  ldims(1) = dims_wrf(7)
463  ldims(2) = dims_wrf(2)
464  ldims(3) = dims_wrf(3)
465 
466  if ( wrf_file_type ) then
467  wrfout = .true.
468  if( io_l ) write(io_fid_log,*) '+++ WRF-ARW FILE-TYPE: WRF History Output'
469  else
470  wrfout = .false.
471  if( io_l ) write(io_fid_log,*) '+++ WRF-ARW FILE-TYPE: WRF Restart'
472  endif
473 
474 
475  if ( .not. allocated(read_xy) ) then
476  allocate( read_xy( ldims(2),ldims(3),1) )
477  end if
478 
479  allocate( read_lzxy(ldims(1),ldims(2),ldims(3),1) )
480  allocate( read_lz(ldims(1),1) )
481 
482  return
483  end subroutine parentlandsetupwrfarw
484 
485  !-----------------------------------------------------------------------------
486  subroutine parentlandinputwrfarw( &
487  tg_org, &
488  sh2o_org, &
489  lst_org, &
490  ust_org, &
491  albg_org, &
492  topo_org, &
493  lmask_org, &
494  llon_org, &
495  llat_org, &
496  lz_org, &
497  basename, &
498  ldims, &
499  use_file_landwater, &
500  it )
501  use scale_const, only: &
502  d2r => const_d2r, &
503  undef => const_undef, &
504  i_lw => const_i_lw, &
505  i_sw => const_i_sw
506  use scale_external_io, only: &
507  externalfileread, &
509  implicit none
510  real(RP), intent(out) :: tg_org(:,:,:)
511  real(RP), intent(out) :: sh2o_org(:,:,:)
512  real(RP), intent(out) :: lst_org(:,:)
513  real(RP), intent(out) :: ust_org(:,:)
514  real(RP), intent(out) :: albg_org(:,:,:)
515  real(RP), intent(out) :: topo_org(:,:)
516  real(RP), intent(out) :: lmask_org(:,:)
517  real(RP), intent(out) :: llon_org(:,:)
518  real(RP), intent(out) :: llat_org(:,:)
519  real(RP), intent(out) :: lz_org(:)
520  character(len=*), intent( in) :: basename
521  integer, intent( in) :: ldims(3)
522  logical, intent( in) :: use_file_landwater ! use land water data from files
523  integer, intent( in) :: it
524 
525 
526  integer :: k, i, j
527 
528  logical :: existence
529 
530  !---------------------------------------------------------------------------
531 
532  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[LandInputWRF]'
533 
534  call externalfileread( read_xy(:,:,:), basename, "XLAT", it, it, myrank, mdlid, single=.true. )
535  llat_org(:,:) = read_xy(:,:,1) * d2r
536 
537  call externalfileread( read_xy(:,:,:), basename, "XLONG", it, it, myrank, mdlid, single=.true. )
538  llon_org(:,:) = read_xy(:,:,1) * d2r
539 
540  call externalfileread( read_xy(:,:,:), basename, "HGT", it, it, myrank, mdlid, single=.true. )
541  topo_org(:,:) = read_xy(:,:,1)
542 
543  ! depth
544  call externalfileread( read_lz(:,:), &
545  basename, "ZS", it, 1, myrank, mdlid, ldims(1), single=.true. )
546  lz_org(:) = read_lz(:,1)
547 
548  ! land mask (1:land, 0:water)
549  call externalfileread( read_xy(:,:,:), &
550  basename, "LANDMASK", it, 1, myrank, mdlid, single=.true. )
551  lmask_org(:,:) = read_xy(:,:,1)
552 
553  ! soil temperature [K]
554  call externalfileread( read_lzxy(:,:,:,:), &
555  basename, "TSLB", it, 1, myrank, mdlid, single=.true., landgrid=.true. )
556  tg_org(:,:,:) = read_lzxy(:,:,:,1)
557 
558  ! soil liquid water [m3 m-3] (no wrfout-default)
559  if( use_file_landwater ) then
560  call externalfilevarexistence( existence, basename, "SH2O", myrank, mdlid, single=.true. )
561  if ( existence ) then
562  call externalfileread( read_lzxy(:,:,:,:), &
563  basename, "SH2O", it, 1, myrank, mdlid, single=.true., landgrid=.true. )
564  sh2o_org(:,:,:) = read_lzxy(:,:,:,1)
565  else
566  sh2o_org(:,:,:) = undef
567  endif
568  endif
569 
570 ! ! surface runoff [mm]
571 ! call ExternalFileRead( read_xy(:,:,:), &
572 ! BASENAME, "SFROFF", it, 1, myrank, mdlid, single=.true. )
573 ! do j = 1, dims(9)
574 ! do i = 1, dims(8)
575 ! org_3D(i,j) = read_xy(i,j,1) * 1000.0_DP * dwatr
576 ! enddo
577 ! enddo
578 
579 
580  ! SURFACE SKIN TEMPERATURE [K]
581  call externalfileread( read_xy(:,:,:), &
582  basename, "TSK", it, 1, myrank, mdlid, single=.true. )
583  lst_org(:,:) = read_xy(:,:,1)
584 
585  ust_org(:,:) = lst_org(:,:)
586 
587  ! ALBEDO [-]
588  call externalfileread( read_xy(:,:,:), &
589  basename, "ALBEDO", it, 1, myrank, mdlid, single=.true. )
590  albg_org(:,:,i_sw) = read_xy(:,:,1)
591 
592  ! SURFACE EMISSIVITY [-]
593  call externalfileread( read_xy(:,:,:), &
594  basename, "EMISS", it, 1, myrank, mdlid, single=.true. )
595  do j = 1, ldims(3)
596  do i = 1, ldims(2)
597  albg_org(i,j,i_lw) = 1.0_dp - read_xy(i,j,1)
598  enddo
599  enddo
600 
601 
602 ! ! SNOW WATER EQUIVALENT [kg m-2] (no wrfout-default)
603 ! call ExternalFileVarExistence( existence, BASENAME, "SNOW", myrank, mdlid, single=.true. )
604 ! if ( existence ) then
605 ! call ExternalFileRead( read_xy(:,:,:), &
606 ! BASENAME, "SNOW", it, 1, myrank, mdlid, single=.true. )
607 ! swowq_org(:,:,:) = read_xy(:,:,1)
608 ! else
609 ! snowq_org(:,:) = UNDEF
610 ! endif
611 
612 ! ! AVERAGE SNOW TEMPERATURE [C] (no wrfout-default)
613 ! call ExternalFileVarExistence( existence, BASENAME, "TSNAV", myrank, mdlid, single=.true. )
614 ! if ( existence ) then
615 ! call ExternalFileRead( read_xy(:,:,:), &
616 ! BASENAME, "TSNAV", it, 1, myrank, mdlid, single=.true. )
617 ! do j = 1, dims(9)
618 ! do i = 1, dims(8)
619 ! snowt_org(i,j,n) = read_xy(i,j,1) + TEM00
620 ! enddo
621 ! enddo
622 ! else
623 ! snowt_org(:,:) = UNDEFF
624 ! endif
625 
626  return
627  end subroutine parentlandinputwrfarw
628 
629  !-----------------------------------------------------------------------------
631  subroutine parentoceansetupwrfarw( &
632  odims, &
633  timelen, &
634  basename_org )
635  use scale_external_io, only: &
636  iwrfarw, &
638  implicit none
639 
640  integer, intent(out) :: odims(2)
641  integer, intent(out) :: timelen
642  character(len=H_LONG), intent(in) :: basename_org
643 
644  logical :: WRF_FILE_TYPE = .false. ! wrf filetype: T=wrfout, F=wrfrst
645 
646  namelist / param_mkinit_real_wrfarw / &
647  wrf_file_type
648 
649  integer :: dims_wrf(7)
650  integer :: ierr
651  !---------------------------------------------------------------------------
652 
653  if( io_l ) write(io_fid_log,*) '+++ Real Case/Ocean Input File Type: WRF-ARW'
654 
655  !--- read namelist
656  rewind(io_fid_conf)
657  read(io_fid_conf,nml=param_mkinit_real_wrfarw,iostat=ierr)
658  if( ierr > 0 ) then
659  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_REAL_WRFARW. Check!'
660  call prc_mpistop
661  endif
662  if( io_lnml ) write(io_fid_log,nml=param_mkinit_real_wrfarw)
663 
664 
665  call externalfilegetshape( dims_wrf, timelen, mdlid, basename_org, myrank, single=.true. )
666  odims(1) = dims_wrf(2)
667  odims(2) = dims_wrf(3)
668 
669  if ( wrf_file_type ) then
670  wrfout = .true.
671  if( io_l ) write(io_fid_log,*) '+++ WRF-ARW FILE-TYPE: WRF History Output'
672  else
673  wrfout = .false.
674  if( io_l ) write(io_fid_log,*) '+++ WRF-ARW FILE-TYPE: WRF Restart'
675  endif
676 
677 
678  if ( .not. allocated(read_xy) ) then
679  allocate( read_xy( odims(1),odims(2),1) )
680  end if
681 
682  return
683  end subroutine parentoceansetupwrfarw
684 
685  !-----------------------------------------------------------------------------
686  subroutine parentoceanopenwrfarw
687  implicit none
688 
689  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanOpenWRFARW]'
690  return
691  end subroutine parentoceanopenwrfarw
692 
693  !-----------------------------------------------------------------------------
694  subroutine parentoceaninputwrfarw( &
695  tw_org, &
696  sst_org, &
697  albw_org, &
698  z0w_org, &
699  omask_org, &
700  olon_org, &
701  olat_org, &
702  basename, &
703  odims, &
704  it )
705  use scale_const, only: &
706  d2r => const_d2r, &
707  undef => const_undef, &
708  i_lw => const_i_lw, &
709  i_sw => const_i_sw
710  use scale_external_io, only: &
711  externalfileread, &
713  implicit none
714  real(RP), intent(out) :: tw_org(:,:)
715  real(RP), intent(out) :: sst_org(:,:)
716  real(RP), intent(out) :: albw_org(:,:,:)
717  real(RP), intent(out) :: z0w_org(:,:)
718  real(RP), intent(out) :: omask_org(:,:)
719  real(RP), intent(out) :: olon_org(:,:)
720  real(RP), intent(out) :: olat_org(:,:)
721  character(len=*), intent( in) :: basename
722  integer, intent( in) :: odims(2)
723  integer, intent( in) :: it
724 
725 
726  integer :: i, j
727 
728  logical :: existence
729 
730  !---------------------------------------------------------------------------
731 
732  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanInputWRF]'
733 
734  call externalfileread( read_xy(:,:,:), basename, "XLAT", it, it, myrank, mdlid, single=.true. )
735  olat_org(:,:) = read_xy(:,:,1) * d2r
736 
737  call externalfileread( read_xy(:,:,:), basename, "XLONG", it, it, myrank, mdlid, single=.true. )
738  olon_org(:,:) = read_xy(:,:,1) * d2r
739 
740  ! land mask (1:land, 0:water)
741  call externalfileread( read_xy(:,:,:), &
742  basename, "LANDMASK", it, 1, myrank, mdlid, single=.true. )
743  omask_org(:,:) = read_xy(:,:,1)
744 
745  ! SEA SURFACE TEMPERATURE [K]
746  call externalfileread( read_xy(:,:,:), &
747  basename, "SST", it, 1, myrank, mdlid, single=.true. )
748  sst_org(:,:) = read_xy(:,:,1)
749 
750  tw_org(:,:) = sst_org(:,:)
751 
752  ! ALBEDO [-]
753  call externalfileread( read_xy(:,:,:), &
754  basename, "ALBEDO", it, 1, myrank, mdlid, single=.true. )
755  albw_org(:,:,i_sw) = read_xy(:,:,1)
756 
757  ! SURFACE EMISSIVITY [-]
758  call externalfileread( read_xy(:,:,:), &
759  basename, "EMISS", it, 1, myrank, mdlid, single=.true. )
760  do j = 1, odims(2)
761  do i = 1, odims(1)
762  albw_org(i,j,i_lw) = 1.0_dp - read_xy(i,j,1)
763  enddo
764  enddo
765 
766  ! TIME-VARYING ROUGHNESS LENGTH [m] (no wrfout-default)
767  call externalfilevarexistence( existence, basename, "ZNT", myrank, mdlid, single=.true. )
768  if ( existence ) then
769  call externalfileread( read_xy(:,:,:), &
770  basename, "ZNT", it, 1, myrank, mdlid, single=.true. )
771  z0w_org(:,:) = read_xy(:,:,1)
772  else
773  z0w_org(:,:) = undef
774  endif
775 
776  return
777  end subroutine parentoceaninputwrfarw
778 
779  !-----------------------------------------------------------------------------
781  !-----------------------------------------------------------------------------
782  subroutine wrf_arwpost_calc_uvmet( &
783  u_latlon, & ! (out)
784  v_latlon, & ! (out)
785  u_on_map, & ! (in)
786  v_on_map, & ! (in)
787  xlon, & ! (in)
788  xlat, & ! (in)
789  basename , & ! (in)
790  k1, i1, j1 ) ! (in)
791  use scale_const, only: &
792  d2r => const_d2r, &
793  pi => const_pi
794  use scale_external_io, only: &
795  externalfilegetglobalattv
796  implicit none
797  real(RP), intent(out) :: u_latlon(:,:,:)
798  real(RP), intent(out) :: v_latlon(:,:,:)
799  real(RP), intent(in ) :: u_on_map(:,:,:)
800  real(RP), intent(in ) :: v_on_map(:,:,:)
801  real(RP), intent(in ) :: xlon(:,:)
802  real(RP), intent(in ) :: xlat(:,:)
803  integer , intent(in ) :: K1, I1, J1
804 
805  character(len=*), intent( in) :: basename
806 
807  real(RP) :: truelat1, truelat2
808  real(RP) :: stand_lon
809  real(RP) :: diff
810  real(RP) :: alpha
811  real(RP) :: sine(i1,j1)
812  real(RP) :: cose(i1,j1)
813  real(RP) :: cone
814  integer :: map_proj
815 
816  real(RP) :: dum_r(1)
817  integer :: dum_i(1)
818 
819 
820  integer :: k, i, j
821  !---------------------------------------------------------------------------
822 
823  call externalfilegetglobalattv( dum_i, iwrfarw, basename, "MAP_PROJ", myrank, single=.true. )
824  map_proj = dum_i(1)
825  call externalfilegetglobalattv( dum_r, iwrfarw, basename, "TRUELAT1", myrank, single=.true. )
826  truelat1 = dum_r(1) * d2r
827  call externalfilegetglobalattv( dum_r, iwrfarw, basename, "TRUELAT2", myrank, single=.true. )
828  truelat2 = dum_r(1) * d2r
829  call externalfilegetglobalattv( dum_r, iwrfarw, basename, "STAND_LON", myrank, single=.true. )
830  stand_lon = dum_r(1) * d2r
831 
832  ! No need to rotate
833  if ( map_proj .ge. 3 ) then
834  u_latlon(:,:,:) = u_on_map(:,:,:)
835  v_latlon(:,:,:) = v_on_map(:,:,:)
836 
837  return
838  endif
839 
840  ! Lambert Conformal mapping
841  cone = 1.0_rp ! PS
842  if ( map_proj .eq. 1 ) then
843  if ( abs(truelat1-truelat2) .gt. 0.1_rp*d2r ) then
844  cone = ( log(cos(truelat1)) - &
845  log(cos(truelat2)) ) / &
846  ( log(tan((pi*0.5_rp-abs(truelat1))*0.5_rp )) - &
847  log(tan((pi*0.5_rp-abs(truelat2))*0.5_rp )) )
848  else
849  cone = sin( abs(truelat1) )
850  endif
851  endif
852 
853  do j = 1, j1
854  do i = 1, i1
855  diff = xlon(i,j) - stand_lon
856  if ( diff .gt. pi ) then
857  diff = diff - pi*2.0_rp
858  endif
859  if ( diff .lt. -pi ) then
860  diff = diff + pi*2.0_rp
861  endif
862  alpha = diff * cone * sign(1.0_rp, xlat(i,j))
863  sine(i,j) = sin( alpha )
864  cose(i,j) = cos( alpha )
865  enddo
866  enddo
867 
868  do j = 1, j1
869  do i = 1, i1
870  do k = 1, k1
871  u_latlon(k,i,j) = v_on_map(k,i,j)*sine(i,j) + u_on_map(k,i,j)*cose(i,j)
872  v_latlon(k,i,j) = v_on_map(k,i,j)*cose(i,j) - u_on_map(k,i,j)*sine(i,j)
873  enddo
874  enddo
875  enddo
876 
877  return
878  end subroutine wrf_arwpost_calc_uvmet
879 
880 end module mod_realinput_wrfarw
881 !-------------------------------------------------------------------------------
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
subroutine, public parentatomopenwrfarw
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:98
subroutine, public externalfilevarexistence(existence, basename, varname, myrank, mdlid, single)
Check Existence of a Variable.
subroutine, public prc_mpistop
Abort MPI.
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)
module REAL input WRF-ARW
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine, public parentlandsetupwrfarw(ldims, basename_land)
Land Setup.
integer, public i_ng
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:60
subroutine, public parentoceaninputwrfarw(tw_org, sst_org, albw_org, z0w_org, omask_org, olon_org, olat_org, basename, odims, it)
integer, parameter, public iwrfarw
module FILE I/O (netcdf)
real(rp), public const_undef
Definition: scale_const.F90:43
subroutine, public parentoceansetupwrfarw(odims, timelen, basename_org)
Ocean Setup.
module TRACER
integer, public i_ni
integer, public i_qv
subroutine, public externalfilegetshape(dims, timelen, mdlid, basename, myrank, single)
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
integer, public i_nc
module PROCESS
subroutine, public log(type, message)
Definition: dc_log.f90:133
subroutine, public parentatominputwrfarw(velz_org, llvelx_org, llvely_org, pres_org, temp_org, qtrc_org, lon_org, lat_org, cz_org, basename, mptype_parent, dims, it)
subroutine, public parentatomsetupwrfarw(dims, timelen, basename_org)
Atmos Setup.
module CONSTANT
Definition: scale_const.F90:14
integer, public i_nr
integer, public prc_myrank
process num in local communicator
integer, public i_qs
integer, public i_qi
module ATMOSPHERE / Thermodynamics
integer, public i_ns
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public i_qg
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
module PRECISION
real(rp), public const_pi
pi
Definition: scale_const.F90:34
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
subroutine, public parentoceanopenwrfarw
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public rp
integer, public i_qr
integer, public i_qc