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