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  flg_bin, &
164  flg_intrp, &
165  basename_org, &
166  mptype_parent, &
167  dims, &
168  it ) ! (in)
169  use scale_const, only: &
170  p00 => const_pre00, &
171  cpdry => const_cpdry, &
172  rdry => const_rdry
173  use scale_comm, only: &
174  comm_vars8, &
175  comm_wait
176  use gtool_file, only: &
177  fileread
178  use scale_atmos_hydrostatic, only: &
179  hydrostatic_buildrho_real => atmos_hydrostatic_buildrho_real
180  use scale_atmos_thermodyn, only: &
181  thermodyn_temp_pres => atmos_thermodyn_temp_pres, &
182  thermodyn_pott => atmos_thermodyn_pott
183  use scale_gridtrans, only: &
184  rotc => gtrans_rotc
185  use scale_topography, only: &
186  topo => topo_zsfc
187  use scale_interpolation_nest, only: &
191  use scale_atmos_phy_mp_convert, only: &
193  implicit none
194 
195  real(RP), intent(out) :: velz_org(:,:,:)
196  real(RP), intent(out) :: velx_org(:,:,:)
197  real(RP), intent(out) :: vely_org(:,:,:)
198  real(RP), intent(out) :: pres_org(:,:,:)
199  real(RP), intent(out) :: dens_org(:,:,:)
200  real(RP), intent(out) :: pott_org(:,:,:)
201  real(RP), intent(out) :: qtrc_org(:,:,:,:)
202  logical, intent(in) :: flg_bin ! flag for SBM(S10)
203  logical, intent(in) :: flg_intrp ! flag for interpolation of SBM
204  character(len=*), intent(in) :: basename_org
205  integer, intent(in) :: mptype_parent
206  integer, intent(in) :: dims(6)
207  integer, intent(in) :: it
208 
209  ! work
210 
211  real(RP) :: momz_org(dims(1)+2,dims(2),dims(3))
212  real(RP) :: momx_org(dims(1)+2,dims(2),dims(3))
213  real(RP) :: momy_org(dims(1)+2,dims(2),dims(3))
214  real(RP) :: rhot_org(dims(1)+2,dims(2),dims(3))
215  real(RP) :: tsfc_org( dims(2),dims(3))
216  real(RP) :: temp_org
217 
218  integer :: xs, xe
219  integer :: ys, ye
220  integer :: xloc, yloc
221  integer :: rank
222 
223  integer :: k, i, j, iq
224  logical :: lack_of_val
225  !---------------------------------------------------------------------------
226 
227 
228  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[AtomInputSCALE]'
229 
230  do i = 1, size( nest_tile_id(:) )
231  ! read data from split files
232  rank = nest_tile_id(i)
233 
234  xloc = mod( i-1, nest_tile_num_x ) + 1
235  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
236 
237  xs = parent_imax(handle) * (xloc-1) + 1
238  xe = parent_imax(handle) * xloc
239  ys = parent_jmax(handle) * (yloc-1) + 1
240  ye = parent_jmax(handle) * yloc
241 
242  call fileread( read2d(:,:), basename_org, "T2", it, rank )
243  tsfc_org(xs:xe,ys:ye) = read2d(:,:)
244 
245  call fileread( read2d(:,:), basename_org, "MSLP", it, rank )
246  pres_org(1,xs:xe,ys:ye) = read2d(:,:)
247 
248  call fileread( read3d(:,:,:), basename_org, "DENS", it, rank )
249  do k = 1, dims(1)
250  dens_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
251  end do
252 
253  call fileread( read3d(:,:,:), basename_org, "MOMZ", it, rank )
254  do k = 1, dims(1)
255  momz_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
256  end do
257 
258  call fileread( read3d(:,:,:), basename_org, "MOMX", it, rank )
259  do k = 1, dims(1)
260  momx_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
261  end do
262 
263  call fileread( read3d(:,:,:), basename_org, "MOMY", it, rank )
264  do k = 1, dims(1)
265  momy_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
266  end do
267 
268  call fileread( read3d(:,:,:), basename_org, "RHOT", it, rank )
269  do k = 1, dims(1)
270  rhot_org(k+2,xs:xe,ys:ye) = read3d(:,:,k)
271  end do
272 
273  if( flg_bin .and. flg_intrp ) then !--- tracers created from parent domain by interpolation
274 
275  if( io_l ) write(io_fid_log,*) '+++ SDF of SBM(S10) is interpolated from Qxx and Nxx'
276  if( io_l ) write(io_fid_log,*) '+++ of outer bulk MP model'
277 
278  call atmos_phy_mp_bulk2bin( xs, xe, & ! [IN]
279  ys, ye, & ! [IN]
280  dims, & ! [IN]
281  it, & ! [IN]
282  rank, & ! [IN]
283  handle, & ! [IN]
284  basename_org, & ! [IN]
285  dens_org, & ! [IN]
286  qtrc_org ) ! [INOUT]
287  do iq = 1, qa
288  qtrc_org(2,xs:xe,ys:ye,iq) = qtrc_org(3,xs:xe,ys:ye,iq)
289  enddo
290 
291  else !--- tracers of paremt domain directly used
292 
293  do iq = 1, mptype_parent
294  call fileread( read3d(:,:,:), basename_org, aq_name(iq), it, rank )
295  do k = 1, dims(1)
296  qtrc_org(k+2,xs:xe,ys:ye,iq) = read3d(:,:,k)
297  end do
298  qtrc_org(2,xs:xe,ys:ye,iq) = qtrc_org(3,xs:xe,ys:ye,iq)
299  end do
300  do iq = mptype_parent+1, qa
301  qtrc_org(:,xs:xe,ys:ye,iq) = 0.0_rp
302  end do
303 
304  endif
305 
306 ! call FileRead( read2D(:,:), BASENAME_ORG, "Q2", it, rank )
307 ! qtrc_org(2,xs:xe,ys:ye,I_QV) = read2D(:,:)
308 
309  end do
310 
311  ! convert from momentum to velocity
312  do j = 1, dims(3)
313  do i = 1, dims(2)
314  do k = 4, dims(1)+2
315  velz_org(k,i,j) = ( momz_org(k-1,i,j) + momz_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
316  end do
317  end do
318  end do
319  do j = 1, dims(3)
320  do i = 1, dims(2)
321  velz_org(1:3 ,i,j) = 0.0_rp
322  velz_org(dims(1)+2,i,j) = 0.0_rp
323  end do
324  end do
325 
326  ! convert from momentum to velocity
327  do j = 1, dims(3)
328  do i = 2, dims(2)
329  do k = 3, dims(1)+2
330  velx_org(k,i,j) = ( momx_org(k,i-1,j) + momx_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
331  end do
332  end do
333  end do
334  do j = 1, dims(3)
335  do k = 3, dims(1)+2
336  velx_org(k,1,j) = momx_org(k,1,j) / dens_org(k,1,j)
337  end do
338  end do
339  velx_org(1:2,:,:) = 0.0_rp
340 
341  ! convert from momentum to velocity
342  do j = 2, dims(3)
343  do i = 1, dims(2)
344  do k = 3, dims(1)+2
345  vely_org(k,i,j) = ( momy_org(k,i,j-1) + momy_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
346  end do
347  end do
348  end do
349  do i = 1, dims(2)
350  do k = 3, dims(1)+2
351  vely_org(k,i,1) = momy_org(k,i,1) / dens_org(k,i,1)
352  end do
353  end do
354  vely_org(1:2,:,:) = 0.0_rp
355 
356 
357  !!! must be rotate !!!
358 
359 
360  do j = 1, dims(3)
361  do i = 1, dims(2)
362  do k = 3, dims(1)+2
363  ! diagnose temp and pres
364  call thermodyn_temp_pres( temp_org, & ! [OUT]
365  pres_org(k,i,j), & ! [OUT]
366  dens_org(k,i,j), & ! [IN]
367  rhot_org(k,i,j), & ! [IN]
368  qtrc_org(k,i,j,:) ) ! [IN]
369  pott_org(k,i,j) = rhot_org(k,i,j) / dens_org(k,i,j)
370  end do
371  pott_org(1:2,i,j) = pott_org(3,i,j)
372  pres_org(2,i,j) = p00 * ( tsfc_org(i,j) / pott_org(2,i,j) )**(cpdry/rdry)
373  dens_org(1,i,j) = p00 * ( p00/pres_org(1,i,j) )**(rdry/cpdry-1.0_rp)
374  dens_org(2,i,j) = pres_org(2,i,j) / ( tsfc_org(i,j) * rdry )
375  qtrc_org(1,i,j,:) = qtrc_org(2,i,j,:)
376  end do
377  end do
378 
379  return
380  end subroutine parentatominputscale
381 
382  !-----------------------------------------------------------------------------
384  subroutine parentlandsetupscale( &
385  ldims )
386  implicit none
387 
388  integer, intent(out) :: ldims(3)
389 
390  integer :: i
391  !---------------------------------------------------------------------------
392 
393  if( io_l ) write(io_fid_log,*) '+++ Real Case/Land Input File Type: SCALE-RM'
394  ldims(1) = parent_lkmax(handle)
395  ldims(2) = parent_imax(handle) * nest_tile_num_x
396  ldims(3) = parent_jmax(handle) * nest_tile_num_y
397 
398  if ( .not. allocated(read2d) ) then
399  allocate( read2d( parent_imax(handle), parent_jmax(handle) ) )
400  end if
401  allocate( read3dl( parent_imax(handle), parent_jmax(handle), ldims(1) ) )
402 
403  return
404  end subroutine parentlandsetupscale
405 
406  !-----------------------------------------------------------------------------
407  subroutine parentlandinputscale( &
408  tg_org, &
409  strg_org, &
410  lst_org, &
411  ust_org, &
412  albg_org, &
413  topo_org, &
414  lmask_org, &
415  llon_org, &
416  llat_org, &
417  lz_org, &
418  basename_land, &
419  ldims, &
420  use_file_landwater, &
421  it )
422  use gtool_file, only: &
423  fileread
424  use scale_const, only: &
425  d2r => const_d2r
426  implicit none
427 
428  real(RP), intent(out) :: tg_org(:,:,:)
429  real(RP), intent(out) :: strg_org(:,:,:)
430  real(RP), intent(out) :: lst_org(:,:)
431  real(RP), intent(out) :: ust_org(:,:)
432  real(RP), intent(out) :: albg_org(:,:,:)
433  real(RP), intent(out) :: topo_org(:,:)
434  real(RP), intent(out) :: lmask_org(:,:)
435  real(RP), intent(out) :: llon_org(:,:)
436  real(RP), intent(out) :: llat_org(:,:)
437  real(RP), intent(out) :: lz_org(:)
438 
439  character(len=*), intent(in) :: basename_land
440  integer, intent(in) :: ldims(3)
441  logical, intent(in) :: use_file_landwater ! use land water data from files
442  integer, intent(in) :: it
443 
444  integer :: rank
445 
446  integer :: k, i, j, n
447  integer :: xloc, yloc
448  integer :: xs, xe
449  integer :: ys, ye
450  !---------------------------------------------------------------------------
451 
452  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[LandInputSCALE]'
453 
454  do i = 1, size( nest_tile_id(:) )
455  ! read data from split files
456  rank = nest_tile_id(i)
457 
458  xloc = mod( i-1, nest_tile_num_x ) + 1
459  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
460 
461  xs = parent_imax(handle) * (xloc-1) + 1
462  xe = parent_imax(handle) * xloc
463  ys = parent_jmax(handle) * (yloc-1) + 1
464  ye = parent_jmax(handle) * yloc
465 
466  call fileread( read3dl(:,:,:), basename_land, "LAND_TEMP", it, rank )
467  do k = 1, ldims(1)
468  tg_org(k,xs:xe,ys:ye) = read3dl(:,:,k)
469  end do
470 
471  if( use_file_landwater )then
472  call fileread( read3dl(:,:,:), basename_land, "LAND_WATER", it, rank )
473  do k = 1, ldims(1)
474  strg_org(k,xs:xe,ys:ye) = read3dl(:,:,k)
475  end do
476  endif
477 
478  call fileread( read2d(:,:), basename_land, "lon", 1, rank )
479  llon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
480 
481  call fileread( read2d(:,:), basename_land, "lat", 1, rank )
482  llat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
483 
484  call fileread( read2d(:,:), basename_land, "LAND_SFC_TEMP", it, rank )
485  lst_org(xs:xe,ys:ye) = read2d(:,:)
486 
487  call fileread( read2d(:,:), basename_land, "URBAN_SFC_TEMP", it, rank )
488  ust_org(xs:xe,ys:ye) = read2d(:,:)
489 
490  call fileread( read2d(:,:), basename_land, "LAND_ALB_LW", it, rank )
491  albg_org(xs:xe,ys:ye,1) = read2d(:,:)
492 
493  call fileread( read2d(:,:), basename_land, "LAND_ALB_SW", it, rank )
494  albg_org(xs:xe,ys:ye,2) = read2d(:,:)
495 
496  call fileread( read2d(:,:), basename_land, "topo", it, rank )
497  topo_org(xs:xe,ys:ye) = read2d(:,:)
498 
499  call fileread( read2d(:,:), basename_land, "lsmask", it, rank )
500  lmask_org(xs:xe,ys:ye) = read2d(:,:)
501 
502  end do
503 
504  call fileread( lz_org(:), basename_land, "lz", 1, rank )
505 
506  return
507  end subroutine parentlandinputscale
508 
509  !-----------------------------------------------------------------------------
511  subroutine parentoceansetupscale( &
512  odims )
513  implicit none
514 
515  integer, intent(out) :: odims(2)
516 
517  integer :: i
518  !---------------------------------------------------------------------------
519 
520  if( io_l ) write(io_fid_log,*) '+++ Real Case/Ocean Input File Type: SCALE-RM'
521  odims(1) = parent_imax(handle) * nest_tile_num_x
522  odims(2) = parent_jmax(handle) * nest_tile_num_y
523 
524  if ( .not. allocated(read2d) ) then
525  allocate( read2d( parent_imax(handle), parent_jmax(handle) ) )
526  end if
527 
528  return
529  end subroutine parentoceansetupscale
530 
531  !-----------------------------------------------------------------------------
532  subroutine parentoceanopenscale( &
533  olon_org, &
534  olat_org, &
535  omask_org, &
536  basename_ocean, &
537  odims )
538  use scale_const, only: &
539  d2r => const_d2r
540  use gtool_file, only: &
541  fileread
542  implicit none
543  real(RP), intent(out) :: olon_org (:,:)
544  real(RP), intent(out) :: olat_org (:,:)
545  real(RP), intent(out) :: omask_org(:,:)
546  character(len=*), intent(in) :: basename_ocean
547  integer, intent(in) :: odims(2)
548 
549  integer :: rank
550  integer :: xloc, yloc
551  integer :: xs, xe, ys, ye
552 
553  integer :: i, k
554 
555  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanOpenSCALE]'
556 
557  do i = 1, size( nest_tile_id(:) )
558  ! read data from split files
559  rank = nest_tile_id(i)
560 
561  xloc = mod( i-1, nest_tile_num_x ) + 1
562  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
563 
564  xs = parent_imax(handle) * (xloc-1) + 1
565  xe = parent_imax(handle) * xloc
566  ys = parent_jmax(handle) * (yloc-1) + 1
567  ye = parent_jmax(handle) * yloc
568 
569  call fileread( read2d(:,:), basename_ocean, "lon", 1, rank )
570  olon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
571 
572  call fileread( read2d(:,:), basename_ocean, "lat", 1, rank )
573  olat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
574 
575  call fileread( read2d(:,:), basename_ocean, "lsmask", 1, rank )
576  omask_org(xs:xe,ys:ye) = read2d(:,:)
577 
578  end do
579 
580  return
581  end subroutine parentoceanopenscale
582 
583  !-----------------------------------------------------------------------------
584  subroutine parentoceaninputscale( &
585  tw_org, &
586  sst_org, &
587  albw_org, &
588  z0w_org, &
589  omask_org, &
590  basename_ocean, &
591  odims, &
592  it )
593  use gtool_file, only: &
594  fileread
595  implicit none
596 
597  real(RP), intent(out) :: tw_org(:,:)
598  real(RP), intent(out) :: sst_org(:,:)
599  real(RP), intent(out) :: albw_org(:,:,:)
600  real(RP), intent(out) :: z0w_org(:,:)
601  real(RP), intent(out) :: omask_org(:,:)
602 
603  character(len=*), intent(in) :: basename_ocean
604  integer, intent(in) :: odims(2)
605  integer, intent(in) :: it
606 
607  integer :: rank
608 
609  integer :: i, j, n
610  integer :: xloc, yloc
611  integer :: xs, xe
612  integer :: ys, ye
613  !---------------------------------------------------------------------------
614 
615  if( io_l ) write(io_fid_log,*) '+++ ScaleLib/IO[realinput]/Categ[OceanInputSCALE]'
616 
617  do i = 1, size( nest_tile_id(:) )
618  ! read data from split files
619  rank = nest_tile_id(i)
620 
621  xloc = mod( i-1, nest_tile_num_x ) + 1
622  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
623 
624  xs = parent_imax(handle) * (xloc-1) + 1
625  xe = parent_imax(handle) * xloc
626  ys = parent_jmax(handle) * (yloc-1) + 1
627  ye = parent_jmax(handle) * yloc
628 
629  call fileread( read2d(:,:), basename_ocean, "OCEAN_TEMP", it, rank )
630  tw_org(xs:xe,ys:ye) = read2d(:,:)
631 
632  call fileread( read2d(:,:), basename_ocean, "OCEAN_SFC_TEMP", it, rank )
633  sst_org(xs:xe,ys:ye) = read2d(:,:)
634 
635  call fileread( read2d(:,:), basename_ocean, "OCEAN_ALB_LW", it, rank )
636  albw_org(xs:xe,ys:ye,1) = read2d(:,:)
637 
638  call fileread( read2d(:,:), basename_ocean, "OCEAN_ALB_SW", it, rank )
639  albw_org(xs:xe,ys:ye,2) = read2d(:,:)
640 
641  call fileread( read2d(:,:), basename_ocean, "OCEAN_SFC_Z0M", it, rank )
642  z0w_org(xs:xe,ys:ye) = read2d(:,:)
643 
644  call fileread( read2d(:,:), basename_ocean, "lsmask", it, rank )
645  omask_org(xs:xe,ys:ye) = read2d(:,:)
646 
647  end do
648 
649  return
650  end subroutine parentoceaninputscale
651 
652 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)
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine, public parentoceanopenscale(olon_org, olat_org, omask_org, basename_ocean, odims)
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
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:93
module REAL input SCALE
character(len=h_short), dimension(:), allocatable, public aq_name
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
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)
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
subroutine, public parentatominputscale(velz_org, velx_org, vely_org, pres_org, dens_org, pott_org, qtrc_org, flg_bin, flg_intrp, basename_org, mptype_parent, dims, it)