SCALE-RM
Functions/Subroutines
mod_realinput_scale Module Reference

module REAL input SCALE More...

Functions/Subroutines

subroutine, public parentatmossetupscale (dims)
 Atmos Setup. More...
 
subroutine, public parentatmosopenscale (lon_org, lat_org, cz_org, basename_org, dims)
 
subroutine, public parentatmosinputscale (velz_org, velx_org, vely_org, pres_org, dens_org, pott_org, qv_org, qtrc_org, cz_org, basename_org, sfc_diagnoses, same_mptype, dims, it)
 
subroutine, public parentlandsetupscale (ldims)
 Land Setup. More...
 
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)
 
subroutine, public parentoceansetupscale (odims)
 Ocean Setup. More...
 
subroutine, public parentoceanopenscale (olon_org, olat_org, omask_org, basename_ocean)
 
subroutine, public parentoceaninputscale (tw_org, sst_org, albw_org, z0w_org, omask_org, basename_ocean, it)
 

Detailed Description

module REAL input SCALE

Description
read data from SCALE file for real atmospheric simulations
Author
Team SCALE

Function/Subroutine Documentation

◆ parentatmossetupscale()

subroutine, public mod_realinput_scale::parentatmossetupscale ( integer, dimension(6), intent(out)  dims)

Atmos Setup.

Definition at line 64 of file mod_realinput_scale.F90.

64  use scale_comm_cartesc_nest, only: &
65  parent_kmax, &
66  parent_imax, &
67  parent_jmax, &
68  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
69  nest_tile_num_y => comm_cartesc_nest_tile_num_y
70  implicit none
71 
72  integer, intent(out) :: dims(6)
73  !---------------------------------------------------------------------------
74 
75  log_info("ParentAtmosSetupSCALE",*) 'Real Case/Atmos Input File Type: SCALE-RM'
76 
77  ! full level
78  dims(1) = parent_kmax(handle)
79  dims(2) = parent_imax(handle) * nest_tile_num_x
80  dims(3) = parent_jmax(handle) * nest_tile_num_y
81  ! half level
82  dims(4) = dims(1)
83  dims(5) = dims(2)
84  dims(6) = dims(3)
85 
86  allocate( read2d( parent_imax(handle),parent_jmax(handle)) )
87  allocate( read3d(dims(1),parent_imax(handle),parent_jmax(handle)) )
88 
89  allocate( rotc_cos(dims(2),dims(3)), rotc_sin(dims(2),dims(3)) )
90 
91  return

References scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_y, scale_comm_cartesc_nest::parent_imax, scale_comm_cartesc_nest::parent_jmax, and scale_comm_cartesc_nest::parent_kmax.

Referenced by mod_realinput::realinput_surface().

Here is the caller graph for this function:

◆ parentatmosopenscale()

subroutine, public mod_realinput_scale::parentatmosopenscale ( real(rp), dimension(:,:), intent(out)  lon_org,
real(rp), dimension(:,:), intent(out)  lat_org,
real(rp), dimension (:,:,:), intent(out)  cz_org,
character(len=*), intent(in)  basename_org,
integer, dimension(6), intent(in)  dims 
)

Definition at line 101 of file mod_realinput_scale.F90.

101  use scale_const, only: &
102  d2r => const_d2r
103  use scale_prc, only: &
104  prc_abort
105  use scale_comm_cartesc_nest, only: &
106  parent_imax, &
107  parent_jmax, &
108  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
109  nest_tile_id => comm_cartesc_nest_tile_id
110  use scale_file, only: &
111  file_open, &
112  file_get_attribute
113  use scale_file_cartesc, only: &
114  file_cartesc_read
115  use scale_mapprojection, only: &
116  mappinginfo, &
117  mappingparam, &
119  mapprojection_rotcoef
120  implicit none
121 
122  real(RP), intent(out) :: lon_org(:,:)
123  real(RP), intent(out) :: lat_org(:,:)
124  real(RP), intent(out) :: cz_org (:,:,:)
125  character(len=*), intent(in) :: basename_org
126  integer, intent(in) :: dims(6)
127 
128  ! mapping information
129  character(len=H_SHORT) :: vname
130  type(mappinginfo) :: mapping_info
131  type(mappingparam) :: mapping_param
132 
133  integer :: rank
134  integer :: xloc, yloc
135  integer :: xs, xe, ys, ye
136 
137  logical :: existed
138 
139  integer :: fid
140  integer :: i, j
141  !---------------------------------------------------------------------------
142 
143  do i = 1, size( nest_tile_id(:) )
144  ! read data from split files
145  rank = nest_tile_id(i)
146 
147  xloc = mod( i-1, nest_tile_num_x ) + 1
148  yloc = int( real(i-1) / real(nest_tile_num_x) ) + 1
149 
150  xs = parent_imax(handle) * (xloc-1) + 1
151  xe = parent_imax(handle) * xloc
152  ys = parent_jmax(handle) * (yloc-1) + 1
153  ye = parent_jmax(handle) * yloc
154 
155  call file_open( basename_org, & ! (in)
156  fid, & ! (out)
157  aggregate=.false., rankid=rank ) ! (in)
158 
159  call file_cartesc_read( fid, "lon", read2d(:,:) )
160  lon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
161 
162  call file_cartesc_read( fid, "lat", read2d(:,:) )
163  lat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
164 
165  call file_cartesc_read( fid, "height", read3d(:,:,:) )
166  cz_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
167 
168  call file_cartesc_read( fid, "topo", read2d(:,:), existed=existed )
169  cz_org(2,xs:xe,ys:ye) = read2d(:,:)
170 
171  if ( i == 1 ) then
172  call file_get_attribute( fid, "DENS", "grid_mapping", vname, existed=existed )
173  if ( existed ) then
174  call file_get_attribute( fid, vname, "grid_mapping_name", mapping_info%mapping_name )
175  call file_get_attribute( fid, vname, "false_easting", mapping_info%false_easting, existed=existed )
176  call file_get_attribute( fid, vname, "false_northing", mapping_info%false_northing, existed=existed )
177  call file_get_attribute( fid, vname, "longitude_of_central_meridian", mapping_info%longitude_of_central_meridian, existed=existed )
178  call file_get_attribute( fid, vname, "longitude_of_projection_origin", mapping_info%longitude_of_projection_origin, existed=existed )
179  call file_get_attribute( fid, vname, "latitude_of_projection_origin", mapping_info%latitude_of_projection_origin, existed=existed )
180  call file_get_attribute( fid, vname, "straight_vertical_longitude_from_pole", mapping_info%straight_vertical_longitude_from_pole, existed=existed )
181  call file_get_attribute( fid, vname, "standard_parallel", mapping_info%standard_parallel(:), existed=existed )
182  if ( .not. existed ) then
183  call file_get_attribute( fid, vname, "standard_parallel", mapping_info%standard_parallel(1), existed=existed )
184  end if
185  call file_get_attribute( fid, vname, "rotation", mapping_info%rotation, existed=existed )
186  if ( .not. existed ) mapping_info%rotation = 0.0_dp
187  else
188  mapping_info%mapping_name = ""
189  end if
190  end if
191 
192  end do
193 
194  cz_org(1,:,:) = 0.0_rp
195 
196  if ( mapping_info%mapping_name .ne. "" ) then
197  call mapprojection_get_param( mapping_info, & ! (in)
198  mapping_param ) ! (out)
199  call mapprojection_rotcoef( dims(2), 1, dims(2), dims(3), 1, dims(3), &
200  lon_org(:,:), lat_org(:,:), & ! (in)
201  mapping_info%mapping_name, & ! (in)
202  mapping_param, & ! (in)
203  rotc_cos(:,:), rotc_sin(:,:) ) ! (out)
204  else
205 !OCL XFILL
206  !$omp parallel do
207  do j = 1, dims(3)
208  do i = 1, dims(2)
209  rotc_cos(i,j) = 1.0_rp
210  rotc_sin(i,j) = 0.0_rp
211  end do
212  end do
213  end if
214 
215  return

