SCALE-RM
mod_realinput_scale.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_tracer
20  !-----------------------------------------------------------------------------
21  implicit none
22  private
23  !-----------------------------------------------------------------------------
24  !
25  !++ Public procedure
26  !
27  public :: parentatmossetupscale
28  public :: parentatmosopenscale
29  public :: parentatmosinputscale
30  public :: parentlandsetupscale
31  public :: parentlandinputscale
32  public :: parentoceansetupscale
33  public :: parentoceanopenscale
34  public :: parentoceaninputscale
35 
36  !-----------------------------------------------------------------------------
37  !
38  !++ Public parameters & variables
39  !
40  !-----------------------------------------------------------------------------
41  !
42  !++ Private procedure
43  !
44  !-----------------------------------------------------------------------------
45  !
46  !++ Private parameters & variables
47  !
48  integer, private, parameter :: handle = 1
49 
50  real(RP), private, allocatable :: read2d (:,:)
51  real(RP), private, allocatable :: read3d (:,:,:)
52  real(RP), private, allocatable :: read3dl(:,:,:)
53  real(RP), private, allocatable :: read3do(:,:,:)
54 
55  !-----------------------------------------------------------------------------
56 contains
57  !-----------------------------------------------------------------------------
59  subroutine parentatmossetupscale( &
60  dims )
62  parent_kmax, &
63  parent_imax, &
64  parent_jmax, &
65  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
66  nest_tile_num_y => comm_cartesc_nest_tile_num_y
67  implicit none
68 
69  integer, intent(out) :: dims(6)
70  !---------------------------------------------------------------------------
71 
72  log_info("ParentAtmosSetupSCALE",*) 'Real Case/Atmos Input File Type: SCALE-RM'
73 
74  ! full level
75  dims(1) = parent_kmax(handle)
76  dims(2) = parent_imax(handle) * nest_tile_num_x
77  dims(3) = parent_jmax(handle) * nest_tile_num_y
78  ! half level
79  dims(4) = dims(1)
80  dims(5) = dims(2)
81  dims(6) = dims(3)
82 
83  allocate( read2d( parent_imax(handle),parent_jmax(handle)) )
84  allocate( read3d(dims(1),parent_imax(handle),parent_jmax(handle)) )
85 
86  return
87  end subroutine parentatmossetupscale
88 
89  !-----------------------------------------------------------------------------
90  subroutine parentatmosopenscale( &
91  lon_org, &
92  lat_org, &
93  cz_org, &
94  basename_org, &
95  dims )
96  use scale_const, only: &
97  d2r => const_d2r
98  use scale_comm_cartesc_nest, only: &
99  parent_imax, &
100  parent_jmax, &
101  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
102  nest_tile_id => comm_cartesc_nest_tile_id
103  use scale_file, only: &
104  file_open
105  use scale_file_cartesc, only: &
106  file_cartesc_read
107  implicit none
108 
109  real(RP), intent(out) :: lon_org(:,:)
110  real(RP), intent(out) :: lat_org(:,:)
111  real(RP), intent(out) :: cz_org (:,:,:)
112  character(len=*), intent(in) :: basename_org
113  integer, intent(in) :: dims(6)
114 
115  integer :: rank
116  integer :: xloc, yloc
117  integer :: xs, xe, ys, ye
118 
119  logical :: existed
120 
121  integer :: fid
122  integer :: i
123  !---------------------------------------------------------------------------
124 
125  do i = 1, size( nest_tile_id(:) )
126  ! read data from split files
127  rank = nest_tile_id(i)
128 
129  xloc = mod( i-1, nest_tile_num_x ) + 1
130  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
131 
132  xs = parent_imax(handle) * (xloc-1) + 1
133  xe = parent_imax(handle) * xloc
134  ys = parent_jmax(handle) * (yloc-1) + 1
135  ye = parent_jmax(handle) * yloc
136 
137  call file_open( basename_org, & ! (in)
138  fid, & ! (out)
139  aggregate=.false., rankid=rank ) ! (in)
140 
141  call file_cartesc_read( fid, "lon", read2d(:,:) )
142  lon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
143 
144  call file_cartesc_read( fid, "lat", read2d(:,:) )
145  lat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
146 
147  call file_cartesc_read( fid, "height", read3d(:,:,:) )
148  cz_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
149 
150  call file_cartesc_read( fid, "topo", read2d(:,:), existed=existed )
151  cz_org(2,xs:xe,ys:ye) = read2d(:,:)
152 
153  end do
154 
155  cz_org(1,:,:) = 0.0_rp
156 
157  return
158  end subroutine parentatmosopenscale
159 
160  !-----------------------------------------------------------------------------
161  subroutine parentatmosinputscale( &
162  velz_org, &
163  velx_org, &
164  vely_org, &
165  pres_org, &
166  dens_org, &
167  pott_org, &
168  qv_org, &
169  qtrc_org, &
170  cz_org, &
171  basename_org, &
172  same_mptype, &
173  dims, &
174  it )
175  use scale_const, only: &
176  p00 => const_pre00, &
177  cpdry => const_cpdry, &
178  rdry => const_rdry, &
179  grav => const_grav, &
180  laps => const_laps
181  use scale_comm_cartesc, only: &
182  comm_vars8, &
183  comm_wait
184  use scale_comm_cartesc_nest, only: &
185  parent_imax, &
186  parent_jmax, &
187  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
188  nest_tile_id => comm_cartesc_nest_tile_id
189  use scale_file, only: &
190  file_open
191  use scale_file_cartesc, only: &
192  file_cartesc_read
195  use scale_atmos_thermodyn, only: &
196  thermodyn_specific_heat => atmos_thermodyn_specific_heat, &
197  thermodyn_rhot2temp_pres => atmos_thermodyn_rhot2temp_pres
198  use scale_atmos_hydrostatic, only: &
199  hydrostatic_buildrho_real => atmos_hydrostatic_buildrho_real
200  use scale_atmos_hydrometeor, only: &
201  n_hyd, &
202  hyd_name, &
203  num_name
204  use mod_atmos_phy_mp_driver, only: &
206  use mod_atmos_phy_mp_vars, only: &
207  qa_mp, &
208  qs_mp, &
209  qe_mp
210  implicit none
211 
212  real(RP), intent(out) :: velz_org(:,:,:)
213  real(RP), intent(out) :: velx_org(:,:,:)
214  real(RP), intent(out) :: vely_org(:,:,:)
215  real(RP), intent(out) :: pres_org(:,:,:)
216  real(RP), intent(out) :: dens_org(:,:,:)
217  real(RP), intent(out) :: pott_org(:,:,:)
218  real(RP), intent(out) :: qv_org (:,:,:)
219  real(RP), intent(out) :: qtrc_org(:,:,:,:)
220  real(RP), intent(in) :: cz_org (:,:,:)
221  character(len=*), intent(in) :: basename_org
222  logical, intent(in) :: same_mptype
223  integer, intent(in) :: dims(6)
224  integer, intent(in) :: it
225 
226  real(RP) :: momz_org(dims(1)+2,dims(2),dims(3))
227  real(RP) :: momx_org(dims(1)+2,dims(2),dims(3))
228  real(RP) :: momy_org(dims(1)+2,dims(2),dims(3))
229  real(RP) :: rhot_org(dims(1)+2,dims(2),dims(3))
230  real(RP) :: tsfc_org( dims(2),dims(3))
231  real(RP) :: qhyd_org(dims(1)+2,dims(2),dims(3),n_hyd)
232  real(RP) :: qnum_org(dims(1)+2,dims(2),dims(3),n_hyd)
233  real(RP) :: qdry (dims(1)+2)
234  real(RP) :: Rtot (dims(1)+2)
235  real(RP) :: CVtot (dims(1)+2)
236  real(RP) :: CPtot (dims(1)+2)
237  real(RP) :: temp_org(dims(1)+2)
238  real(RP) :: dz
239 
240  integer :: xs, xe
241  integer :: ys, ye
242  integer :: nx, ny
243  integer :: xloc, yloc
244  integer :: rank
245 
246  integer :: fid
247  logical :: existed, existed_t2, existed_mslp
248  integer :: k, i, j, iq
249  !---------------------------------------------------------------------------
250 
251  do i = 1, size( nest_tile_id(:) )
252  ! read data from split files
253  rank = nest_tile_id(i)
254 
255  xloc = mod( i-1, nest_tile_num_x ) + 1
256  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
257 
258  xs = parent_imax(handle) * (xloc-1) + 1
259  xe = parent_imax(handle) * xloc
260  ys = parent_jmax(handle) * (yloc-1) + 1
261  ye = parent_jmax(handle) * yloc
262 
263  nx = xe - xs + 1
264  ny = ye - ys + 1
265 
266  call file_open( basename_org, & ! (in)
267  fid, & ! (out)
268  aggregate=.false., rankid=rank ) ! (in)
269 
270  call file_cartesc_read( fid, "T2", read2d(:,:), step=it, existed=existed_t2 )
271 !OCL XFILL
272  if ( existed_t2 ) tsfc_org(xs:xe,ys:ye) = read2d(:,:)
273 
274  call file_cartesc_read( fid, "MSLP", read2d(:,:), step=it, existed=existed_mslp )
275 !OCL XFILL
276  if ( existed_mslp ) pres_org(1,xs:xe,ys:ye) = read2d(:,:)
277 
278  call file_cartesc_read( fid, "DENS", read3d(:,:,:), step=it )
279 !OCL XFILL
280  dens_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
281 
282  call file_cartesc_read( fid, "MOMZ", read3d(:,:,:), step=it )
283 !OCL XFILL
284  momz_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
285 
286  call file_cartesc_read( fid, "MOMX", read3d(:,:,:), step=it )
287 !OCL XFILL
288  momx_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
289 
290  call file_cartesc_read( fid, "MOMY", read3d(:,:,:), step=it )
291 !OCL XFILL
292  momy_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
293 
294 
295  call file_cartesc_read( fid, "RHOT", read3d(:,:,:), step=it )
296 !OCL XFILL
297  rhot_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
298 
299 !OCL XFILL
300  do iq = 1, n_hyd
301  qhyd_org(:,xs:xe,ys:ye,iq) = 0.0_rp
302  qnum_org(:,xs:xe,ys:ye,iq) = 0.0_rp
303  end do
304 
305  if ( same_mptype ) then
306 
307  do iq = qs_mp, qe_mp
308  call file_cartesc_read( fid, tracer_name(iq), read3d(:,:,:), step=it )
309 !OCL XFILL
310  qtrc_org(3:dims(1)+2,xs:xe,ys:ye,iq) = read3d(:,:,:)
311 !OCL XFILL
312  qtrc_org(2,xs:xe,ys:ye,iq) = qtrc_org(3,xs:xe,ys:ye,iq)
313 !OCL XFILL
314  qtrc_org(1,xs:xe,ys:ye,iq) = qtrc_org(3,xs:xe,ys:ye,iq)
315  enddo
316 
317  else
318 
319  call file_cartesc_read( fid, "QV", read3d(:,:,:), step=it, existed=existed )
320  if ( existed ) then
321 !OCL XFILL
322  qv_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
323 !OCL XFILL
324  qv_org(2,xs:xe,ys:ye) = qv_org(3,xs:xe,ys:ye)
325 !OCL XFILL
326  qv_org(1,xs:xe,ys:ye) = qv_org(3,xs:xe,ys:ye)
327  else
328 !OCL XFILL
329  qv_org(:,:,:) = 0.0_rp
330  end if
331 
332  do iq = 1, n_hyd
333  call file_cartesc_read( fid, hyd_name(iq), read3d(:,:,:), step=it, existed=existed )
334  if ( existed ) then
335 !OCL XFILL
336  qhyd_org(3:dims(1)+2,xs:xe,ys:ye,iq) = read3d(:,:,:)
337 !OCL XFILL
338  qhyd_org(2,xs:xe,ys:ye,iq) = qhyd_org(3,xs:xe,ys:ye,iq)
339 !OCL XFILL
340  qhyd_org(1,xs:xe,ys:ye,iq) = qhyd_org(3,xs:xe,ys:ye,iq)
341  else
342 !OCL XFILL
343  qhyd_org(:,:,:,iq) = 0.0_rp
344  end if
345 
346  ! number density
347  call file_cartesc_read( fid, num_name(iq), read3d(:,:,:), step=it, existed=existed )
348  if ( existed ) then
349 !OCL XFILL
350  qnum_org(3:dims(1)+2,xs:xe,ys:ye,iq) = read3d(:,:,:)
351 !OCL XFILL
352  qnum_org(2,xs:xe,ys:ye,iq) = qnum_org(3,xs:xe,ys:ye,iq)
353 !OCL XFILL
354  qnum_org(1,xs:xe,ys:ye,iq) = qnum_org(3,xs:xe,ys:ye,iq)
355  else
356 !OCL XFILL
357  qnum_org(:,:,:,iq) = 0.0_rp
358  end if
359  end do
360 
361  endif
362 
363 ! call FILE_CARTESC_read( fid, "Q2", read2D(:,:), step=it )
364 ! qv_org(2,xs:xe,ys:ye) = read2D(:,:)
365 
366  end do
367 
368  ! convert from momentum to velocity
369 !OCL XFILL
370  do j = 1, dims(3)
371  do i = 1, dims(2)
372  do k = 4, dims(1)+2
373  velz_org(k,i,j) = ( momz_org(k-1,i,j) + momz_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
374  end do
375  end do
376  end do
377 !OCL XFILL
378  do j = 1, dims(3)
379  do i = 1, dims(2)
380  velz_org(1:3 ,i,j) = 0.0_rp
381  velz_org(dims(1)+2,i,j) = 0.0_rp
382  end do
383  end do
384 
385  ! convert from momentum to velocity
386 !OCL XFILL
387  do j = 1, dims(3)
388  do i = 2, dims(2)
389  do k = 3, dims(1)+2
390  velx_org(k,i,j) = ( momx_org(k,i-1,j) + momx_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
391  end do
392  end do
393  end do
394 !OCL XFILL
395  do j = 1, dims(3)
396  do k = 3, dims(1)+2
397  velx_org(k,1,j) = momx_org(k,1,j) / dens_org(k,1,j)
398  end do
399  end do
400 !OCL XFILL
401  velx_org(1:2,:,:) = 0.0_rp
402 
403  ! convert from momentum to velocity
404 !OCL XFILL
405  do j = 2, dims(3)
406  do i = 1, dims(2)
407  do k = 3, dims(1)+2
408  vely_org(k,i,j) = ( momy_org(k,i,j-1) + momy_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
409  end do
410  end do
411  end do
412 !OCL XFILL
413  do i = 1, dims(2)
414  do k = 3, dims(1)+2
415  vely_org(k,i,1) = momy_org(k,i,1) / dens_org(k,i,1)
416  end do
417  end do
418 !OCL XFILL
419  vely_org(1:2,:,:) = 0.0_rp
420 
421 
422  !!! must be rotate !!!
423 
424 
425  if ( qa_mp > 0 .AND. .NOT. same_mptype ) then
426  call atmos_phy_mp_driver_qhyd2qtrc( dims(1)+2, 1, dims(1)+2, dims(2), 1, dims(2), dims(3), 1, dims(3), &
427  qv_org(:,:,:), qhyd_org(:,:,:,:), & ! [IN]
428  qtrc_org(:,:,:,qs_mp:qe_mp), & ! [OUT]
429  qnum=qnum_org(:,:,:,:) ) ! [IN]
430  end if
431 
432  !$omp parallel do default(none) &
433  !$omp private(dz,temp_org,qdry,rtot,cvtot,cptot) &
434  !$omp shared(dims,QA,LAPS,P00,Rdry,CPdry,GRAV,TRACER_MASS,TRACER_R,TRACER_CV,TRACER_CP, &
435  !$omp dens_org,rhot_org,qtrc_org,pres_org,pott_org,tsfc_org,cz_org, &
436  !$omp existed_t2,existed_mslp)
437  do j = 1, dims(3)
438  do i = 1, dims(2)
439 
440  ! diagnose temp and pres
441  do k = 3, dims(1)+2
442 
443  call thermodyn_specific_heat( dims(1)+2, 3, dims(1)+2, qa, & ! [IN]
444  qtrc_org(:,i,j,:), & ! [IN]
445  tracer_mass(:), tracer_r(:), tracer_cv(:), tracer_cp(:), & ! [IN]
446  qdry(:), rtot(:), cvtot(:), cptot(:) ) ! [OUT]
447 
448  call thermodyn_rhot2temp_pres( dims(1)+2, 3, dims(1)+2, & ! [IN]
449  dens_org(:,i,j), rhot_org(:,i,j), rtot, cvtot, cptot, & ! [IN]
450  temp_org(:), pres_org(:,i,j) ) ! [OUT]
451 
452  end do
453 
454 !OCL XFILL
455  do k = 3, dims(1)+2
456  pott_org(k,i,j) = rhot_org(k,i,j) / dens_org(k,i,j)
457  end do
458 
459  ! at the surface
460  dz = cz_org(3,i,j) - cz_org(2,i,j)
461 
462  if ( .not. existed_t2 ) then
463  tsfc_org(i,j) = temp_org(3) + laps * dz
464  end if
465  dens_org(2,i,j) = ( pres_org(3,i,j) + grav * dens_org(3,i,j) * dz * 0.5_rp ) &
466  / ( rdry * tsfc_org(i,j) - grav * dz * 0.5_rp )
467  pres_org(2,i,j) = dens_org(2,i,j) * rdry * tsfc_org(i,j)
468  pott_org(2,i,j) = tsfc_org(i,j) * ( p00 / pres_org(2,i,j) )**(rdry/cpdry)
469 
470  ! at the sea-level
471  temp_org(1) = tsfc_org(i,j) + laps * cz_org(2,i,j)
472  if ( .not. existed_mslp ) then
473  dens_org(1,i,j) = ( pres_org(2,i,j) + grav * dens_org(2,i,j) * cz_org(2,i,j) * 0.5_rp ) &
474  / ( rdry * temp_org(1) - grav * cz_org(2,i,j) * 0.5_rp )
475  pres_org(1,i,j) = dens_org(1,i,j) * rdry * temp_org(1)
476  pott_org(1,i,j) = temp_org(1) * ( p00 / pres_org(1,i,j) )**(rdry/cpdry)
477  else
478  pott_org(1,i,j) = temp_org(1) * ( p00 / pres_org(1,i,j) )**(rdry/cpdry)
479  dens_org(1,i,j) = pres_org(1,i,j) / ( rdry * temp_org(1) )
480  end if
481 
482  end do
483  end do
484 
485  return
486  end subroutine parentatmosinputscale
487 
488  !-----------------------------------------------------------------------------
490  subroutine parentlandsetupscale( &
491  ldims )
493  parent_imax, &
494  parent_jmax, &
495  parent_lkmax, &
496  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
497  nest_tile_num_y => comm_cartesc_nest_tile_num_y
498  implicit none
499 
500  integer, intent(out) :: ldims(3)
501  !---------------------------------------------------------------------------
502 
503  log_info("ParentLandSetupSCALE",*) 'Real Case/Land Input File Type: SCALE-RM'
504 
505  ldims(1) = parent_lkmax(handle)
506  ldims(2) = parent_imax(handle) * nest_tile_num_x
507  ldims(3) = parent_jmax(handle) * nest_tile_num_y
508 
509  if ( .NOT. allocated(read2d) ) then
510  allocate( read2d(parent_imax(handle),parent_jmax(handle)) )
511  endif
512 
513  allocate( read3dl(ldims(1),parent_imax(handle),parent_jmax(handle)) )
514 
515  return
516  end subroutine parentlandsetupscale
517 
518  !-----------------------------------------------------------------------------
519  subroutine parentlandinputscale( &
520  tg_org, &
521  strg_org, &
522  lst_org, &
523  ust_org, &
524  albg_org, &
525  topo_org, &
526  lmask_org, &
527  llon_org, &
528  llat_org, &
529  lz_org, &
530  basename_land, &
531  ldims, &
532  use_file_landwater, &
533  it )
534  use scale_const, only: &
535  undef => const_undef, &
536  d2r => const_d2r
537  use scale_comm_cartesc_nest, only: &
538  parent_imax, &
539  parent_jmax, &
540  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
541  nest_tile_id => comm_cartesc_nest_tile_id
542  use scale_file, only: &
543  file_open, &
544  file_read, &
545  file_get_datainfo
546  use scale_file_cartesc, only: &
547  file_cartesc_read
548  implicit none
549 
550  real(RP), intent(out) :: tg_org (:,:,:)
551  real(RP), intent(out) :: strg_org (:,:,:)
552  real(RP), intent(out) :: lst_org (:,:)
553  real(RP), intent(out) :: ust_org (:,:)
554  real(RP), intent(out) :: albg_org (:,:,:,:)
555  real(RP), intent(out) :: topo_org (:,:)
556  real(RP), intent(out) :: lmask_org(:,:)
557  real(RP), intent(out) :: llon_org (:,:)
558  real(RP), intent(out) :: llat_org (:,:)
559  real(RP), intent(out) :: lz_org (:)
560  character(len=*), intent(in) :: basename_land
561  integer, intent(in) :: ldims(3)
562  logical, intent(in) :: use_file_landwater ! use land water data from files
563  integer, intent(in) :: it
564 
565  integer :: rank
566  integer :: fid
567  integer :: xloc, yloc
568  integer :: xs, xe
569  integer :: ys, ye
570  logical :: existed
571 
572  integer :: i
573  !---------------------------------------------------------------------------
574 
575  do i = 1, size( nest_tile_id(:) )
576  ! read data from split files
577  rank = nest_tile_id(i)
578 
579  xloc = mod( i-1, nest_tile_num_x ) + 1
580  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
581 
582  xs = parent_imax(handle) * (xloc-1) + 1
583  xe = parent_imax(handle) * xloc
584  ys = parent_jmax(handle) * (yloc-1) + 1
585  ye = parent_jmax(handle) * yloc
586 
587  call file_open( basename_land, & ! (in)
588  fid, & ! (out)
589  aggregate=.false., rankid=rank ) ! (in)
590 
591  call file_cartesc_read( fid, "LAND_TEMP", read3dl(:,:,:), step=it )
592  tg_org(1:ldims(1),xs:xe,ys:ye) = read3dl(:,:,:)
593 
594  if( use_file_landwater )then
595  call file_cartesc_read( fid, "LAND_WATER", read3dl(:,:,:), step=it )
596  strg_org(1:ldims(1),xs:xe,ys:ye) = read3dl(:,:,:)
597  endif
598 
599  call file_cartesc_read( fid, "lon", read2d(:,:) )
600  llon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
601 
602  call file_cartesc_read( fid, "lat", read2d(:,:) )
603  llat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
604 
605  call file_cartesc_read( fid, "LAND_SFC_TEMP", read2d(:,:), step=it )
606  lst_org(xs:xe,ys:ye) = read2d(:,:)
607 
608  call file_get_datainfo( fid, "URBAN_SFC_TEMP", istep=it, existed=existed )
609  if ( existed ) then
610  call file_cartesc_read( fid, "URBAN_SFC_TEMP", read2d(:,:), step=it )
611  ust_org(xs:xe,ys:ye) = read2d(:,:)
612  else
613  ust_org(xs:xe,ys:ye) = undef
614  end if
615 
616  call file_cartesc_read( fid, "LAND_SFC_ALB_IR_dir", read2d(:,:), step=it )
617  albg_org(xs:xe,ys:ye,i_r_direct ,i_r_ir ) = read2d(:,:)
618 
619  call file_cartesc_read( fid, "LAND_SFC_ALB_IR_dif", read2d(:,:), step=it )
620  albg_org(xs:xe,ys:ye,i_r_diffuse,i_r_ir ) = read2d(:,:)
621 
622  call file_cartesc_read( fid, "LAND_SFC_ALB_NIR_dir", read2d(:,:), step=it )
623  albg_org(xs:xe,ys:ye,i_r_direct ,i_r_nir) = read2d(:,:)
624 
625  call file_cartesc_read( fid, "LAND_SFC_ALB_NIR_dif", read2d(:,:), step=it )
626  albg_org(xs:xe,ys:ye,i_r_diffuse,i_r_nir) = read2d(:,:)
627 
628  call file_cartesc_read( fid, "LAND_SFC_ALB_VIS_dir", read2d(:,:), step=it )
629  albg_org(xs:xe,ys:ye,i_r_direct ,i_r_vis) = read2d(:,:)
630 
631  call file_cartesc_read( fid, "LAND_SFC_ALB_VIS_dif", read2d(:,:), step=it )
632  albg_org(xs:xe,ys:ye,i_r_diffuse,i_r_vis) = read2d(:,:)
633 
634  call file_cartesc_read( fid, "topo", read2d(:,:) )
635  topo_org(xs:xe,ys:ye) = read2d(:,:)
636 
637  call file_cartesc_read( fid, "lsmask", read2d(:,:) )
638  lmask_org(xs:xe,ys:ye) = read2d(:,:)
639 
640  end do
641 
642  call file_read( fid, "lz", lz_org(:) )
643 
644  return
645  end subroutine parentlandinputscale
646 
647  !-----------------------------------------------------------------------------
649  subroutine parentoceansetupscale( &
650  odims )
652  parent_imax, &
653  parent_jmax, &
654  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
655  nest_tile_num_y => comm_cartesc_nest_tile_num_y
656  implicit none
657 
658  integer, intent(out) :: odims(2)
659  !---------------------------------------------------------------------------
660 
661  log_info("ParentOceanSetupSCALE",*) 'Real Case/Ocean Input File Type: SCALE-RM'
662 
663  odims(1) = parent_imax(handle) * nest_tile_num_x
664  odims(2) = parent_jmax(handle) * nest_tile_num_y
665 
666  if ( .NOT. allocated(read2d) ) then
667  allocate( read2d(parent_imax(handle),parent_jmax(handle)) )
668  end if
669 
670  allocate( read3do(1,parent_imax(handle),parent_jmax(handle)) )
671 
672  return
673  end subroutine parentoceansetupscale
674 
675  !-----------------------------------------------------------------------------
676  subroutine parentoceanopenscale( &
677  olon_org, &
678  olat_org, &
679  omask_org, &
680  basename_ocean, &
681  odims )
682  use scale_const, only: &
683  d2r => const_d2r
684  use scale_comm_cartesc_nest, only: &
685  parent_imax, &
686  parent_jmax, &
687  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
688  nest_tile_id => comm_cartesc_nest_tile_id
689  use scale_file, only: &
690  file_open
691  use scale_file_cartesc, only: &
692  file_cartesc_read
693  implicit none
694 
695  real(RP), intent(out) :: olon_org (:,:)
696  real(RP), intent(out) :: olat_org (:,:)
697  real(RP), intent(out) :: omask_org(:,:)
698  character(len=*), intent(in) :: basename_ocean
699  integer, intent(in) :: odims(2)
700 
701  integer :: rank
702  integer :: xloc, yloc
703  integer :: xs, xe, ys, ye
704 
705  integer :: fid
706  integer :: i
707  !---------------------------------------------------------------------------
708 
709  do i = 1, size( nest_tile_id(:) )
710  ! read data from split files
711  rank = nest_tile_id(i)
712 
713  xloc = mod( i-1, nest_tile_num_x ) + 1
714  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
715 
716  xs = parent_imax(handle) * (xloc-1) + 1
717  xe = parent_imax(handle) * xloc
718  ys = parent_jmax(handle) * (yloc-1) + 1
719  ye = parent_jmax(handle) * yloc
720 
721  call file_open( basename_ocean, & ! (in)
722  fid, & ! (out)
723  aggregate=.false., rankid=rank ) ! (in)
724 
725  call file_cartesc_read( fid, "lon", read2d(:,:) )
726  olon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
727 
728  call file_cartesc_read( fid, "lat", read2d(:,:) )
729  olat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
730 
731  call file_cartesc_read( fid, "lsmask", read2d(:,:) )
732  omask_org(xs:xe,ys:ye) = read2d(:,:)
733 
734  end do
735 
736  return
737  end subroutine parentoceanopenscale
738 
739  !-----------------------------------------------------------------------------
740  subroutine parentoceaninputscale( &
741  tw_org, &
742  sst_org, &
743  albw_org, &
744  z0w_org, &
745  omask_org, &
746  basename_ocean, &
747  odims, &
748  it )
750  parent_imax, &
751  parent_jmax, &
752  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
753  nest_tile_id => comm_cartesc_nest_tile_id
754  use scale_file, only: &
755  file_open, &
756  file_get_datainfo
757  use scale_file_cartesc, only: &
758  file_cartesc_read
759  implicit none
760 
761  real(RP), intent(out) :: tw_org (:,:)
762  real(RP), intent(out) :: sst_org (:,:)
763  real(RP), intent(out) :: albw_org (:,:,:,:)
764  real(RP), intent(out) :: z0w_org (:,:)
765  real(RP), intent(out) :: omask_org(:,:)
766  character(len=*), intent(in) :: basename_ocean
767  integer, intent(in) :: odims(2)
768  integer, intent(in) :: it
769 
770  integer :: rank
771 
772  integer :: fid
773  integer :: ndim
774  integer :: i
775  integer :: xloc, yloc
776  integer :: xs, xe
777  integer :: ys, ye
778  !---------------------------------------------------------------------------
779 
780  do i = 1, size( nest_tile_id(:) )
781  ! read data from split files
782  rank = nest_tile_id(i)
783 
784  xloc = mod( i-1, nest_tile_num_x ) + 1
785  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
786 
787  xs = parent_imax(handle) * (xloc-1) + 1
788  xe = parent_imax(handle) * xloc
789  ys = parent_jmax(handle) * (yloc-1) + 1
790  ye = parent_jmax(handle) * yloc
791 
792  call file_open( basename_ocean, & ! (in)
793  fid, & ! (out)
794  aggregate=.false., rankid=rank ) ! (in)
795 
796  call file_get_datainfo( fid, "OCEAN_TEMP", dim_rank=ndim )
797  select case (ndim)
798  case ( 2 ) ! for old file
799  call file_cartesc_read( fid, "OCEAN_TEMP", read2d(:,:), step=it )
800  tw_org(xs:xe,ys:ye) = read2d(:,:)
801  case ( 3 )
802  call file_cartesc_read( fid, "OCEAN_TEMP", read3do(:,:,:), step=it )
803  tw_org(xs:xe,ys:ye) = read3do(1,:,:)
804  end select
805 
806  call file_cartesc_read( fid, "OCEAN_SFC_TEMP", read2d(:,:), step=it )
807  sst_org(xs:xe,ys:ye) = read2d(:,:)
808 
809  call file_cartesc_read( fid, "OCEAN_SFC_ALB_IR_dir", read2d(:,:), step=it )
810  albw_org(xs:xe,ys:ye,i_r_direct ,i_r_ir ) = read2d(:,:)
811 
812  call file_cartesc_read( fid, "OCEAN_SFC_ALB_IR_dif", read2d(:,:), step=it )
813  albw_org(xs:xe,ys:ye,i_r_diffuse,i_r_ir ) = read2d(:,:)
814 
815  call file_cartesc_read( fid, "OCEAN_SFC_ALB_NIR_dir", read2d(:,:), step=it )
816  albw_org(xs:xe,ys:ye,i_r_direct ,i_r_nir) = read2d(:,:)
817 
818  call file_cartesc_read( fid, "OCEAN_SFC_ALB_NIR_dif", read2d(:,:), step=it )
819  albw_org(xs:xe,ys:ye,i_r_diffuse,i_r_nir) = read2d(:,:)
820 
821  call file_cartesc_read( fid, "OCEAN_SFC_ALB_VIS_dir", read2d(:,:), step=it )
822  albw_org(xs:xe,ys:ye,i_r_direct ,i_r_vis) = read2d(:,:)
823 
824  call file_cartesc_read( fid, "OCEAN_SFC_ALB_VIS_dif", read2d(:,:), step=it )
825  albw_org(xs:xe,ys:ye,i_r_diffuse,i_r_vis) = read2d(:,:)
826 
827  call file_cartesc_read( fid, "OCEAN_SFC_Z0M", read2d(:,:), step=it )
828  z0w_org(xs:xe,ys:ye) = read2d(:,:)
829 
830  call file_cartesc_read( fid, "lsmask", read2d(:,:) )
831  omask_org(xs:xe,ys:ye) = read2d(:,:)
832 
833  end do
834 
835  return
836  end subroutine parentoceaninputscale
837 
838 end module mod_realinput_scale
integer, dimension(2), public parent_kmax
parent max number in z-direction
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
character(len=h_short), dimension(n_hyd), parameter, public hyd_name
module coupler / surface-atmospehre
integer, public comm_cartesc_nest_tile_num_y
parent tile number in y-direction
real(rp), dimension(qa_max), public tracer_r
module Atmosphere / Physics Cloud Microphysics
subroutine, public parentoceanopenscale(olon_org, olat_org, omask_org, basename_ocean, odims)
integer, parameter, public i_r_vis
module Atmosphere Grid CartesianC metirc
integer, public qa
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
subroutine, public parentlandsetupscale(ldims)
Land Setup.
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:58
character(len=h_short), dimension(qa_max), public tracer_name
real(rp), dimension(qa_max), public tracer_cv
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
real(rp), public const_undef
Definition: scale_const.F90:41
module COMMUNICATION
real(rp), dimension(qa_max), public tracer_cp
module file
Definition: scale_file.F90:15
module TRACER
module atmosphere / hydrometeor
integer, dimension(2), public parent_jmax
parent max number in y-direction
subroutine, public atmos_phy_mp_driver_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
subroutine, public parentoceaninputscale(tw_org, sst_org, albw_org, z0w_org, omask_org, basename_ocean, odims, it)
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
module REAL input SCALE
subroutine, public parentlandinputscale(tg_org, strg_org, lst_org, ust_org, albg_org, topo_org, lmask_org, llon_org, llat_org, lz_org, basename_land, ldims, use_file_landwater, it)
module atmosphere / hydrostatic barance
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
subroutine, public parentatmossetupscale(dims)
Atmos Setup.
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
module CONSTANT
Definition: scale_const.F90:11
module Communication CartesianC nesting
integer, parameter, public i_r_direct
integer, dimension(2), public parent_lkmax
parent max number in lz-direction
subroutine, public parentatmosinputscale(velz_org, velx_org, vely_org, pres_org, dens_org, pott_org, qv_org, qtrc_org, cz_org, basename_org, same_mptype, dims, it)
integer, dimension(:), allocatable, public comm_cartesc_nest_tile_id
parent tile real id
integer, parameter, public i_r_nir
module atmosphere / thermodyn
subroutine, public parentoceansetupscale(odims)
Ocean Setup.
module PRECISION
module file / cartesianC
integer, dimension(2), public parent_imax
parent max number in x-direction
subroutine, public parentatmosopenscale(lon_org, lat_org, cz_org, basename_org, dims)
integer, parameter, public i_r_ir
module STDIO
Definition: scale_io.F90:10
character(len=h_short), dimension(n_hyd), parameter, public num_name
integer, parameter, public i_r_diffuse
integer, parameter, public n_hyd
module atmosphere / physics / cloud microphysics
integer, public comm_cartesc_nest_tile_num_x
parent tile number in x-direction
real(rp), dimension(qa_max), public tracer_mass