SCALE-RM
mod_realinput_scale.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_grid_nest, only: &
25  parent_kmax, &
26  parent_imax, &
27  parent_jmax, &
28  parent_lkmax, &
32 
33  !-----------------------------------------------------------------------------
34  implicit none
35  private
36  !-----------------------------------------------------------------------------
37  !
38  !++ Public procedure
39  !
40  public :: parentatomsetupscale
41  public :: parentatomopenscale
42  public :: parentatominputscale
43  public :: parentlandsetupscale
44  public :: parentlandinputscale
45  public :: parentoceansetupscale
46  public :: parentoceanopenscale
47  public :: parentoceaninputscale
48 
49  !-----------------------------------------------------------------------------
50  !
51  !++ Public parameters & variables
52  !
53  !-----------------------------------------------------------------------------
54  !
55  !++ Private procedure
56  !
57  !-----------------------------------------------------------------------------
58  !
59  !++ Private parameters & variables
60  !
61  integer, parameter :: handle = 1
62 
63  real(RP), allocatable :: read2D(:,:)
64  real(RP), allocatable :: read3D(:,:,:)
65  real(RP), allocatable :: read3DL(:,:,:)
66 
67  !-----------------------------------------------------------------------------
68 contains
69  !-----------------------------------------------------------------------------
71  subroutine parentatomsetupscale( &
72  dims )
73  implicit none
74 
75  integer, intent(out) :: dims(6)
76 
77  integer :: i
78  !---------------------------------------------------------------------------
79 
80  if( io_l ) write(io_fid_log,*) '+++ Real Case/Atom Input File Type: SCALE-RM'
81  ! full level
82  dims(1) = parent_kmax(handle)
83  dims(2) = parent_imax(handle) * nest_tile_num_x
84  dims(3) = parent_jmax(handle) * nest_tile_num_y
85  ! half level
86  dims(4) = dims(1)
87  dims(5) = dims(2)
88  dims(6) = dims(3)
89 
90  allocate( read2d( parent_imax(handle), parent_jmax(handle) ) )
91  allocate( read3d( parent_imax(handle), parent_jmax(handle), dims(1) ) )
92 
93  return
94  end subroutine parentatomsetupscale
95 
96  !-----------------------------------------------------------------------------
97  subroutine parentatomopenscale( &
98  lon_org, &
99  lat_org, &
100  cz_org, &
101  basename_org, &
102  dims )
103  use scale_const, only: &
104  d2r => const_d2r
105  use gtool_file, only: &
106  fileread
107  implicit none
108  real(RP), intent(out) :: lon_org(:,:)
109  real(RP), intent(out) :: lat_org(:,:)
110  real(RP), intent(out) :: cz_org (:,:,:)
111  character(len=*), intent(in) :: basename_org
112  integer, intent(in) :: dims(6)
113 
114  integer :: rank
115  integer :: xloc, yloc
116  integer :: xs, xe, ys, ye
117 
118  integer :: i, k
119 
120  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[AtomOpenSCALE]'
121 
122  do i = 1, size( nest_tile_id(:) )
123  ! read data from split files
124  rank = nest_tile_id(i)
125 
126  xloc = mod( i-1, nest_tile_num_x ) + 1
127  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
128 
129  xs = parent_imax(handle) * (xloc-1) + 1
130  xe = parent_imax(handle) * xloc
131  ys = parent_jmax(handle) * (yloc-1) + 1
132  ye = parent_jmax(handle) * yloc
133 
134  call fileread( read2d(:,:), basename_org, "lon", 1, rank )
135  lon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
136 
137  call fileread( read2d(:,:), basename_org, "lat", 1, rank )
138  lat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
139 
140  call fileread( read3d(:,:,:), basename_org, "height", 1, rank )
141  do k = 1, dims(1)
142  cz_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
143  end do
144 
145  call fileread( read2d(:,:), basename_org, "topo", 1, rank )
146  cz_org(2,xs:xe,ys:ye) = read2d(:,:)
147 
148  end do
149 
150  cz_org(1,:,:) = 0.0_rp
151  return
152  end subroutine parentatomopenscale
153 
154  !-----------------------------------------------------------------------------
155  subroutine parentatominputscale( &
156  velz_org, &
157  velx_org, &
158  vely_org, &
159  pres_org, &
160  dens_org, &
161  pott_org, &
162  qtrc_org, &
163  cz_org, &
164  flg_bin, &
165  flg_intrp, &
166  basename_org, &
167  mptype_parent, &
168  dims, &
169  it ) ! (in)
170  use scale_const, only: &
171  p00 => const_pre00, &
172  cpdry => const_cpdry, &
173  rdry => const_rdry, &
174  grav => const_grav, &
175  laps => const_laps
176  use scale_comm, only: &
177  comm_vars8, &
178  comm_wait
179  use gtool_file, only: &
180  fileread
181  use scale_atmos_hydrostatic, only: &
182  hydrostatic_buildrho_real => atmos_hydrostatic_buildrho_real
183  use scale_atmos_thermodyn, only: &
184  thermodyn_temp_pres => atmos_thermodyn_temp_pres, &
185  thermodyn_pott => atmos_thermodyn_pott
186  use scale_atmos_phy_mp, only: &
187  qs_mp
188  use scale_gridtrans, only: &
189  rotc => gtrans_rotc
190  use scale_topography, only: &
191  topo => topo_zsfc
192  use scale_interpolation_nest, only: &
196  use scale_atmos_phy_mp_convert, only: &
198  implicit none
199 
200  real(RP), intent(out) :: velz_org(:,:,:)
201  real(RP), intent(out) :: velx_org(:,:,:)
202  real(RP), intent(out) :: vely_org(:,:,:)
203  real(RP), intent(out) :: pres_org(:,:,:)
204  real(RP), intent(out) :: dens_org(:,:,:)
205  real(RP), intent(out) :: pott_org(:,:,:)
206  real(RP), intent(out) :: qtrc_org(:,:,:,:)
207  real(RP), intent(in) :: cz_org(:,:,:)
208  logical, intent(in) :: flg_bin ! flag for SBM(S10)
209  logical, intent(in) :: flg_intrp ! flag for interpolation of SBM
210  character(len=*), intent(in) :: basename_org
211  integer, intent(in) :: mptype_parent
212  integer, intent(in) :: dims(6)
213  integer, intent(in) :: it
214 
215  ! work
216 
217  real(RP) :: momz_org(dims(1)+2,dims(2),dims(3))
218  real(RP) :: momx_org(dims(1)+2,dims(2),dims(3))
219  real(RP) :: momy_org(dims(1)+2,dims(2),dims(3))
220  real(RP) :: rhot_org(dims(1)+2,dims(2),dims(3))
221  real(RP) :: tsfc_org( dims(2),dims(3))
222  real(RP) :: temp_org
223  real(RP) :: dz
224 
225  integer :: xs, xe
226  integer :: ys, ye
227  integer :: xloc, yloc
228  integer :: rank
229 
230  integer :: k, i, j, iq, iqa
231  !---------------------------------------------------------------------------
232 
233 
234  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[AtomInputSCALE]'
235 
236  do i = 1, size( nest_tile_id(:) )
237  ! read data from split files
238  rank = nest_tile_id(i)
239 
240  xloc = mod( i-1, nest_tile_num_x ) + 1
241  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
242 
243  xs = parent_imax(handle) * (xloc-1) + 1
244  xe = parent_imax(handle) * xloc
245  ys = parent_jmax(handle) * (yloc-1) + 1
246  ye = parent_jmax(handle) * yloc
247 
248  call fileread( read2d(:,:), basename_org, "T2", it, rank )
249  tsfc_org(xs:xe,ys:ye) = read2d(:,:)
250 
251  call fileread( read2d(:,:), basename_org, "MSLP", it, rank )
252  pres_org(1,xs:xe,ys:ye) = read2d(:,:)
253 
254  call fileread( read3d(:,:,:), basename_org, "DENS", it, rank )
255  do k = 1, dims(1)
256  dens_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
257  end do
258 
259  call fileread( read3d(:,:,:), basename_org, "MOMZ", it, rank )
260  do k = 1, dims(1)
261  momz_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
262  end do
263 
264  call fileread( read3d(:,:,:), basename_org, "MOMX", it, rank )
265  do k = 1, dims(1)
266  momx_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
267  end do
268 
269  call fileread( read3d(:,:,:), basename_org, "MOMY", it, rank )
270  do k = 1, dims(1)
271  momy_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
272  end do
273 
274  call fileread( read3d(:,:,:), basename_org, "RHOT", it, rank )
275  do k = 1, dims(1)
276  rhot_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
277  end do
278 
279  do iq = 1, qa
280  qtrc_org(:,xs:xe,ys:ye,iq) = 0.0_rp
281  end do
282 
283  if( flg_bin .and. flg_intrp ) then !--- tracers created from parent domain by interpolation
284 
285  if( io_l ) write(io_fid_log,*) '+++ SDF of SBM(S10) is interpolated from Qxx and Nxx'
286  if( io_l ) write(io_fid_log,*) '+++ of outer bulk MP model'
287 
288  call atmos_phy_mp_bulk2bin( xs, xe, & ! [IN]
289  ys, ye, & ! [IN]
290  dims, & ! [IN]
291  it, & ! [IN]
292  rank, & ! [IN]
293  handle, & ! [IN]
294  basename_org, & ! [IN]
295  dens_org, & ! [IN]
296  qtrc_org ) ! [INOUT]
297  do iq = 1, qa
298  qtrc_org(2,xs:xe,ys:ye,iq) = qtrc_org(3,xs:xe,ys:ye,iq)
299  enddo
300 
301  else !--- tracers of paremt domain directly used
302 
303  do iq = 1, mptype_parent
304  iqa = qs_mp + iq - 1
305  call fileread( read3d(:,:,:), basename_org, tracer_name(iq), it, rank )
306  do k = 1, dims(1)
307  qtrc_org(k+2,xs:xe,ys:ye,iqa) = read3d(:,:,k)
308  end do
309  qtrc_org(2,xs:xe,ys:ye,iqa) = qtrc_org(3,xs:xe,ys:ye,iqa)
310  qtrc_org(1,xs:xe,ys:ye,iqa) = qtrc_org(3,xs:xe,ys:ye,iqa)
311  end do
312 
313  endif
314 
315 ! call FileRead( read2D(:,:), BASENAME_ORG, "Q2", it, rank )
316 ! qtrc_org(2,xs:xe,ys:ye,I_QV) = read2D(:,:)
317 
318  end do
319 
320  ! convert from momentum to velocity
321  do j = 1, dims(3)
322  do i = 1, dims(2)
323  do k = 4, dims(1)+2
324  velz_org(k,i,j) = ( momz_org(k-1,i,j) + momz_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
325  end do
326  end do
327  end do
328  do j = 1, dims(3)
329  do i = 1, dims(2)
330  velz_org(1:3 ,i,j) = 0.0_rp
331  velz_org(dims(1)+2,i,j) = 0.0_rp
332  end do
333  end do
334 
335  ! convert from momentum to velocity
336  do j = 1, dims(3)
337  do i = 2, dims(2)
338  do k = 3, dims(1)+2
339  velx_org(k,i,j) = ( momx_org(k,i-1,j) + momx_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
340  end do
341  end do
342  end do
343  do j = 1, dims(3)
344  do k = 3, dims(1)+2
345  velx_org(k,1,j) = momx_org(k,1,j) / dens_org(k,1,j)
346  end do
347  end do
348  velx_org(1:2,:,:) = 0.0_rp
349 
350  ! convert from momentum to velocity
351  do j = 2, dims(3)
352  do i = 1, dims(2)
353  do k = 3, dims(1)+2
354  vely_org(k,i,j) = ( momy_org(k,i,j-1) + momy_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
355  end do
356  end do
357  end do
358  do i = 1, dims(2)
359  do k = 3, dims(1)+2
360  vely_org(k,i,1) = momy_org(k,i,1) / dens_org(k,i,1)
361  end do
362  end do
363  vely_org(1:2,:,:) = 0.0_rp
364 
365 
366  !!! must be rotate !!!
367 
368 
369  do j = 1, dims(3)
370  do i = 1, dims(2)
371  do k = 3, dims(1)+2
372  ! diagnose temp and pres
373  call thermodyn_temp_pres( temp_org, & ! [OUT]
374  pres_org(k,i,j), & ! [OUT]
375  dens_org(k,i,j), & ! [IN]
376  rhot_org(k,i,j), & ! [IN]
377  qtrc_org(k,i,j,:), & ! [IN]
378  tracer_cv(:), & ! [IN]
379  tracer_r(:), & ! [IN]
380  tracer_mass(:) ) ! [IN]
381  pott_org(k,i,j) = rhot_org(k,i,j) / dens_org(k,i,j)
382  end do
383  dz = cz_org(3,i,j) - cz_org(2,i,j)
384  dens_org(2,i,j) = ( pres_org(3,i,j) + grav * dens_org(3,i,j) * dz * 0.5_rp ) &
385  / ( rdry * tsfc_org(i,j) - grav * dz * 0.5_rp )
386  pres_org(2,i,j) = dens_org(2,i,j) * rdry * tsfc_org(i,j)
387  pott_org(2,i,j) = tsfc_org(i,j) * ( p00 / pres_org(2,i,j) )**(rdry/cpdry)
388  temp_org = tsfc_org(i,j) + laps * cz_org(2,i,j)
389  pott_org(1,i,j) = temp_org * ( p00 / pres_org(1,i,j) )**(rdry/cpdry)
390  dens_org(1,i,j) = pres_org(1,i,j) / ( rdry * temp_org )
391  end do
392  end do
393 
394  return
395  end subroutine parentatominputscale
396 
397  !-----------------------------------------------------------------------------
399  subroutine parentlandsetupscale( &
400  ldims )
401  implicit none
402 
403  integer, intent(out) :: ldims(3)
404 
405  integer :: i
406  !---------------------------------------------------------------------------
407 
408  if( io_l ) write(io_fid_log,*) '+++ Real Case/Land Input File Type: SCALE-RM'
409  ldims(1) = parent_lkmax(handle)
410  ldims(2) = parent_imax(handle) * nest_tile_num_x
411  ldims(3) = parent_jmax(handle) * nest_tile_num_y
412 
413  if ( .not. allocated(read2d) ) then
414  allocate( read2d( parent_imax(handle), parent_jmax(handle) ) )
415  end if
416  allocate( read3dl( parent_imax(handle), parent_jmax(handle), ldims(1) ) )
417 
418  return
419  end subroutine parentlandsetupscale
420 
421  !-----------------------------------------------------------------------------
422  subroutine parentlandinputscale( &
423  tg_org, &
424  strg_org, &
425  lst_org, &
426  ust_org, &
427  albg_org, &
428  topo_org, &
429  lmask_org, &
430  llon_org, &
431  llat_org, &
432  lz_org, &
433  basename_land, &
434  ldims, &
435  use_file_landwater, &
436  it )
437  use gtool_file, only: &
438  fileread
439  use scale_const, only: &
440  d2r => const_d2r
441  implicit none
442 
443  real(RP), intent(out) :: tg_org(:,:,:)
444  real(RP), intent(out) :: strg_org(:,:,:)
445  real(RP), intent(out) :: lst_org(:,:)
446  real(RP), intent(out) :: ust_org(:,:)
447  real(RP), intent(out) :: albg_org(:,:,:)
448  real(RP), intent(out) :: topo_org(:,:)
449  real(RP), intent(out) :: lmask_org(:,:)
450  real(RP), intent(out) :: llon_org(:,:)
451  real(RP), intent(out) :: llat_org(:,:)
452  real(RP), intent(out) :: lz_org(:)
453 
454  character(len=*), intent(in) :: basename_land
455  integer, intent(in) :: ldims(3)
456  logical, intent(in) :: use_file_landwater ! use land water data from files
457  integer, intent(in) :: it
458 
459  integer :: rank
460 
461  integer :: k, i, j, n
462  integer :: xloc, yloc
463  integer :: xs, xe
464  integer :: ys, ye
465  !---------------------------------------------------------------------------
466 
467  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[LandInputSCALE]'
468 
469  do i = 1, size( nest_tile_id(:) )
470  ! read data from split files
471  rank = nest_tile_id(i)
472 
473  xloc = mod( i-1, nest_tile_num_x ) + 1
474  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
475 
476  xs = parent_imax(handle) * (xloc-1) + 1
477  xe = parent_imax(handle) * xloc
478  ys = parent_jmax(handle) * (yloc-1) + 1
479  ye = parent_jmax(handle) * yloc
480 
481  call fileread( read3dl(:,:,:), basename_land, "LAND_TEMP", it, rank )
482  do k = 1, ldims(1)
483  tg_org(k,xs:xe,ys:ye) = read3dl(:,:,k)
484  end do
485 
486  if( use_file_landwater )then
487  call fileread( read3dl(:,:,:), basename_land, "LAND_WATER", it, rank )
488  do k = 1, ldims(1)
489  strg_org(k,xs:xe,ys:ye) = read3dl(:,:,k)
490  end do
491  endif
492 
493  call fileread( read2d(:,:), basename_land, "lon", 1, rank )
494  llon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
495 
496  call fileread( read2d(:,:), basename_land, "lat", 1, rank )
497  llat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
498 
499  call fileread( read2d(:,:), basename_land, "LAND_SFC_TEMP", it, rank )
500  lst_org(xs:xe,ys:ye) = read2d(:,:)
501 
502  call fileread( read2d(:,:), basename_land, "URBAN_SFC_TEMP", it, rank )
503  ust_org(xs:xe,ys:ye) = read2d(:,:)
504 
505  call fileread( read2d(:,:), basename_land, "LAND_ALB_LW", it, rank )
506  albg_org(xs:xe,ys:ye,1) = read2d(:,:)
507 
508  call fileread( read2d(:,:), basename_land, "LAND_ALB_SW", it, rank )
509  albg_org(xs:xe,ys:ye,2) = read2d(:,:)
510 
511  call fileread( read2d(:,:), basename_land, "topo", it, rank )
512  topo_org(xs:xe,ys:ye) = read2d(:,:)
513 
514  call fileread( read2d(:,:), basename_land, "lsmask", it, rank )
515  lmask_org(xs:xe,ys:ye) = read2d(:,:)
516 
517  end do
518 
519  call fileread( lz_org(:), basename_land, "lz", 1, rank )
520 
521  return
522  end subroutine parentlandinputscale
523 
524  !-----------------------------------------------------------------------------
526  subroutine parentoceansetupscale( &
527  odims )
528  implicit none
529 
530  integer, intent(out) :: odims(2)
531 
532  integer :: i
533  !---------------------------------------------------------------------------
534 
535  if( io_l ) write(io_fid_log,*) '+++ Real Case/Ocean Input File Type: SCALE-RM'
536  odims(1) = parent_imax(handle) * nest_tile_num_x
537  odims(2) = parent_jmax(handle) * nest_tile_num_y
538 
539  if ( .not. allocated(read2d) ) then
540  allocate( read2d( parent_imax(handle), parent_jmax(handle) ) )
541  end if
542 
543  return
544  end subroutine parentoceansetupscale
545 
546  !-----------------------------------------------------------------------------
547  subroutine parentoceanopenscale( &
548  olon_org, &
549  olat_org, &
550  omask_org, &
551  basename_ocean, &
552  odims )
553  use scale_const, only: &
554  d2r => const_d2r
555  use gtool_file, only: &
556  fileread
557  implicit none
558  real(RP), intent(out) :: olon_org (:,:)
559  real(RP), intent(out) :: olat_org (:,:)
560  real(RP), intent(out) :: omask_org(:,:)
561  character(len=*), intent(in) :: basename_ocean
562  integer, intent(in) :: odims(2)
563 
564  integer :: rank
565  integer :: xloc, yloc
566  integer :: xs, xe, ys, ye
567 
568  integer :: i, k
569 
570  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanOpenSCALE]'
571 
572  do i = 1, size( nest_tile_id(:) )
573  ! read data from split files
574  rank = nest_tile_id(i)
575 
576  xloc = mod( i-1, nest_tile_num_x ) + 1
577  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
578 
579  xs = parent_imax(handle) * (xloc-1) + 1
580  xe = parent_imax(handle) * xloc
581  ys = parent_jmax(handle) * (yloc-1) + 1
582  ye = parent_jmax(handle) * yloc
583 
584  call fileread( read2d(:,:), basename_ocean, "lon", 1, rank )
585  olon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
586 
587  call fileread( read2d(:,:), basename_ocean, "lat", 1, rank )
588  olat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
589 
590  call fileread( read2d(:,:), basename_ocean, "lsmask", 1, rank )
591  omask_org(xs:xe,ys:ye) = read2d(:,:)
592 
593  end do
594 
595  return
596  end subroutine parentoceanopenscale
597 
598  !-----------------------------------------------------------------------------
599  subroutine parentoceaninputscale( &
600  tw_org, &
601  sst_org, &
602  albw_org, &
603  z0w_org, &
604  omask_org, &
605  basename_ocean, &
606  odims, &
607  it )
608  use gtool_file, only: &
609  fileread
610  implicit none
611 
612  real(RP), intent(out) :: tw_org(:,:)
613  real(RP), intent(out) :: sst_org(:,:)
614  real(RP), intent(out) :: albw_org(:,:,:)
615  real(RP), intent(out) :: z0w_org(:,:)
616  real(RP), intent(out) :: omask_org(:,:)
617 
618  character(len=*), intent(in) :: basename_ocean
619  integer, intent(in) :: odims(2)
620  integer, intent(in) :: it
621 
622  integer :: rank
623 
624  integer :: i, j, n
625  integer :: xloc, yloc
626  integer :: xs, xe
627  integer :: ys, ye
628  !---------------------------------------------------------------------------
629 
630  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanInputSCALE]'
631 
632  do i = 1, size( nest_tile_id(:) )
633  ! read data from split files
634  rank = nest_tile_id(i)
635 
636  xloc = mod( i-1, nest_tile_num_x ) + 1
637  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
638 
639  xs = parent_imax(handle) * (xloc-1) + 1
640  xe = parent_imax(handle) * xloc
641  ys = parent_jmax(handle) * (yloc-1) + 1
642  ye = parent_jmax(handle) * yloc
643 
644  call fileread( read2d(:,:), basename_ocean, "OCEAN_TEMP", it, rank )
645  tw_org(xs:xe,ys:ye) = read2d(:,:)
646 
647  call fileread( read2d(:,:), basename_ocean, "OCEAN_SFC_TEMP", it, rank )
648  sst_org(xs:xe,ys:ye) = read2d(:,:)
649 
650  call fileread( read2d(:,:), basename_ocean, "OCEAN_ALB_LW", it, rank )
651  albw_org(xs:xe,ys:ye,1) = read2d(:,:)
652 
653  call fileread( read2d(:,:), basename_ocean, "OCEAN_ALB_SW", it, rank )
654  albw_org(xs:xe,ys:ye,2) = read2d(:,:)
655 
656  call fileread( read2d(:,:), basename_ocean, "OCEAN_SFC_Z0M", it, rank )
657  z0w_org(xs:xe,ys:ye) = read2d(:,:)
658 
659  call fileread( read2d(:,:), basename_ocean, "lsmask", it, rank )
660  omask_org(xs:xe,ys:ye) = read2d(:,:)
661 
662  end do
663 
664  return
665  end subroutine parentoceaninputscale
666 
667 end module mod_realinput_scale
module GTOOL_FILE
Definition: gtool_file.f90:17
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
subroutine, public prc_mpistop
Abort MPI.
module ATMOSPHERE / Physics Cloud Microphysics - Convert
integer, dimension(2), public parent_jmax
parent max number in y-direction
module GRID (nesting system)
real(rp), dimension(qa_max), public tracer_r
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
subroutine, public parentoceanopenscale(olon_org, olat_org, omask_org, basename_ocean, odims)
module ATMOSPHERE / Physics Cloud Microphysics
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
subroutine, public atmos_phy_mp_bulk2bin(xs, xe, ys, ye, dims, it, rank, handle, basename_org, dens_org, qtrc_org)
Bulk to Bin.
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
procedure(intrpnest_intfc_interp_2d), pointer, public intrpnest_interp_2d
subroutine, public parentlandsetupscale(ldims)
Land Setup.
integer, dimension(:), allocatable, public nest_tile_id
parent tile real id
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:60
character(len=h_short), dimension(qa_max), public tracer_name
real(rp), dimension(qa_max), public tracer_cv
subroutine, public intrpnest_interp_fact_llz(hfact, vfact, kgrd, igrd, jgrd, ncopy, myhgt, mylat, mylon, myKS, myKE, myIA, myJA, inhgt, inlat, inlon, inKA, inIA, inJA, landgrid)
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
subroutine, public parentatomopenscale(lon_org, lat_org, cz_org, basename_org, dims)
module TRACER
module GRIDTRANS
subroutine, public parentoceaninputscale(tw_org, sst_org, albw_org, z0w_org, omask_org, basename_ocean, odims, it)
subroutine, public parentatomsetupscale(dims)
Atmos Setup.
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
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 COMMUNICATION
Definition: scale_comm.F90:23
module ATMOSPHERE / Hydrostatic barance
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
module PROCESS
real(rp), dimension(:,:,:), allocatable, public gtrans_rotc
rotation coefficient
integer, dimension(2), public parent_lkmax
parent max number in lz-direction
procedure(intrpnest_intfc_interp_3d), pointer, public intrpnest_interp_3d
integer, public nest_tile_num_y
parent tile number in y-direction
module CONSTANT
Definition: scale_const.F90:14
integer, dimension(2), public parent_kmax
parent max number in z-direction
integer, public prc_myrank
process num in local communicator
module INTERPOLATION (nesting system)
subroutine, public parentatominputscale(velz_org, velx_org, vely_org, pres_org, dens_org, pott_org, qtrc_org, cz_org, flg_bin, flg_intrp, basename_org, mptype_parent, dims, it)
integer, dimension(2), public parent_imax
parent max number in x-direction
module ATMOSPHERE / Thermodynamics
integer, public nest_tile_num_x
parent tile number in x-direction
subroutine, public parentoceansetupscale(odims)
Ocean Setup.
module PRECISION
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
module TOPOGRAPHY
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
real(rp), dimension(qa_max), public tracer_mass