References scale_comm_cartesc_nest::comm_cartesc_nest_tile_id, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_const::const_d2r, scale_file::file_open(), scale_mapprojection::mapprojection_get_param(), scale_comm_cartesc_nest::parent_imax, scale_comm_cartesc_nest::parent_jmax, and scale_prc::prc_abort().

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentatmosinputscale()

subroutine, public mod_realinput_scale::parentatmosinputscale ( real(rp), dimension(:,:,:), intent(out)  velz_org,
real(rp), dimension(:,:,:), intent(out)  velx_org,
real(rp), dimension(:,:,:), intent(out)  vely_org,
real(rp), dimension(:,:,:), intent(out)  pres_org,
real(rp), dimension(:,:,:), intent(out)  dens_org,
real(rp), dimension(:,:,:), intent(out)  pott_org,
real(rp), dimension (:,:,:), intent(out)  qv_org,
real(rp), dimension(:,:,:,:), intent(out)  qtrc_org,
real(rp), dimension (:,:,:), intent(in)  cz_org,
character(len=*), intent(in)  basename_org,
logical, intent(in)  sfc_diagnoses,
logical, intent(in)  same_mptype,
integer, dimension(6), intent(in)  dims,
integer, intent(in)  it 
)

Definition at line 234 of file mod_realinput_scale.F90.

234  use scale_const, only: &
235  undef => const_undef, &
236  p00 => const_pre00, &
237  cpdry => const_cpdry, &
238  rdry => const_rdry, &
239  grav => const_grav, &
240  laps => const_laps
241  use scale_comm_cartesc, only: &
242  comm_vars8, &
243  comm_wait
244  use scale_comm_cartesc_nest, only: &
245  parent_imax, &
246  parent_jmax, &
247  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
248  nest_tile_id => comm_cartesc_nest_tile_id
249  use scale_file, only: &
250  file_open
251  use scale_file_cartesc, only: &
252  file_cartesc_read
255  use scale_atmos_thermodyn, only: &
256  thermodyn_specific_heat => atmos_thermodyn_specific_heat, &
257  thermodyn_rhot2temp_pres => atmos_thermodyn_rhot2temp_pres
258  use scale_atmos_hydrometeor, only: &
259  n_hyd, &
260  hyd_name, &
261  num_name
262  use mod_atmos_phy_mp_driver, only: &
264  use mod_atmos_phy_mp_vars, only: &
265  qa_mp, &
266  qs_mp, &
267  qe_mp
268  implicit none
269 
270  real(RP), intent(out) :: velz_org(:,:,:)
271  real(RP), intent(out) :: velx_org(:,:,:)
272  real(RP), intent(out) :: vely_org(:,:,:)
273  real(RP), intent(out) :: pres_org(:,:,:)
274  real(RP), intent(out) :: dens_org(:,:,:)
275  real(RP), intent(out) :: pott_org(:,:,:)
276  real(RP), intent(out) :: qv_org (:,:,:)
277  real(RP), intent(out) :: qtrc_org(:,:,:,:)
278  real(RP), intent(in) :: cz_org (:,:,:)
279  character(len=*), intent(in) :: basename_org
280  logical, intent(in) :: sfc_diagnoses
281  logical, intent(in) :: same_mptype
282  integer, intent(in) :: dims(6)
283  integer, intent(in) :: it
284 
285  real(RP) :: momz_org(dims(1)+2,dims(2),dims(3))
286  real(RP) :: momx_org(dims(1)+2,dims(2),dims(3))
287  real(RP) :: momy_org(dims(1)+2,dims(2),dims(3))
288  real(RP) :: rhot_org(dims(1)+2,dims(2),dims(3))
289  real(RP) :: tsfc_org( dims(2),dims(3))
290  real(RP) :: qhyd_org(dims(1)+2,dims(2),dims(3),N_HYD)
291  real(RP) :: qnum_org(dims(1)+2,dims(2),dims(3),N_HYD)
292  real(RP) :: qdry (dims(1)+2)
293  real(RP) :: Rtot (dims(1)+2)
294  real(RP) :: CVtot (dims(1)+2)
295  real(RP) :: CPtot (dims(1)+2)
296  real(RP) :: temp_org(dims(1)+2)
297  real(RP) :: dz
298 
299  real(RP) :: u, v
300 
301  integer :: xs, xe
302  integer :: ys, ye
303  integer :: nx, ny
304  integer :: xloc, yloc
305  integer :: rank
306 
307  logical :: qnum_flag
308 
309  integer :: fid
310  logical :: existed, existed_t2, existed_mslp
311  integer :: k, i, j, iq
312  !---------------------------------------------------------------------------
313 
314  qnum_flag = .false.
315 
316  do i = 1, size( nest_tile_id(:) )
317  ! read data from split files
318  rank = nest_tile_id(i)
319 
320  xloc = mod( i-1, nest_tile_num_x ) + 1
321  yloc = int( real(i-1) / real(nest_tile_num_x) ) + 1
322 
323  xs = parent_imax(handle) * (xloc-1) + 1
324  xe = parent_imax(handle) * xloc
325  ys = parent_jmax(handle) * (yloc-1) + 1
326  ye = parent_jmax(handle) * yloc
327 
328  nx = xe - xs + 1
329  ny = ye - ys + 1
330 
331  call file_open( basename_org, & ! (in)
332  fid, & ! (out)
333  aggregate=.false., rankid=rank ) ! (in)
334 
335  if ( sfc_diagnoses ) then
336  call file_cartesc_read( fid, "T2", read2d(:,:), step=it, existed=existed_t2 )
337  if ( existed_t2 ) then
338 !OCL XFILL
339  tsfc_org(xs:xe,ys:ye) = read2d(:,:)
340  end if
341 
342  call file_cartesc_read( fid, "MSLP", read2d(:,:), step=it, existed=existed_mslp )
343  if ( existed_mslp ) then
344 !OCL XFILL
345  pres_org(1,xs:xe,ys:ye) = read2d(:,:)
346  end if
347  end if
348 
349  call file_cartesc_read( fid, "DENS", read3d(:,:,:), step=it )
350 !OCL XFILL
351  dens_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
352 
353  call file_cartesc_read( fid, "MOMZ", read3d(:,:,:), step=it )
354 !OCL XFILL
355  momz_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
356 
357  call file_cartesc_read( fid, "MOMX", read3d(:,:,:), step=it )
358 !OCL XFILL
359  momx_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
360 
361  call file_cartesc_read( fid, "MOMY", read3d(:,:,:), step=it )
362 !OCL XFILL
363  momy_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
364 
365 
366  call file_cartesc_read( fid, "RHOT", read3d(:,:,:), step=it )
367 !OCL XFILL
368  rhot_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
369 
370 
371  if ( same_mptype ) then
372 
373  do iq = qs_mp, qe_mp
374  call file_cartesc_read( fid, tracer_name(iq), read3d(:,:,:), step=it )
375 !OCL XFILL
376  qtrc_org(3:dims(1)+2,xs:xe,ys:ye,iq) = read3d(:,:,:)
377 !OCL XFILL
378  qtrc_org(2,xs:xe,ys:ye,iq) = qtrc_org(3,xs:xe,ys:ye,iq)
379 !OCL XFILL
380  qtrc_org(1,xs:xe,ys:ye,iq) = qtrc_org(3,xs:xe,ys:ye,iq)
381  enddo
382 
383  else
384 
385 !OCL XFILL
386  do iq = 1, n_hyd
387  qhyd_org(:,xs:xe,ys:ye,iq) = 0.0_rp
388  qnum_org(:,xs:xe,ys:ye,iq) = 0.0_rp
389  end do
390 
391  call file_cartesc_read( fid, "QV", read3d(:,:,:), step=it, existed=existed )
392  if ( existed ) then
393 !OCL XFILL
394  qv_org(3:dims(1)+2,xs:xe,ys:ye) = read3d(:,:,:)
395 !OCL XFILL
396  qv_org(2,xs:xe,ys:ye) = qv_org(3,xs:xe,ys:ye)
397 !OCL XFILL
398  qv_org(1,xs:xe,ys:ye) = qv_org(3,xs:xe,ys:ye)
399  else
400 !OCL XFILL
401  qv_org(:,:,:) = 0.0_rp
402  end if
403 
404  do iq = 1, n_hyd
405  call file_cartesc_read( fid, hyd_name(iq), read3d(:,:,:), step=it, existed=existed )
406  if ( existed ) then
407 !OCL XFILL
408  qhyd_org(3:dims(1)+2,xs:xe,ys:ye,iq) = read3d(:,:,:)
409 !OCL XFILL
410  qhyd_org(2,xs:xe,ys:ye,iq) = qhyd_org(3,xs:xe,ys:ye,iq)
411 !OCL XFILL
412  qhyd_org(1,xs:xe,ys:ye,iq) = qhyd_org(3,xs:xe,ys:ye,iq)
413  else
414 !OCL XFILL
415  qhyd_org(:,:,:,iq) = 0.0_rp
416  end if
417 
418  ! number density
419  call file_cartesc_read( fid, num_name(iq), read3d(:,:,:), step=it, existed=existed )
420  if ( existed ) then
421 !OCL XFILL
422  qnum_org(3:dims(1)+2,xs:xe,ys:ye,iq) = read3d(:,:,:)
423 !OCL XFILL
424  qnum_org(2,xs:xe,ys:ye,iq) = qnum_org(3,xs:xe,ys:ye,iq)
425 !OCL XFILL
426  qnum_org(1,xs:xe,ys:ye,iq) = qnum_org(3,xs:xe,ys:ye,iq)
427  qnum_flag = .true.
428  else
429 !OCL XFILL
430  qnum_org(:,:,:,iq) = undef
431  end if
432  end do
433 
434  endif
435 
436 ! call FILE_CARTESC_read( fid, "Q2", read2D(:,:), step=it )
437 ! qv_org(2,xs:xe,ys:ye) = read2D(:,:)
438 
439  end do
440 
441  ! convert from momentum to velocity
442 !OCL XFILL
443  do j = 1, dims(3)
444  do i = 1, dims(2)
445  do k = 4, dims(1)+2
446  velz_org(k,i,j) = ( momz_org(k-1,i,j) + momz_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
447  end do
448  end do
449  end do
450 !OCL XFILL
451  do j = 1, dims(3)
452  do i = 1, dims(2)
453  velz_org(1:3 ,i,j) = 0.0_rp
454  velz_org(dims(1)+2,i,j) = 0.0_rp
455  end do
456  end do
457 
458  ! convert from momentum to velocity
459 !OCL XFILL
460  do j = 1, dims(3)
461  do i = 2, dims(2)
462  do k = 3, dims(1)+2
463  velx_org(k,i,j) = ( momx_org(k,i-1,j) + momx_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
464  end do
465  end do
466  end do
467 !OCL XFILL
468  do j = 1, dims(3)
469  do k = 3, dims(1)+2
470  velx_org(k,1,j) = momx_org(k,1,j) / dens_org(k,1,j)
471  end do
472  end do
473  if ( sfc_diagnoses ) then
474 !OCL XFILL
475  velx_org(1:2,:,:) = 0.0_rp
476  else
477 !OCL XFILL
478  velx_org(1:2,:,:) = undef
479  end if
480 
481  ! convert from momentum to velocity
482 !OCL XFILL
483  do j = 2, dims(3)
484  do i = 1, dims(2)
485  do k = 3, dims(1)+2
486  vely_org(k,i,j) = ( momy_org(k,i,j-1) + momy_org(k,i,j) ) / dens_org(k,i,j) * 0.5_rp
487  end do
488  end do
489  end do
490 !OCL XFILL
491  do i = 1, dims(2)
492  do k = 3, dims(1)+2
493  vely_org(k,i,1) = momy_org(k,i,1) / dens_org(k,i,1)
494  end do
495  end do
496  if ( sfc_diagnoses ) then
497 !OCL XFILL
498  vely_org(1:2,:,:) = 0.0_rp
499  else
500 !OCL XFILL
501  vely_org(1:2,:,:) = undef
502  end if
503 
504 
505  ! rotation
506  !$omp parallel do &
507  !$omp private(u,v)
508  do j = 1, dims(3)
509  do i = 1, dims(2)
510  do k = 1, dims(1)+2
511  u = velx_org(k,i,j)
512  v = vely_org(k,i,j)
513  if ( u .ne. undef .and. v .ne.undef ) then
514  velx_org(k,i,j) = u * rotc_cos(i,j) - v * rotc_sin(i,j)
515  vely_org(k,i,j) = u * rotc_sin(i,j) + v * rotc_cos(i,j)
516  else
517  velx_org(k,i,j) = undef
518  vely_org(k,i,j) = undef
519  end if
520  end do
521  end do
522  end do
523 
524 
525  do iq = 1, size(qtrc_org,4)
526  if ( iq >= qs_mp .and. iq <= qe_mp ) cycle
527 !OCL XFILL
528  qtrc_org(:,:,:,iq) = 0.0_rp
529  end do
530 
531  if ( qa_mp > 0 .AND. .NOT. same_mptype ) then
532  if ( qnum_flag ) then
533  call atmos_phy_mp_driver_qhyd2qtrc( dims(1)+2, 1, dims(1)+2, dims(2), 1, dims(2), dims(3), 1, dims(3), &
534  qv_org(:,:,:), qhyd_org(:,:,:,:), & ! [IN]
535  qtrc_org(:,:,:,qs_mp:qe_mp), & ! [OUT]
536  qnum=qnum_org(:,:,:,:) ) ! [IN]
537  else
538  call atmos_phy_mp_driver_qhyd2qtrc( dims(1)+2, 1, dims(1)+2, dims(2), 1, dims(2), dims(3), 1, dims(3), &
539  qv_org(:,:,:), qhyd_org(:,:,:,:), & ! [IN]
540  qtrc_org(:,:,:,qs_mp:qe_mp) ) ! [OUT]
541  end if
542  end if
543 
544  !$omp parallel do default(none) &
545  !$omp private(dz,temp_org,qdry,rtot,cvtot,cptot) &
546  !$omp shared(dims,QA,UNDEF,LAPS,P00,Rdry,CPdry,GRAV,TRACER_MASS,TRACER_R,TRACER_CV,TRACER_CP, &
547  !$omp dens_org,rhot_org,qtrc_org,pres_org,pott_org,tsfc_org,cz_org, &
548  !$omp existed_t2,existed_mslp,sfc_diagnoses)
549  do j = 1, dims(3)
550  do i = 1, dims(2)
551 
552  ! diagnose temp and pres
553  call thermodyn_specific_heat( dims(1)+2, 3, dims(1)+2, qa, & ! [IN]
554  qtrc_org(:,i,j,:), & ! [IN]
555  tracer_mass(:), tracer_r(:), tracer_cv(:), tracer_cp(:), & ! [IN]
556  qdry(:), rtot(:), cvtot(:), cptot(:) ) ! [OUT]
557 
558  call thermodyn_rhot2temp_pres( dims(1), 1, dims(1), & ! [IN]
559  dens_org(3:dims(1)+2,i,j), rhot_org(3:dims(1)+2,i,j), & ! [IN]
560  rtot(3:dims(1)+2), cvtot(3:dims(1)+2), cptot(3:dims(1)+2), & ! [IN]
561  temp_org(3:dims(1)+2), pres_org(3:dims(1)+2,i,j) ) ! [OUT]
562 
563 !OCL XFILL
564  do k = 3, dims(1)+2
565  pott_org(k,i,j) = rhot_org(k,i,j) / dens_org(k,i,j)
566  end do
567 
568  if ( sfc_diagnoses ) then
569  ! at the surface
570  dz = cz_org(3,i,j) - cz_org(2,i,j)
571 
572  if ( .not. existed_t2 ) then
573  tsfc_org(i,j) = temp_org(3) + laps * dz
574  end if
575  dens_org(2,i,j) = ( pres_org(3,i,j) + grav * dens_org(3,i,j) * dz * 0.5_rp ) &
576  / ( rdry * tsfc_org(i,j) - grav * dz * 0.5_rp )
577  pres_org(2,i,j) = dens_org(2,i,j) * rdry * tsfc_org(i,j)
578  pott_org(2,i,j) = tsfc_org(i,j) * ( p00 / pres_org(2,i,j) )**(rdry/cpdry)
579 
580  ! at the sea-level
581  temp_org(1) = tsfc_org(i,j) + laps * cz_org(2,i,j)
582  if ( existed_mslp ) then
583  pott_org(1,i,j) = temp_org(1) * ( p00 / pres_org(1,i,j) )**(rdry/cpdry)
584  dens_org(1,i,j) = pres_org(1,i,j) / ( rdry * temp_org(1) )
585  else
586  dens_org(1,i,j) = ( pres_org(2,i,j) + grav * dens_org(2,i,j) * cz_org(2,i,j) * 0.5_rp ) &
587  / ( rdry * temp_org(1) - grav * cz_org(2,i,j) * 0.5_rp )
588  pres_org(1,i,j) = dens_org(1,i,j) * rdry * temp_org(1)
589  pott_org(1,i,j) = temp_org(1) * ( p00 / pres_org(1,i,j) )**(rdry/cpdry)
590  end if
591  else
592  dens_org(1:2,i,j) = undef
593  pres_org(1:2,i,j) = undef
594  pott_org(1:2,i,j) = undef
595  end if
596 
597  end do
598  end do
599 
600  return

References scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc, mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc(), scale_comm_cartesc_nest::comm_cartesc_nest_tile_id, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_const::const_cpdry, scale_const::const_grav, scale_const::const_laps, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_undef, scale_file::file_open(), scale_atmos_hydrometeor::hyd_name, scale_tracer::k, scale_atmos_hydrometeor::n_hyd, scale_atmos_hydrometeor::num_name, scale_comm_cartesc_nest::parent_imax, scale_comm_cartesc_nest::parent_jmax, scale_tracer::qa, mod_atmos_phy_mp_vars::qa_mp, mod_atmos_phy_mp_vars::qe_mp, mod_atmos_phy_mp_vars::qs_mp, scale_tracer::tracer_cp, scale_tracer::tracer_cv, scale_tracer::tracer_mass, scale_tracer::tracer_name, and scale_tracer::tracer_r.

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentlandsetupscale()

subroutine, public mod_realinput_scale::parentlandsetupscale ( integer, dimension(3), intent(out)  ldims)

Land Setup.

Definition at line 607 of file mod_realinput_scale.F90.

607  use scale_comm_cartesc_nest, only: &
608  parent_imax, &
609  parent_jmax, &
610  parent_lkmax, &
611  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
612  nest_tile_num_y => comm_cartesc_nest_tile_num_y
613  implicit none
614 
615  integer, intent(out) :: ldims(3)
616  !---------------------------------------------------------------------------
617 
618  log_info("ParentLandSetupSCALE",*) 'Real Case/Land Input File Type: SCALE-RM'
619 
620  ldims(1) = parent_lkmax(handle)
621  ldims(2) = parent_imax(handle) * nest_tile_num_x
622  ldims(3) = parent_jmax(handle) * nest_tile_num_y
623 
624  if ( .NOT. allocated(read2d) ) then
625  allocate( read2d(parent_imax(handle),parent_jmax(handle)) )
626  endif
627 
628  allocate( read3dl(ldims(1),parent_imax(handle),parent_jmax(handle)) )
629 
630  return

References scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_y, scale_comm_cartesc_nest::parent_imax, scale_comm_cartesc_nest::parent_jmax, and scale_comm_cartesc_nest::parent_lkmax.

Referenced by mod_realinput::realinput_surface().

Here is the caller graph for this function:

◆ parentlandinputscale()

subroutine, public mod_realinput_scale::parentlandinputscale ( real(rp), dimension (:,:,:), intent(out)  tg_org,
real(rp), dimension (:,:,:), intent(out)  strg_org,
real(rp), dimension (:,:), intent(out)  lst_org,
real(rp), dimension (:,:), intent(out)  ust_org,
real(rp), dimension (:,:,:,:), intent(out)  albg_org,
real(rp), dimension (:,:), intent(out)  topo_org,
real(rp), dimension(:,:), intent(out)  lmask_org,
real(rp), dimension (:,:), intent(out)  llon_org,
real(rp), dimension (:,:), intent(out)  llat_org,
real(rp), dimension (:), intent(out)  lz_org,
character(len=*), intent(in)  basename_land,
integer, dimension(3), intent(in)  ldims,
logical, intent(in)  use_file_landwater,
integer, intent(in)  it 
)

Definition at line 649 of file mod_realinput_scale.F90.

649  use scale_const, only: &
650  undef => const_undef, &
651  d2r => const_d2r
652  use scale_comm_cartesc_nest, only: &
653  parent_imax, &
654  parent_jmax, &
655  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
656  nest_tile_id => comm_cartesc_nest_tile_id
657  use scale_file, only: &
658  file_open, &
659  file_read, &
660  file_get_datainfo
661  use scale_file_cartesc, only: &
662  file_cartesc_read
663  implicit none
664 
665  real(RP), intent(out) :: tg_org (:,:,:)
666  real(RP), intent(out) :: strg_org (:,:,:)
667  real(RP), intent(out) :: lst_org (:,:)
668  real(RP), intent(out) :: ust_org (:,:)
669  real(RP), intent(out) :: albg_org (:,:,:,:)
670  real(RP), intent(out) :: topo_org (:,:)
671  real(RP), intent(out) :: lmask_org(:,:)
672  real(RP), intent(out) :: llon_org (:,:)
673  real(RP), intent(out) :: llat_org (:,:)
674  real(RP), intent(out) :: lz_org (:)
675  character(len=*), intent(in) :: basename_land
676  integer, intent(in) :: ldims(3)
677  logical, intent(in) :: use_file_landwater ! use land water data from files
678  integer, intent(in) :: it
679 
680  integer :: rank
681  integer :: fid
682  integer :: xloc, yloc
683  integer :: xs, xe
684  integer :: ys, ye
685  logical :: existed
686 
687  integer :: i
688  !---------------------------------------------------------------------------
689 
690  do i = 1, size( nest_tile_id(:) )
691  ! read data from split files
692  rank = nest_tile_id(i)
693 
694  xloc = mod( i-1, nest_tile_num_x ) + 1
695  yloc = int( real(i-1) / real(nest_tile_num_x) ) + 1
696 
697  xs = parent_imax(handle) * (xloc-1) + 1
698  xe = parent_imax(handle) * xloc
699  ys = parent_jmax(handle) * (yloc-1) + 1
700  ye = parent_jmax(handle) * yloc
701 
702  call file_open( basename_land, & ! (in)
703  fid, & ! (out)
704  aggregate=.false., rankid=rank ) ! (in)
705 
706  call file_cartesc_read( fid, "LAND_TEMP", read3dl(:,:,:), step=it )
707  tg_org(1:ldims(1),xs:xe,ys:ye) = read3dl(:,:,:)
708 
709  if( use_file_landwater )then
710  call file_cartesc_read( fid, "LAND_WATER", read3dl(:,:,:), step=it )
711  strg_org(1:ldims(1),xs:xe,ys:ye) = read3dl(:,:,:)
712  endif
713 
714  call file_cartesc_read( fid, "lon", read2d(:,:) )
715  llon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
716 
717  call file_cartesc_read( fid, "lat", read2d(:,:) )
718  llat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
719 
720  call file_cartesc_read( fid, "LAND_SFC_TEMP", read2d(:,:), step=it )
721  lst_org(xs:xe,ys:ye) = read2d(:,:)
722 
723  call file_get_datainfo( fid, "URBAN_SFC_TEMP", istep=it, existed=existed )
724  if ( existed ) then
725  call file_cartesc_read( fid, "URBAN_SFC_TEMP", read2d(:,:), step=it )
726  ust_org(xs:xe,ys:ye) = read2d(:,:)
727  else
728  ust_org(xs:xe,ys:ye) = undef
729  end if
730 
731  call file_cartesc_read( fid, "LAND_SFC_ALB_IR_dir", read2d(:,:), step=it )
732  albg_org(xs:xe,ys:ye,i_r_direct ,i_r_ir ) = read2d(:,:)
733 
734  call file_cartesc_read( fid, "LAND_SFC_ALB_IR_dif", read2d(:,:), step=it )
735  albg_org(xs:xe,ys:ye,i_r_diffuse,i_r_ir ) = read2d(:,:)
736 
737  call file_cartesc_read( fid, "LAND_SFC_ALB_NIR_dir", read2d(:,:), step=it )
738  albg_org(xs:xe,ys:ye,i_r_direct ,i_r_nir) = read2d(:,:)
739 
740  call file_cartesc_read( fid, "LAND_SFC_ALB_NIR_dif", read2d(:,:), step=it )
741  albg_org(xs:xe,ys:ye,i_r_diffuse,i_r_nir) = read2d(:,:)
742 
743  call file_cartesc_read( fid, "LAND_SFC_ALB_VIS_dir", read2d(:,:), step=it )
744  albg_org(xs:xe,ys:ye,i_r_direct ,i_r_vis) = read2d(:,:)
745 
746  call file_cartesc_read( fid, "LAND_SFC_ALB_VIS_dif", read2d(:,:), step=it )
747  albg_org(xs:xe,ys:ye,i_r_diffuse,i_r_vis) = read2d(:,:)
748 
749  call file_cartesc_read( fid, "topo", read2d(:,:) )
750  topo_org(xs:xe,ys:ye) = read2d(:,:)
751 
752  call file_cartesc_read( fid, "lsmask", read2d(:,:) )
753  lmask_org(xs:xe,ys:ye) = read2d(:,:)
754 
755  end do
756 
757  call file_read( fid, "lz", lz_org(:) )
758 
759  return

References scale_comm_cartesc_nest::comm_cartesc_nest_tile_id, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_const::const_d2r, scale_const::const_undef, scale_file::file_open(), scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_comm_cartesc_nest::parent_imax, and scale_comm_cartesc_nest::parent_jmax.

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentoceansetupscale()

subroutine, public mod_realinput_scale::parentoceansetupscale ( integer, dimension(2), intent(out)  odims)

Ocean Setup.

Definition at line 766 of file mod_realinput_scale.F90.

766  use scale_comm_cartesc_nest, only: &
767  parent_imax, &
768  parent_jmax, &
769  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
770  nest_tile_num_y => comm_cartesc_nest_tile_num_y
771  implicit none
772 
773  integer, intent(out) :: odims(2)
774  !---------------------------------------------------------------------------
775 
776  log_info("ParentOceanSetupSCALE",*) 'Real Case/Ocean Input File Type: SCALE-RM'
777 
778  odims(1) = parent_imax(handle) * nest_tile_num_x
779  odims(2) = parent_jmax(handle) * nest_tile_num_y
780 
781  if ( .NOT. allocated(read2d) ) then
782  allocate( read2d(parent_imax(handle),parent_jmax(handle)) )
783  end if
784 
785  allocate( read3do(1,parent_imax(handle),parent_jmax(handle)) )
786 
787  return

References scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_y, scale_comm_cartesc_nest::parent_imax, and scale_comm_cartesc_nest::parent_jmax.

Referenced by mod_realinput::realinput_surface().

Here is the caller graph for this function:

◆ parentoceanopenscale()

subroutine, public mod_realinput_scale::parentoceanopenscale ( real(rp), dimension (:,:), intent(out)  olon_org,
real(rp), dimension (:,:), intent(out)  olat_org,
real(rp), dimension(:,:), intent(out)  omask_org,
character(len=*), intent(in)  basename_ocean 
)

Definition at line 796 of file mod_realinput_scale.F90.

796  use scale_const, only: &
797  d2r => const_d2r
798  use scale_comm_cartesc_nest, only: &
799  parent_imax, &
800  parent_jmax, &
801  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
802  nest_tile_id => comm_cartesc_nest_tile_id
803  use scale_file, only: &
804  file_open
805  use scale_file_cartesc, only: &
806  file_cartesc_read
807  implicit none
808 
809  real(RP), intent(out) :: olon_org (:,:)
810  real(RP), intent(out) :: olat_org (:,:)
811  real(RP), intent(out) :: omask_org(:,:)
812  character(len=*), intent(in) :: basename_ocean
813 
814  integer :: rank
815  integer :: xloc, yloc
816  integer :: xs, xe, ys, ye
817 
818  integer :: fid
819  integer :: i
820  !---------------------------------------------------------------------------
821 
822  do i = 1, size( nest_tile_id(:) )
823  ! read data from split files
824  rank = nest_tile_id(i)
825 
826  xloc = mod( i-1, nest_tile_num_x ) + 1
827  yloc = int( real(i-1) / real(nest_tile_num_x) ) + 1
828 
829  xs = parent_imax(handle) * (xloc-1) + 1
830  xe = parent_imax(handle) * xloc
831  ys = parent_jmax(handle) * (yloc-1) + 1
832  ye = parent_jmax(handle) * yloc
833 
834  call file_open( basename_ocean, & ! (in)
835  fid, & ! (out)
836  aggregate=.false., rankid=rank ) ! (in)
837 
838  call file_cartesc_read( fid, "lon", read2d(:,:) )
839  olon_org(xs:xe,ys:ye) = read2d(:,:) * d2r
840 
841  call file_cartesc_read( fid, "lat", read2d(:,:) )
842  olat_org(xs:xe,ys:ye) = read2d(:,:) * d2r
843 
844  call file_cartesc_read( fid, "lsmask", read2d(:,:) )
845  omask_org(xs:xe,ys:ye) = read2d(:,:)
846 
847  end do
848 
849  return

References scale_comm_cartesc_nest::comm_cartesc_nest_tile_id, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_const::const_d2r, scale_file::file_open(), scale_comm_cartesc_nest::parent_imax, and scale_comm_cartesc_nest::parent_jmax.

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ parentoceaninputscale()

subroutine, public mod_realinput_scale::parentoceaninputscale ( real(rp), dimension (:,:), intent(out)  tw_org,
real(rp), dimension (:,:), intent(out)  sst_org,
real(rp), dimension (:,:,:,:), intent(out)  albw_org,
real(rp), dimension (:,:), intent(out)  z0w_org,
real(rp), dimension(:,:), intent(out)  omask_org,
character(len=*), intent(in)  basename_ocean,
integer, intent(in)  it 
)

Definition at line 861 of file mod_realinput_scale.F90.

861  use scale_comm_cartesc_nest, only: &
862  parent_imax, &
863  parent_jmax, &
864  nest_tile_num_x => comm_cartesc_nest_tile_num_x, &
865  nest_tile_id => comm_cartesc_nest_tile_id
866  use scale_file, only: &
867  file_open, &
868  file_get_datainfo
869  use scale_file_cartesc, only: &
870  file_cartesc_read
871  implicit none
872 
873  real(RP), intent(out) :: tw_org (:,:)
874  real(RP), intent(out) :: sst_org (:,:)
875  real(RP), intent(out) :: albw_org (:,:,:,:)
876  real(RP), intent(out) :: z0w_org (:,:)
877  real(RP), intent(out) :: omask_org(:,:)
878  character(len=*), intent(in) :: basename_ocean
879  integer, intent(in) :: it
880 
881  integer :: rank
882 
883  integer :: fid
884  integer :: ndim
885  integer :: i
886  integer :: xloc, yloc
887  integer :: xs, xe
888  integer :: ys, ye
889  !---------------------------------------------------------------------------
890 
891  do i = 1, size( nest_tile_id(:) )
892  ! read data from split files
893  rank = nest_tile_id(i)
894 
895  xloc = mod( i-1, nest_tile_num_x ) + 1
896  yloc = int( real(i-1) / real(nest_tile_num_x) ) + 1
897 
898  xs = parent_imax(handle) * (xloc-1) + 1
899  xe = parent_imax(handle) * xloc
900  ys = parent_jmax(handle) * (yloc-1) + 1
901  ye = parent_jmax(handle) * yloc
902 
903  call file_open( basename_ocean, & ! (in)
904  fid, & ! (out)
905  aggregate=.false., rankid=rank ) ! (in)
906 
907  call file_get_datainfo( fid, "OCEAN_TEMP", dim_rank=ndim )
908  select case (ndim)
909  case ( 2 ) ! for old file
910  call file_cartesc_read( fid, "OCEAN_TEMP", read2d(:,:), step=it )
911  tw_org(xs:xe,ys:ye) = read2d(:,:)
912  case ( 3 )
913  call file_cartesc_read( fid, "OCEAN_TEMP", read3do(:,:,:), step=it )
914  tw_org(xs:xe,ys:ye) = read3do(1,:,:)
915  end select
916 
917  call file_cartesc_read( fid, "OCEAN_SFC_TEMP", read2d(:,:), step=it )
918  sst_org(xs:xe,ys:ye) = read2d(:,:)
919 
920  call file_cartesc_read( fid, "OCEAN_SFC_ALB_IR_dir", read2d(:,:), step=it )
921  albw_org(xs:xe,ys:ye,i_r_direct ,i_r_ir ) = read2d(:,:)
922 
923  call file_cartesc_read( fid, "OCEAN_SFC_ALB_IR_dif", read2d(:,:), step=it )
924  albw_org(xs:xe,ys:ye,i_r_diffuse,i_r_ir ) = read2d(:,:)
925 
926  call file_cartesc_read( fid, "OCEAN_SFC_ALB_NIR_dir", read2d(:,:), step=it )
927  albw_org(xs:xe,ys:ye,i_r_direct ,i_r_nir) = read2d(:,:)
928 
929  call file_cartesc_read( fid, "OCEAN_SFC_ALB_NIR_dif", read2d(:,:), step=it )
930  albw_org(xs:xe,ys:ye,i_r_diffuse,i_r_nir) = read2d(:,:)
931 
932  call file_cartesc_read( fid, "OCEAN_SFC_ALB_VIS_dir", read2d(:,:), step=it )
933  albw_org(xs:xe,ys:ye,i_r_direct ,i_r_vis) = read2d(:,:)
934 
935  call file_cartesc_read( fid, "OCEAN_SFC_ALB_VIS_dif", read2d(:,:), step=it )
936  albw_org(xs:xe,ys:ye,i_r_diffuse,i_r_vis) = read2d(:,:)
937 
938  call file_cartesc_read( fid, "OCEAN_SFC_Z0M", read2d(:,:), step=it )
939  z0w_org(xs:xe,ys:ye) = read2d(:,:)
940 
941  call file_cartesc_read( fid, "lsmask", read2d(:,:) )
942  omask_org(xs:xe,ys:ye) = read2d(:,:)
943 
944  end do
945 
946  return

References scale_comm_cartesc_nest::comm_cartesc_nest_tile_id, scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x, scale_file::file_open(), scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_comm_cartesc_nest::parent_imax, and scale_comm_cartesc_nest::parent_jmax.

Referenced by mod_realinput::realinput_surface().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
scale_atmos_hydrometeor::num_name
character(len=h_short), dimension(n_hyd), parameter, public num_name
Definition: scale_atmos_hydrometeor.F90:92
mod_atmos_phy_mp_vars
module Atmosphere / Physics Cloud Microphysics
Definition: mod_atmos_phy_mp_vars.F90:12
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
mod_atmos_phy_mp_vars::qa_mp
integer, public qa_mp
Definition: mod_atmos_phy_mp_vars.F90:77
scale_atmos_hydrometeor::hyd_name
character(len=h_short), dimension(n_hyd), parameter, public hyd_name
Definition: scale_atmos_hydrometeor.F90:88
scale_comm_cartesc_nest::comm_cartesc_nest_tile_id
integer, dimension(:), allocatable, public comm_cartesc_nest_tile_id
parent tile real id
Definition: scale_comm_cartesC_nest.F90:52
mod_atmos_phy_mp_vars::qs_mp
integer, public qs_mp
Definition: mod_atmos_phy_mp_vars.F90:78
scale_comm_cartesc_nest
module Communication CartesianC nesting
Definition: scale_comm_cartesC_nest.F90:12
scale_atmos_grid_cartesc_metric
module Atmosphere Grid CartesianC metirc
Definition: scale_atmos_grid_cartesC_metric.F90:12
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_file
module file
Definition: scale_file.F90:15
scale_comm_cartesc_nest::parent_jmax
integer, dimension(2), public parent_jmax
parent max number in y-direction
Definition: scale_comm_cartesC_nest.F90:56
mod_atmos_phy_mp_driver
module atmosphere / physics / cloud microphysics
Definition: mod_atmos_phy_mp_driver.F90:12
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_comm_cartesc_nest::parent_lkmax
integer, dimension(2), public parent_lkmax
parent max number in lz-direction
Definition: scale_comm_cartesC_nest.F90:61
scale_file::file_open
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:487
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_y
integer, public comm_cartesc_nest_tile_num_y
parent tile number in y-direction
Definition: scale_comm_cartesC_nest.F90:51
scale_comm_cartesc_nest::parent_kmax
integer, dimension(2), public parent_kmax
parent max number in z-direction
Definition: scale_comm_cartesC_nest.F90:54
scale_comm_cartesc_nest::parent_imax
integer, dimension(2), public parent_imax
parent max number in x-direction
Definition: scale_comm_cartesC_nest.F90:55
scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x
integer, public comm_cartesc_nest_tile_num_x
parent tile number in x-direction
Definition: scale_comm_cartesC_nest.F90:50
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_mapprojection::mapprojection_get_param
subroutine, public mapprojection_get_param(info, param)
Definition: scale_mapprojection.F90:333
scale_atmos_thermodyn
module atmosphere / thermodyn
Definition: scale_atmos_thermodyn.F90:11
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
mod_atmos_phy_mp_vars::qe_mp
integer, public qe_mp
Definition: mod_atmos_phy_mp_vars.F90:79
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
scale_const::const_laps
real(rp), public const_laps
lapse rate of ISA [K/m]
Definition: scale_const.F90:58
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc
subroutine, public atmos_phy_mp_driver_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
Definition: mod_atmos_phy_mp_driver.F90:1361
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
scale_mapprojection::mappingparam
Definition: scale_mapprojection.F90:106
scale_file_cartesc
module file / cartesianC
Definition: scale_file_cartesC.F90:11
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:79
scale_mapprojection::mappinginfo
Definition: scale_mapprojection.F90:93
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
Definition: scale_atmos_grid_cartesC_metric.F90:35