SCALE-RM
Functions/Subroutines | Variables
mod_mkinit Module Reference

module INITIAL More...

Functions/Subroutines

subroutine, public mkinit_setup
 Setup. More...
 
subroutine, public mkinit
 Driver. More...
 
subroutine rect_setup
 Bubble. More...
 
subroutine aerosol_activation_kajino13 (c_kappa, super, temp_k, ia_m0, ia_m2, ia_m3, n_atr, n_siz_max, n_kap_max, n_ctg, n_siz, n_kap, d_ct, aerosol_procs, aerosol_activ)
 
subroutine diag_ds (m0, m2, m3, dg, sg, dm2)
 
real(rp) function faero (f0, r0, x, alpha, rhoa)
 
subroutine flux_setup
 flux setup More...
 
subroutine land_setup
 Land setup. More...
 
subroutine ocean_setup
 Ocean setup. More...
 
subroutine urban_setup
 Urban setup. More...
 
subroutine read_sounding (DENS, VELX, VELY, POTT, QV)
 Read sounding data from file. More...
 
subroutine interporation_fact (fact0, fact1, idx0, idx1, x, x_org, nx, loop)
 

Variables

integer, public mkinit_type = -1
 
integer, parameter, public i_ignore = 0
 
integer, parameter, public i_planestate = 1
 
integer, parameter, public i_tracerbubble = 2
 
integer, parameter, public i_coldbubble = 3
 
integer, parameter, public i_lambwave = 4
 
integer, parameter, public i_gravitywave = 5
 
integer, parameter, public i_khwave = 6
 
integer, parameter, public i_turbulence = 7
 
integer, parameter, public i_mountainwave = 8
 
integer, parameter, public i_warmbubble = 9
 
integer, parameter, public i_supercell = 10
 
integer, parameter, public i_squallline = 11
 
integer, parameter, public i_wk1982 = 12
 
integer, parameter, public i_dycoms2_rf01 = 13
 
integer, parameter, public i_dycoms2_rf02 = 14
 
integer, parameter, public i_rico = 15
 
integer, parameter, public i_interporation = 16
 
integer, parameter, public i_landcouple = 17
 
integer, parameter, public i_oceancouple = 18
 
integer, parameter, public i_urbancouple = 19
 
integer, parameter, public i_triplecouple = 20
 
integer, parameter, public i_bubblecouple = 21
 
integer, parameter, public i_seabreeze = 22
 
integer, parameter, public i_heatisland = 23
 
integer, parameter, public i_dycoms2_rf02_dns = 24
 
integer, parameter, public i_real = 25
 
integer, parameter, public i_grayzone = 26
 
integer, parameter, public i_boxaero = 27
 
integer, parameter, public i_warmbubbleaero = 28
 
integer, parameter, public i_cavityflow = 29
 

Detailed Description

module INITIAL

Description
subroutines for preparing initial data
Author
Team SCALE
History
  • 2011-11-11 (H.Yashiro) [new] Imported from SCALE-RM ver.2
  • 2012-01-31 (Y.Miyamoto) [add] Lamb wave test
  • 2012-01-31 (Y.Miyamoto) [add] KH wave test
  • 2012-01-31 (Y.Miyamoto) [add] turbulence test
  • 2011-02-01 (H.Yashiro) [add] supercell test, follow the supercell test of WRF
  • 2012-02-06 (Y.Miyamoto) [add] advection test
  • 2012-02-16 (Y.Miyamoto) [mod] added hydrostatic balance calculation
  • 2012-03-27 (H.Yashiro) [mod] change subroutines into one module
  • 2012-04-04 (Y.Miyamoto) [new] SQUALLINE test, for GCSS model comparison (Redelsperger et al. 2000)
  • 2012-04-06 (H.Yashiro) [new] uniform state test
  • 2012-04-08 (H.Yashiro) [mod] merge all init programs
  • 2012-06-13 (Y.Sato) [mod] add hbinw option (***HBINW)
  • 2013-02-25 (H.Yashiro) [mod] ISA profile
  • 2014-03-27 (A.Noda) [mod] add DYCOMS2_RF02_DNS
  • 2014-04-28 (R.Yoshida) [add] real case experiment
  • 2014-08-26 (A.Noda) [mod] add GRAYZONE
  • 2015-03-27 (Y.Sato) [mod] add Box aero
  • 2015-04-30 (Y.Sato) [mod] add WARMBUBBLE-AERO
NAMELIST
  • PARAM_MKINIT
    nametypedefault valuecomment
    MKINIT_INITNAME character(len=H_SHORT) 'NONE'
    FLG_INTRP logical .false.

  • PARAM_BUBBLE
    nametypedefault valuecomment
    BBL_EACHNODE logical .false. Arrange bubble at each node? [kg/kg]
    BBL_CZ real(RP) 2.E3_RP center location [m]: z
    BBL_CX real(RP) 2.E3_RP center location [m]: x
    BBL_CY real(RP) 2.E3_RP center location [m]: y
    BBL_RZ real(RP) 0.0_RP bubble radius [m]: z
    BBL_RX real(RP) 0.0_RP bubble radius [m]: x
    BBL_RY real(RP) 0.0_RP bubble radius [m]: y

  • PARAM_RECT
    nametypedefault valuecomment
    RCT_EACHNODE logical .false. Arrange rectangle at each node? [kg/kg]
    RCT_CZ real(RP) 2.E3_RP center location [m]: z
    RCT_CX real(RP) 2.E3_RP center location [m]: x
    RCT_CY real(RP) 2.E3_RP center location [m]: y
    RCT_RZ real(RP) 2.E3_RP rectangle z width [m]: z
    RCT_RX real(RP) 2.E3_RP rectangle x width [m]: x
    RCT_RY real(RP) 2.E3_RP rectangle y width [m]: y

  • PARAM_AERO
    nametypedefault valuecomment
    M0_INIT real(RP) 0.0_RP initial total num. conc. of modes (Atk,Acm,Cor) [#/m3]
    DG_INIT real(RP) 80.e-9_RP initial number equivalen diameters of modes [m]
    SG_INIT real(RP) 1.6_RP initial standard deviation [-]
    D_MIN_INP real(RP) d_min_def
    D_MAX_INP real(RP) d_max_def
    K_MIN_INP real(RP) k_min_def
    K_MAX_INP real(RP) k_max_def
    N_KAP_INP integer n_kap_def

  • PARAM_SBMAERO
    nametypedefault valuecomment
    F0_AERO real(RP) 1.E+7_RP
    R0_AERO real(RP) 1.E-7_RP
    R_MAX real(RP) 1.E-06_RP
    R_MIN real(RP) 1.E-08_RP
    A_ALPHA real(RP) 3.0_RP
    RHO_AERO real(RP) 2.25E+03_RP

  • PARAM_MKINIT_FLUX
    nametypedefault valuecomment
    FLX_RAIN real(RP) 0.0_RP surface rain flux [kg/m2/s]
    FLX_SNOW real(RP) 0.0_RP surface snow flux [kg/m2/s]
    FLX_LW_DN real(RP) 0.0_RP surface downwad long-wave radiation flux [J/m2/s]
    FLX_SW_DN real(RP) 0.0_RP surface downwad short-wave radiation flux [J/m2/s]

  • PARAM_MKINIT_LAND
    nametypedefault valuecomment
    LND_TEMP real(RP) soil temperature [K]
    LND_WATER real(RP) 0.15_RP soil moisture [m3/m3]
    SFC_TEMP real(RP) ocean skin temperature [K]
    SFC_ALBEDO_LW real(RP) 0.04_RP ocean surface albedo for LW [0-1]
    SFC_ALBEDO_SW real(RP) 0.05_RP ocean surface albedo for SW [0-1]

  • PARAM_MKINIT_OCEAN
    nametypedefault valuecomment
    OCN_TEMP real(RP) ocean temperature [K]
    SFC_TEMP real(RP) ocean skin temperature [K]
    SFC_ALBEDO_LW real(RP) 0.04_RP ocean surface albedo for LW [0-1]
    SFC_ALBEDO_SW real(RP) 0.05_RP ocean surface albedo for SW [0-1]
    SFC_Z0M real(RP) 1.0e-4_RP ocean surface roughness length (momentum) [m]
    SFC_Z0H real(RP) 1.0e-4_RP ocean surface roughness length (heat) [m]
    SFC_Z0E real(RP) 1.0e-4_RP ocean surface roughness length (vapor) [m]

  • PARAM_MKINIT_URBAN
    nametypedefault valuecomment
    URB_ROOF_TEMP real(RP) Surface temperature of roof [K]
    URB_BLDG_TEMP real(RP) Surface temperature of building [K
    URB_GRND_TEMP real(RP) Surface temperature of ground [K]
    URB_CNPY_TEMP real(RP) Diagnostic canopy air temperature
    URB_CNPY_HMDT real(RP) 0.0_RP Diagnostic canopy humidity [-]
    URB_CNPY_WIND real(RP) 0.0_RP Diagnostic canopy wind [m/s]
    URB_ROOF_LAYER_TEMP real(RP) temperature in layer of roof [K]
    URB_BLDG_LAYER_TEMP real(RP) temperature in layer of building [
    URB_GRND_LAYER_TEMP real(RP) temperature in layer of ground [K]
    URB_ROOF_RAIN real(RP) 0.0_RP temperature in layer of roof [K]
    URB_BLDG_RAIN real(RP) 0.0_RP temperature in layer of building [
    URB_GRND_RAIN real(RP) 0.0_RP temperature in layer of ground [K]
    URB_RUNOFF real(RP) 0.0_RP temperature in layer of ground [K]
    URB_SFC_TEMP real(RP) Grid average of surface temperature [K]
    URB_ALB_LW real(RP) 0.0_RP Grid average of surface albedo for LW [0-1]
    URB_ALB_SW real(RP) 0.0_RP Grid average of surface albedo for SW [0-1]

  • PARAM_MKINIT_SOUNDING
    nametypedefault valuecomment
    ENV_IN_SOUNDING_FILE character(len=H_LONG) ''

  • PARAM_MKINIT_PLANESTATE
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    SFC_RH real(RP) 80.0_RP surface relative humidity [%]
    ENV_THETA real(RP) potential temperature of environment [K]
    ENV_TLAPS real(RP) 4.E-3_RP Lapse rate of THETA [K/m]
    ENV_U real(RP) 0.0_RP velocity u of environment [m/s]
    ENV_V real(RP) 0.0_RP velocity v of environment [m/s]
    ENV_RH real(RP) 80.0_RP Relative Humidity of environment [%]
    RANDOM_THETA real(RP) 0.01_RP
    RANDOM_U real(RP) 0.0_RP amplitude of random disturbance u
    RANDOM_V real(RP) 0.0_RP amplitude of random disturbance v
    RANDOM_RH real(RP) 0.0_RP amplitude of random disturbance RH

  • PARAM_MKINIT_TRACERBUBBLE
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    ENV_THETA real(RP) potential temperature of environment [K]
    ENV_U real(RP) 0.0_RP velocity u of environment [m/s]
    ENV_V real(RP) 0.0_RP velocity v of environment [m/s]
    SHAPE_NC character(len=H_SHORT) 'BUBBLE' BUBBLE or RECT
    BBL_NC real(RP) 0.0_RP extremum of NC in bubble [kg/kg]

  • PARAM_MKINIT_COLDBUBBLE
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    ENV_THETA real(RP) potential temperature of environment [K]
    BBL_TEMP real(RP) -15.0_RP extremum of temperature in bubble [K]

  • PARAM_MKINIT_LAMBWAVE
    nametypedefault valuecomment
    SFC_PRES real(RP) surface pressure [Pa]
    ENV_U real(RP) 0.0_RP velocity u of environment [m/s]
    ENV_V real(RP) 0.0_RP velocity v of environment [m/s]
    ENV_TEMP real(RP) 300.0_RP temperature of environment [K]
    BBL_PRES real(RP) 100._RP extremum of pressure in bubble [Pa]

  • PARAM_MKINIT_GRAVITYWAVE
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    ENV_U real(RP) 0.0_RP velocity u of environment [m/s]
    ENV_V real(RP) 0.0_RP velocity v of environment [m/s]
    ENV_BVF real(RP) 0.01_RP Brunt Vaisala frequencies of environment [1/s]
    BBL_THETA real(RP) 1.0_RP extremum of temperature in bubble [K]

  • PARAM_MKINIT_KHWAVE
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    ENV_L1_ZTOP real(RP) 1.E3_RP top height of the layer1 (constant THETA) [m]
    ENV_L3_ZBOTTOM real(RP) 2100.0_RP bottom height of the layer3 (high THETA) [m]
    ENV_L1_THETA real(RP) 300.0_RP THETA in the layer1 (low THETA) [K]
    ENV_L3_THETA real(RP) 301.0_RP THETA in the layer3 (high THETA) [K]
    ENV_L1_U real(RP) 0.0_RP velocity u in the layer1 (low THETA) [K]
    ENV_L3_U real(RP) 20.0_RP velocity u in the layer3 (high THETA) [K]
    RANDOM_U real(RP) 0.0_RP amplitude of random disturbance u

  • PARAM_MKINIT_TURBULENCE
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    SFC_RH real(RP) 80.0_RP surface relative humidity [%]
    ENV_THETA real(RP) potential temperature of environment [K]
    ENV_TLAPS real(RP) 4.E-3_RP Lapse rate of THETA [K/m]
    ENV_U real(RP) 0.0_RP velocity u of environment [m/s]
    ENV_V real(RP) 0.0_RP velocity v of environment [m/s]
    ENV_RH real(RP) 80.0_RP Relative Humidity of environment [%]
    RANDOM_THETA real(RP) 0.01_RP
    RANDOM_U real(RP) 0.0_RP amplitude of random disturbance u
    RANDOM_V real(RP) 0.0_RP amplitude of random disturbance v
    RANDOM_RH real(RP) 0.0_RP amplitude of random disturbance RH

  • PARAM_MKINIT_CAVITYFLOW
    nametypedefault valuecomment
    REYNOLDS_NUM real(RP) 4.D02
    MACH_NUM real(RP) 3.D-2
    TEMP0 real(RP) 300.D0

  • PARAM_MKINIT_MOUNTAINWAVE
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    ENV_U real(RP) 0.0_RP velocity u of environment [m/s]
    ENV_V real(RP) 0.0_RP velocity v of environment [m/s]
    SCORER real(RP) 2.E-3_RP Scorer parameter (~=N/U) [1/m]
    BBL_NC real(RP) 0.0_RP extremum of NC in bubble [kg/kg]

  • PARAM_MKINIT_WARMBUBBLE
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    ENV_U real(RP) 0.0_RP velocity u of environment [m/s]
    ENV_V real(RP) 0.0_RP velocity v of environment [m/s]
    ENV_RH real(RP) 80.0_RP Relative Humidity of environment [%]
    ENV_L1_ZTOP real(RP) 1.E3_RP top height of the layer1 (constant THETA) [m]
    ENV_L2_ZTOP real(RP) 14.E3_RP top height of the layer2 (small THETA gradient) [m]
    ENV_L2_TLAPS real(RP) 4.E-3_RP Lapse rate of THETA in the layer2 (small THETA gradient) [K/m]
    ENV_L3_TLAPS real(RP) 3.E-2_RP Lapse rate of THETA in the layer3 (large THETA gradient) [K/m]
    BBL_THETA real(RP) 1.0_RP extremum of temperature in bubble [K]

  • PARAM_MKINIT_SUPERCELL
    nametypedefault valuecomment
    BBL_THETA real(RP) 1.0_RP extremum of temperature in bubble [K]

  • PARAM_MKINIT_SQUALLLINE
    nametypedefault valuecomment
    RANDOM_THETA real(RP) 0.01_RP
    OFFSET_VELX real(RP) 12.0_RP
    OFFSET_VELY real(RP) -2.0_RP

  • PARAM_MKINIT_WK1982
    nametypedefault valuecomment
    SFC_THETA real(RP) surface potential temperature [K]
    SFC_PRES real(RP) surface pressure [Pa]
    TR_Z real(RP) 12000.0_RP height of tropopause [m]
    TR_THETA real(RP) 343.0_RP pot. temperature at tropopause [K]
    TR_TEMP real(RP) 213.0_RP temperature at tropopause [K]
    SHEAR_Z real(RP) 3000.0_RP center height of shear layer [m]
    SHEAR_U real(RP) 15.0_RP velocity u over the shear layer [m/s]
    BBL_THETA real(RP) 1.0_RP extremum of temperature in bubble [K]

  • PARAM_MKINIT_RF01
    nametypedefault valuecomment
    PERTURB_AMP real(RP) 0.0_RP
    RANDOM_LIMIT integer 0
    RANDOM_FLAG integer 0 0 -> no perturbation
    USE_LWSET logical .false. use liq. water. static energy temp.?

  • PARAM_MKINIT_RF02
    nametypedefault valuecomment
    PERTURB_AMP real(RP) 0.0_RP
    RANDOM_LIMIT integer 0
    RANDOM_FLAG integer 0 0 -> no perturbation

  • PARAM_MKINIT_RF02_DNS
    nametypedefault valuecomment
    ZB real(RP) 750.0_RP domain bottom
    CONST_U real(RP) 0.0_RP
    CONST_V real(RP) 0.0_RP
    PRES_ZB real(RP) 93060.0_RP
    PERTURB_AMP real(RP) 0.0_RP
    RANDOM_LIMIT integer 0
    RANDOM_FLAG integer 0 0 -> no perturbation

  • PARAM_MKINIT_RICO
    nametypedefault valuecomment
    PERTURB_AMP_PT real(RP) 0.1_RP
    PERTURB_AMP_QV real(RP) 2.5E-5_RP

  • PARAM_MKINIT_INTERPORATION
    nametypedefault valuecomment
    BASENAME_ORG character(len=H_LONG) ''

  • PARAM_MKINIT_GRAYZONE
    nametypedefault valuecomment
    PERTURB_AMP real(RP) 0.0_RP
    RANDOM_LIMIT integer 0
    RANDOM_FLAG integer 0 0 -> no perturbation

  • PARAM_MKINIT_BOXAERO
    nametypedefault valuecomment
    INIT_DENS real(RP) 1.12_RP [kg/m3]
    INIT_TEMP real(RP) 298.18_RP [K]
    INIT_PRES real(RP) 1.E+5_RP [Pa]
    INIT_SSLIQ real(RP) 0.01_RP [%]

History Output
No history output

Function/Subroutine Documentation

◆ mkinit_setup()

subroutine, public mod_mkinit::mkinit_setup ( )

Setup.

Definition at line 221 of file mod_mkinit.f90.

References i_boxaero, i_bubblecouple, i_cavityflow, i_coldbubble, i_dycoms2_rf01, i_dycoms2_rf02, i_dycoms2_rf02_dns, i_gravitywave, i_grayzone, i_heatisland, i_ignore, i_interporation, i_khwave, i_lambwave, i_landcouple, i_mountainwave, i_oceancouple, i_planestate, i_real, i_rico, i_seabreeze, i_squallline, i_supercell, i_tracerbubble, i_triplecouple, i_turbulence, i_urbancouple, i_warmbubble, i_warmbubbleaero, i_wk1982, scale_grid_index::ia, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::ja, scale_grid_index::ka, mkinit_type, and scale_process::prc_mpistop().

Referenced by mod_rm_prep::scalerm_prep().

221  implicit none
222 
223  character(len=H_SHORT) :: mkinit_initname = 'NONE'
224 
225  namelist / param_mkinit / &
226  mkinit_initname, &
227  flg_intrp
228 
229  integer :: ierr
230  !---------------------------------------------------------------------------
231 
232  if( io_l ) write(io_fid_log,*)
233  if( io_l ) write(io_fid_log,*) '++++++ Module[make init] / Categ[preprocess] / Origin[SCALE-RM]'
234 
235  !--- read namelist
236  rewind(io_fid_conf)
237  read(io_fid_conf,nml=param_mkinit,iostat=ierr)
238  if( ierr < 0 ) then !--- missing
239  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
240  elseif( ierr > 0 ) then !--- fatal error
241  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT. Check!'
242  call prc_mpistop
243  endif
244  if( io_lnml ) write(io_fid_log,nml=param_mkinit)
245 
246  allocate( pres(ka,ia,ja) )
247  allocate( temp(ka,ia,ja) )
248  allocate( pott(ka,ia,ja) )
249  allocate( qsat(ka,ia,ja) )
250  allocate( qv(ka,ia,ja) )
251  allocate( qc(ka,ia,ja) )
252  allocate( velx(ka,ia,ja) )
253  allocate( vely(ka,ia,ja) )
254 
255  allocate( pres_sfc(1,ia,ja) )
256  allocate( temp_sfc(1,ia,ja) )
257  allocate( pott_sfc(1,ia,ja) )
258  allocate( qsat_sfc(1,ia,ja) )
259  allocate( qv_sfc(1,ia,ja) )
260  allocate( qc_sfc(1,ia,ja) )
261 
262  allocate( rndm(ka,ia,ja) )
263  allocate( bubble(ka,ia,ja) )
264  allocate( rect(ka,ia,ja) )
265 
266  select case(trim(mkinit_initname))
267  case('NONE')
268  mkinit_type = i_ignore
269  case('PLANESTATE')
270  mkinit_type = i_planestate
271  case('TRACERBUBBLE')
272  mkinit_type = i_tracerbubble
273  case('COLDBUBBLE')
274  mkinit_type = i_coldbubble
275  call bubble_setup
276  case('LAMBWAVE')
277  mkinit_type = i_lambwave
278  call bubble_setup
279  case('GRAVITYWAVE')
280  mkinit_type = i_gravitywave
281  call bubble_setup
282  case('KHWAVE')
283  mkinit_type = i_khwave
284  case('TURBULENCE')
285  mkinit_type = i_turbulence
286  case('MOUNTAINWAVE')
287  mkinit_type = i_mountainwave
288  call bubble_setup
289  case('WARMBUBBLE')
290  mkinit_type = i_warmbubble
291  call bubble_setup
292  case('SUPERCELL')
293  mkinit_type = i_supercell
294  call bubble_setup
295  case('SQUALLLINE')
296  mkinit_type = i_squallline
297  case('WK1982')
298  mkinit_type = i_wk1982
299  call bubble_setup
300  case('DYCOMS2_RF01')
301  mkinit_type = i_dycoms2_rf01
302  case('DYCOMS2_RF02')
303  mkinit_type = i_dycoms2_rf02
304  case('RICO')
305  mkinit_type = i_rico
306  case('INTERPORATION')
307  mkinit_type = i_interporation
308  case('LANDCOUPLE')
309  mkinit_type = i_landcouple
310  case('OCEANCOUPLE')
311  mkinit_type = i_oceancouple
312  case('URBANCOUPLE')
313  mkinit_type = i_urbancouple
314  case('TRIPLECOUPLE')
315  mkinit_type = i_triplecouple
316  case('BUBBLECOUPLE')
317  mkinit_type = i_bubblecouple
318  call bubble_setup
319  case('SEABREEZE')
320  mkinit_type = i_seabreeze
321  case('HEATISLAND')
322  mkinit_type = i_heatisland
323  case('DYCOMS2_RF02_DNS')
324  mkinit_type = i_dycoms2_rf02_dns
325  case('REAL')
326  mkinit_type = i_real
327  case('GRAYZONE')
328  mkinit_type = i_grayzone
329  case('BOXAERO')
330  mkinit_type = i_boxaero
331  case('WARMBUBBLEAERO')
332  mkinit_type = i_warmbubbleaero
333  call bubble_setup
334  case('CAVITYFLOW')
335  mkinit_type = i_cavityflow
336  case default
337  write(*,*) ' xxx Unsupported TYPE:', trim(mkinit_initname)
338  call prc_mpistop
339  endselect
340 
341  return
subroutine, public prc_mpistop
Abort MPI.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ mkinit()

subroutine, public mod_mkinit::mkinit ( )

Driver.

Definition at line 347 of file mod_mkinit.f90.

References mod_admin_restart::admin_restart(), mod_atmos_vars::atmos_restart_output, mod_atmos_driver::atmos_surface_get(), mod_atmos_vars::atmos_vars_restart_write(), scale_const::const_undef8, scale_grid::grid_cx, scale_grid::grid_cy, scale_grid::grid_cz, i_boxaero, i_bubblecouple, i_cavityflow, i_coldbubble, i_dycoms2_rf01, i_dycoms2_rf02, i_dycoms2_rf02_dns, i_gravitywave, i_grayzone, i_heatisland, i_ignore, i_interporation, i_khwave, i_lambwave, i_landcouple, i_mountainwave, i_oceancouple, i_planestate, i_real, i_rico, i_seabreeze, i_squallline, i_supercell, i_tracerbubble, i_triplecouple, i_turbulence, i_urbancouple, i_warmbubble, i_warmbubbleaero, i_wk1982, scale_grid_index::ia, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, mod_land_vars::land_restart_output, mod_land_driver::land_surface_set(), mod_land_vars::land_vars_restart_write(), scale_landuse::landuse_write(), mkinit_type, mod_ocean_vars::ocean_restart_output, mod_ocean_driver::ocean_surface_set(), mod_ocean_vars::ocean_vars_restart_write(), scale_process::prc_mpistop(), scale_tracer::qa, mod_admin_time::time_doatmos_restart, mod_admin_time::time_doland_restart, mod_admin_time::time_doocean_restart, mod_admin_time::time_dourban_restart, scale_tracer::tracer_type, mod_urban_vars::urban_restart_output, mod_urban_driver::urban_surface_set(), and mod_urban_vars::urban_vars_restart_write().

Referenced by mod_rm_prep::scalerm_prep().

347  use scale_const, only: &
349  use scale_landuse, only: &
351  use mod_atmos_driver, only: &
353  use mod_atmos_vars, only: &
354  atmos_sw_restart => atmos_restart_output, &
356  use mod_ocean_driver, only: &
358  use mod_ocean_vars, only: &
359  ocean_sw_restart => ocean_restart_output, &
361  use mod_land_driver, only: &
363  use mod_land_vars, only: &
364  land_sw_restart => land_restart_output, &
366  use mod_urban_driver, only: &
368  use mod_urban_vars, only: &
369  urban_sw_restart => urban_restart_output, &
371  use mod_admin_restart, only: &
373  use mod_admin_time, only: &
378  implicit none
379 
380  integer :: iq
381  !---------------------------------------------------------------------------
382 
383  if ( mkinit_type == i_ignore ) then
384  if( io_l ) write(io_fid_log,*)
385  if( io_l ) write(io_fid_log,*) '++++++ SKIP MAKING INITIAL DATA ++++++'
386  else
387  if( io_l ) write(io_fid_log,*)
388  if( io_l ) write(io_fid_log,*) '++++++ START MAKING INITIAL DATA ++++++'
389 
390  !--- Initialize variables
391 #ifndef DRY
392  do iq = 2, qa
393  qtrc(:,:,:,iq) = 0.0_rp
394  enddo
395 #endif
396 
397  pres(:,:,:) = const_undef8
398  temp(:,:,:) = const_undef8
399  pott(:,:,:) = const_undef8
400  qsat(:,:,:) = const_undef8
401  qv(:,:,:) = const_undef8
402  qc(:,:,:) = const_undef8
403  velx(:,:,:) = const_undef8
404  vely(:,:,:) = const_undef8
405 
406  rndm(:,:,:) = const_undef8
407 
408  pres_sfc(:,:,:) = const_undef8
409  temp_sfc(:,:,:) = const_undef8
410  pott_sfc(:,:,:) = const_undef8
411  qsat_sfc(:,:,:) = const_undef8
412  qv_sfc(:,:,:) = const_undef8
413  qc_sfc(:,:,:) = const_undef8
414 
415  if ( tracer_type == 'SUZUKI10' ) then
416  if( io_l ) write(io_fid_log,*) '*** Aerosols for SBM are included ***'
417  call sbmaero_setup
418  endif
419 
420  select case(mkinit_type)
421  case(i_planestate)
422  call mkinit_planestate
423  case(i_tracerbubble)
424  call mkinit_tracerbubble
425  case(i_coldbubble)
426  call mkinit_coldbubble
427  case(i_lambwave)
428  call mkinit_lambwave
429  case(i_gravitywave)
430  call mkinit_gravitywave
431  case(i_khwave)
432  call mkinit_khwave
433  case(i_turbulence)
434  call mkinit_turbulence
435  case(i_mountainwave)
436  call mkinit_mountainwave
437  case(i_warmbubble)
438  call mkinit_warmbubble
439  case(i_supercell)
440  call mkinit_supercell
441  case(i_squallline)
442  call mkinit_squallline
443  case(i_wk1982)
444  call mkinit_wk1982
445  case(i_dycoms2_rf01)
446  call mkinit_dycoms2_rf01
447  case(i_dycoms2_rf02)
448  call mkinit_dycoms2_rf02
449  case(i_rico)
450  call mkinit_rico
451  case(i_interporation)
452  call mkinit_interporation
453  case(i_oceancouple)
454  call mkinit_planestate
455  call mkinit_oceancouple
456  case(i_landcouple)
457  call mkinit_planestate
458  call mkinit_landcouple
459  case(i_urbancouple)
460  call mkinit_planestate
461  call mkinit_landcouple ! tentative
462  call mkinit_urbancouple
463  case(i_triplecouple)
464  call mkinit_planestate
465  call mkinit_oceancouple
466  call mkinit_landcouple
467  call mkinit_urbancouple
468  case(i_bubblecouple)
469  call mkinit_planestate
470  call mkinit_warmbubble
471  call mkinit_oceancouple
472  call mkinit_landcouple
473  call mkinit_urbancouple
474  case(i_seabreeze)
475  call mkinit_planestate
476  call mkinit_seabreeze
477  case(i_heatisland)
478  call mkinit_planestate
479  call mkinit_heatisland
480  case(i_dycoms2_rf02_dns)
481  call mkinit_dycoms2_rf02_dns
482  case(i_real)
483  call mkinit_real
484  case(i_grayzone)
485  call mkinit_grayzone
486  case(i_boxaero)
487  call mkinit_boxaero
488  case(i_warmbubbleaero)
489  call mkinit_warmbubbleaero
490  case(i_cavityflow)
491  call mkinit_cavityflow
492  case default
493  write(*,*) ' xxx Unsupported TYPE:', mkinit_type
494  call prc_mpistop
495  endselect
496 
497  if( io_l ) write(io_fid_log,*) '++++++ END MAKING INITIAL DATA ++++++'
498 
499  ! setup surface condition
500  call ocean_surface_set( countup = .false. )
501  call land_surface_set ( countup = .false. )
502  call urban_surface_set( countup = .false. )
503  call atmos_surface_get
504 
505  ! output boundary file
506  call landuse_write
507 
508  ! output restart file
509  ! if( ATMOS_sw_restart ) call ATMOS_vars_restart_write
510  ! if( OCEAN_sw_restart ) call OCEAN_vars_restart_write
511  ! if( LAND_sw_restart ) call LAND_vars_restart_write
512  ! if( URBAN_sw_restart ) call URBAN_vars_restart_write
513  time_doocean_restart = .true.
514  time_doland_restart = .true.
515  time_dourban_restart = .true.
516  time_doatmos_restart = .true.
517  call admin_restart
518 
519  endif
520 
521  return
subroutine, public prc_mpistop
Abort MPI.
logical, public land_restart_output
output restart file?
module URBAN driver
subroutine, public ocean_vars_restart_write
Write ocean restart.
module ATMOSPHERIC Variables
subroutine, public atmos_vars_restart_write
Write restart of atmospheric variables.
logical, public time_doatmos_restart
execute atmosphere restart output in this step?
logical, public time_doocean_restart
execute ocean restart output in this step?
module URBAN Variables
logical, public atmos_restart_output
output restart file?
subroutine, public urban_surface_set(countup)
Set surface boundary to other model.
subroutine, public urban_vars_restart_write
Write urban restart.
subroutine, public atmos_surface_get
Get surface boundary condition.
module LANDUSE
module OCEAN driver
subroutine, public ocean_surface_set(countup)
Put surface boundary to other model.
module ATMOSPHERE driver
subroutine, public landuse_write
Write landuse data.
module LAND Variables
module administrator for restart
module CONSTANT
Definition: scale_const.F90:14
logical, public time_dourban_restart
execute urban restart output in this step?
logical, public urban_restart_output
output restart file?
real(dp), parameter, public const_undef8
undefined value (REAL8)
Definition: scale_const.F90:42
logical, public time_doland_restart
execute land restart output in this step?
module ADMIN TIME
subroutine, public admin_restart
logical, public ocean_restart_output
output restart file?
subroutine, public land_vars_restart_write
Write land restart.
module LAND driver
subroutine, public land_surface_set(countup)
Put surface boundary to other model.
module OCEAN Variables
Here is the call graph for this function:
Here is the caller graph for this function:

◆ rect_setup()

subroutine mod_mkinit::rect_setup ( )

Bubble.

Definition at line 612 of file mod_mkinit.f90.

References scale_tracer::ae_ctg, aerosol_activation_kajino13(), scale_atmos_thermodyn::aq_cv, scale_const::const_cvdry, scale_const::const_pi, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_rvap, scale_const::const_undef8, mod_atmos_vars::dens, scale_tracer::gas_ctg, scale_grid::grid_cx, scale_grid::grid_cy, scale_grid::grid_cz, scale_tracer::i_qv, scale_grid_index::ia, scale_tracer::ic_mix, scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, dc_log::log(), scale_tracer::n_atr, scale_tracer::nkap, scale_tracer::nsiz, scale_process::prc_mpistop(), scale_tracer::qaee, scale_tracer::qaes, scale_tracer::qqe, scale_tracer::qqs, mod_atmos_vars::qtrc, and mod_atmos_vars::rhot.

Referenced by read_sounding().

612  use scale_const, only: &
614  implicit none
615 
616  ! Bubble
617  logical :: rct_eachnode = .false. ! Arrange rectangle at each node? [kg/kg]
618  real(RP) :: rct_cz = 2.e3_rp ! center location [m]: z
619  real(RP) :: rct_cx = 2.e3_rp ! center location [m]: x
620  real(RP) :: rct_cy = 2.e3_rp ! center location [m]: y
621  real(RP) :: rct_rz = 2.e3_rp ! rectangle z width [m]: z
622  real(RP) :: rct_rx = 2.e3_rp ! rectangle x width [m]: x
623  real(RP) :: rct_ry = 2.e3_rp ! rectangle y width [m]: y
624 
625  namelist / param_rect / &
626  rct_eachnode, &
627  rct_cz, &
628  rct_cx, &
629  rct_cy, &
630  rct_rz, &
631  rct_rx, &
632  rct_ry
633 
634  real(RP) :: cz_offset
635  real(RP) :: cx_offset
636  real(RP) :: cy_offset
637  real(RP) :: dist
638 
639  integer :: ierr
640  integer :: k, i, j
641  !---------------------------------------------------------------------------
642 
643  if( io_l ) write(io_fid_log,*)
644  if( io_l ) write(io_fid_log,*) '++++++ Module[mkinit rectangle] / Categ[preprocess] / Origin[SCALE-RM]'
645 
646  !--- read namelist
647  rewind(io_fid_conf)
648  read(io_fid_conf,nml=param_rect,iostat=ierr)
649 
650  if( ierr < 0 ) then !--- missing
651  if( io_l ) write(io_fid_log,*) 'xxx Not found namelist. Check!'
652  call prc_mpistop
653  elseif( ierr > 0 ) then !--- fatal error
654  write(*,*) 'xxx Not appropriate names in namelist PARAM_RECT. Check!'
655  call prc_mpistop
656  endif
657  if( io_lnml ) write(io_fid_log,nml=param_rect)
658 
659  rect(:,:,:) = const_undef8
660 
661  if ( rct_eachnode ) then
662  cz_offset = grid_cz(ks)
663  cx_offset = grid_cx(is)
664  cy_offset = grid_cy(js)
665  else
666  cz_offset = 0.0_rp
667  cx_offset = 0.0_rp
668  cy_offset = 0.0_rp
669  endif
670 
671  do j = 1, ja
672  do i = 1, ia
673  do k = ks, ke
674 
675  ! make tracer rectangle
676  dist = 2.0_rp * max( &
677  abs(grid_cz(k) - cz_offset - rct_cz)/rct_rz, &
678  abs(grid_cx(i) - cx_offset - rct_cx)/rct_rx, &
679  abs(grid_cy(j) - cy_offset - rct_cy)/rct_ry &
680  & )
681  if ( dist <= 1.0_rp ) then
682  rect(k,i,j) = 1.0_rp
683  else
684  rect(k,i,j) = 0.0_rp
685  end if
686  enddo
687  enddo
688  enddo
689 
690  return
subroutine, public prc_mpistop
Abort MPI.
module CONSTANT
Definition: scale_const.F90:14
real(dp), parameter, public const_undef8
undefined value (REAL8)
Definition: scale_const.F90:42
Here is the call graph for this function:
Here is the caller graph for this function:

◆ aerosol_activation_kajino13()

subroutine mod_mkinit::aerosol_activation_kajino13 ( real(rp), intent(in)  c_kappa,
real(rp), intent(in)  super,
real(rp), intent(in)  temp_k,
integer, intent(in)  ia_m0,
integer, intent(in)  ia_m2,
integer, intent(in)  ia_m3,
integer, intent(in)  n_atr,
integer, intent(in)  n_siz_max,
integer, intent(in)  n_kap_max,
integer, intent(in)  n_ctg,
integer, dimension(n_ctg), intent(in)  n_siz,
integer, dimension(n_ctg), intent(in)  n_kap,
real(rp), dimension(n_siz_max,n_ctg), intent(in)  d_ct,
real(rp), dimension(n_atr,n_siz_max,n_kap_max,n_ctg)  aerosol_procs,
real(rp), dimension(n_atr,n_siz_max,n_kap_max,n_ctg)  aerosol_activ 
)

Definition at line 982 of file mod_mkinit.f90.

References scale_const::const_mvap, diag_ds(), and dc_log::log().

Referenced by rect_setup().

982  use scale_const, only: &
983  const_mvap, & ! mean molecular weight for water vapor [g/mol]
984  const_dwatr, & ! water density [kg/m3]
985  const_r, & ! universal gas constant [J/mol/K]
986  const_tem00
987  implicit none
988  !i/o variables
989  real(RP),intent(in) :: super ! supersaturation [-]
990  real(RP),intent(in) :: c_kappa ! hygroscopicity of condensable mass [-]
991  real(RP),intent(in) :: temp_k ! temperature
992  integer, intent(in) :: ia_m0, ia_m2, ia_m3
993  integer, intent(in) :: n_atr
994  integer, intent(in) :: n_siz_max
995  integer, intent(in) :: n_kap_max
996  integer, intent(in) :: n_ctg
997  real(RP),intent(in) :: d_ct(n_siz_max,n_ctg)
998  real(RP) :: aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg)
999  real(RP) :: aerosol_activ(n_atr,n_siz_max,n_kap_max,n_ctg)
1000  integer, intent(in) :: n_siz(n_ctg), n_kap(n_ctg)
1001  !local variables
1002  real(RP),parameter :: two3 = 2._rp/3._rp
1003  real(RP),parameter :: rt2 = sqrt(2._rp)
1004  real(RP),parameter :: twort2 = rt2 ! 2/sqrt(2) = sqrt(2)
1005  real(RP),parameter :: thrrt2 = 3._rp/rt2 ! 3/sqrt(2)
1006  real(RP) :: smax_inv ! inverse
1007  real(RP) :: am,scrit_am,aa,tc,st,bb,ac
1008  real(RP) :: m0t,m2t,m3t,dgt,sgt,dm2
1009  real(RP) :: d_crit ! critical diameter
1010  real(RP) :: tmp1, tmp2, tmp3 !
1011  real(RP) :: ccn_frc,cca_frc,ccv_frc ! activated number,area,volume
1012  integer :: is0, ik, ic
1013 
1014  aerosol_activ(:,:,:,:) = 0._rp
1015 
1016  if (super.le.0._rp) return
1017 
1018  smax_inv = 1._rp / super
1019 
1020  !--- surface tension of water
1021  tc = temp_k - const_tem00
1022  if (tc >= 0._rp ) then
1023  st = 75.94_rp-0.1365_rp*tc-0.3827e-3_rp*tc**2._rp !Gittens, JCIS, 69, Table 4.
1024  else !t[deg C].lt.0.
1025  st = 75.93_rp +0.115_rp*tc & !Eq.5-12, pp.130, PK97
1026  + 6.818e-2_rp*tc**2._rp+6.511e-3_rp*tc**3._rp &
1027  + 2.933e-4_rp*tc**4._rp+6.283e-6_rp*tc**5._rp &
1028  + 5.285e-8_rp*tc**6._rp
1029  endif
1030  st = st * 1.e-3_rp ![J/m2]
1031 
1032  !-- Kelvin effect
1033  ! [J m-2] [kg mol-1] [m3 kg-1] [mol K J-1] [K-1]
1034  aa = 2._rp * st * const_mvap * 1.e-3_rp / (const_dwatr * const_r * temp_k ) ![m] Eq.5 in AR98
1035 
1036  do ic = 1, n_ctg
1037  do ik = 1, n_kap(ic)
1038  do is0 = 1, n_siz(ic)
1039  m0t = aerosol_procs(ia_m0,is0,ik,ic)
1040  m2t = aerosol_procs(ia_m2,is0,ik,ic)
1041  m3t = aerosol_procs(ia_m3,is0,ik,ic)
1042  call diag_ds(m0t,m2t,m3t,dgt,sgt,dm2)
1043  if (dgt <= 0._rp) dgt = d_ct(is0,ic)
1044  if (sgt <= 0._rp) sgt = 1.3_rp
1045  am = dgt * 0.5_rp !geometric dry mean radius [m]
1046  bb = c_kappa
1047  if (bb > 0._rp .and. am > 0._rp ) then
1048  scrit_am = 2._rp/sqrt(bb)*(aa/(3._rp*am))**1.5_rp !AR00 Eq.9
1049  else
1050  scrit_am = 0._rp
1051  endif
1052  ac = am * (scrit_am * smax_inv)**two3 !AR00 Eq.12
1053  d_crit = ac * 2._rp
1054  tmp1 = log(d_crit) - log(dgt)
1055  tmp2 = 1._rp/(rt2*log(sgt))
1056  tmp3 = log(sgt)
1057  ccn_frc= 0.5_rp*(1._rp-erf(tmp1*tmp2))
1058  cca_frc= 0.5_rp*(1._rp-erf(tmp1*tmp2-twort2*tmp3))
1059  ccv_frc= 0.5_rp*(1._rp-erf(tmp1*tmp2-thrrt2*tmp3))
1060  aerosol_activ(ia_m0,is0,ik,ic) = ccn_frc * aerosol_procs(ia_m0,is0,ik,ic)
1061  aerosol_activ(ia_m2,is0,ik,ic) = cca_frc * aerosol_procs(ia_m2,is0,ik,ic)
1062  aerosol_activ(ia_m3,is0,ik,ic) = ccv_frc * aerosol_procs(ia_m3,is0,ik,ic)
1063  enddo !is(1:n_siz(ic))
1064  enddo !ik(1:n_kap(ic))
1065  enddo !ic(1:n_ctg)
1066 
1067  return
real(rp), public const_mvap
mass weight (water vapor) [g/mol]
Definition: scale_const.F90:64
subroutine, public log(type, message)
Definition: dc_log.f90:133
module CONSTANT
Definition: scale_const.F90:14
Here is the call graph for this function:
Here is the caller graph for this function:

◆ diag_ds()

subroutine mod_mkinit::diag_ds ( real(rp)  m0,
real(rp)  m2,
real(rp)  m3,
real(rp)  dg,
real(rp)  sg,
real(rp)  dm2 
)

Definition at line 1072 of file mod_mkinit.f90.

References scale_const::const_pi, faero(), scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, dc_log::log(), scale_process::prc_mpistop(), scale_tracer::qa, scale_tracer::qqa, and mod_atmos_vars::qtrc.

Referenced by aerosol_activation_kajino13().

1072  implicit none
1073  real(RP) :: m0,m2,m3,dg,sg,m3_bar,m2_bar
1074  real(RP) :: m2_new,m2_old,dm2
1075  real(RP), parameter :: sgmax=2.5_rp
1076  real(RP), parameter :: rk1=2._rp
1077  real(RP), parameter :: rk2=3._rp
1078  real(RP), parameter :: ratio =rk1/rk2
1079  real(RP), parameter :: rk1_hat=1._rp/(ratio*(rk2-rk1))
1080  real(RP), parameter :: rk2_hat=ratio/(rk1-rk2)
1081 
1082  dm2=0._rp
1083 
1084  if (m0 <= 0._rp .or. m2 <= 0._rp .or. m3 <= 0._rp) then
1085  m0=0._rp
1086  m2=0._rp
1087  m3=0._rp
1088  dg=-1._rp
1089  sg=-1._rp
1090  return
1091  endif
1092 
1093  m2_old = m2
1094  m3_bar = m3/m0
1095  m2_bar = m2/m0
1096  dg = m2_bar**rk1_hat*m3_bar**rk2_hat
1097 
1098  if (m2_bar/m3_bar**ratio < 1._rp) then !stdev > 1.
1099  sg = exp(sqrt(2._rp/(rk1*(rk1-rk2)) &
1100  * log(m2_bar/m3_bar**ratio) ))
1101  endif
1102 
1103  if (sg > sgmax) then
1104  print *,'sg=',sg
1105  sg = sgmax
1106  call diag_d2(m0,m3,sg, & !i
1107  m2,dg ) !o
1108  ! print *,'warning! modified as sg exceeded sgmax (diag_ds)'
1109  endif
1110 
1111  if (m2_bar/m3_bar**ratio >= 1._rp) then !stdev < 1.
1112  sg = 1._rp
1113  call diag_d2(m0,m3,sg, & !i
1114  m2,dg ) !o
1115  ! print *,'warning! modified as sg < 1. (diag_ds)'
1116  endif
1117 
1118  m2_new = m2
1119  dm2 = m2_old - m2_new !m2_pres - m2_diag
1120 
1121  return
subroutine, public log(type, message)
Definition: dc_log.f90:133
Here is the call graph for this function:
Here is the caller graph for this function:

◆ faero()

real(rp) function mod_mkinit::faero ( real(rp), intent(in)  f0,
real(rp), intent(in)  r0,
real(rp), intent(in)  x,
real(rp), intent(in)  alpha,
real(rp), intent(in)  rhoa 
)

Definition at line 1236 of file mod_mkinit.f90.

References scale_const::const_pi.

Referenced by diag_ds().

1236  use scale_const, only: &
1237  pi => const_pi
1238  implicit none
1239 
1240  real(RP), intent(in) :: x, f0, r0, alpha, rhoa
1241  real(RP) :: faero
1242  real(RP) :: rad
1243  !---------------------------------------------------------------------------
1244 
1245  rad = ( exp(x) * 3.0_rp / 4.0_rp / pi / rhoa )**(1.0_rp/3.0_rp)
1246 
1247  faero = f0 * (rad/r0)**(-alpha)
1248 
1249  return
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_pi
pi
Definition: scale_const.F90:34
Here is the caller graph for this function:

◆ flux_setup()

subroutine mod_mkinit::flux_setup ( )

flux setup

Definition at line 1255 of file mod_mkinit.f90.

References mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_rain, mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_snow, mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_downall, mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_lw_dn, mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_lw_up, mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_sw_dn, mod_atmos_phy_rd_vars::atmos_phy_rd_sflx_sw_up, mod_atmos_phy_rd_vars::atmos_phy_rd_toaflx_lw_dn, mod_atmos_phy_rd_vars::atmos_phy_rd_toaflx_lw_up, mod_atmos_phy_rd_vars::atmos_phy_rd_toaflx_sw_dn, mod_atmos_phy_rd_vars::atmos_phy_rd_toaflx_sw_up, scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, and scale_process::prc_mpistop().

Referenced by interporation_fact(), and read_sounding().

1255  use mod_atmos_phy_mp_vars, only: &
1256  sflx_rain => atmos_phy_mp_sflx_rain, &
1257  sflx_snow => atmos_phy_mp_sflx_snow
1258  use mod_atmos_phy_rd_vars, only: &
1259  sflx_lw_up => atmos_phy_rd_sflx_lw_up, &
1260  sflx_lw_dn => atmos_phy_rd_sflx_lw_dn, &
1261  sflx_sw_up => atmos_phy_rd_sflx_sw_up, &
1262  sflx_sw_dn => atmos_phy_rd_sflx_sw_dn, &
1263  toaflx_lw_up => atmos_phy_rd_toaflx_lw_up, &
1264  toaflx_lw_dn => atmos_phy_rd_toaflx_lw_dn, &
1265  toaflx_sw_up => atmos_phy_rd_toaflx_sw_up, &
1266  toaflx_sw_dn => atmos_phy_rd_toaflx_sw_dn, &
1267  sflx_rad_dn => atmos_phy_rd_sflx_downall
1268  implicit none
1269  ! Flux from Atmosphere
1270  real(RP) :: flx_rain = 0.0_rp ! surface rain flux [kg/m2/s]
1271  real(RP) :: flx_snow = 0.0_rp ! surface snow flux [kg/m2/s]
1272  real(RP) :: flx_lw_dn = 0.0_rp ! surface downwad long-wave radiation flux [J/m2/s]
1273  real(RP) :: flx_sw_dn = 0.0_rp ! surface downwad short-wave radiation flux [J/m2/s]
1274 
1275  namelist / param_mkinit_flux / &
1276  flx_rain, &
1277  flx_snow, &
1278  flx_lw_dn, &
1279  flx_sw_dn
1280 
1281  integer :: i, j
1282  integer :: ierr
1283  !---------------------------------------------------------------------------
1284 
1285  !--- read namelist
1286  rewind(io_fid_conf)
1287  read(io_fid_conf,nml=param_mkinit_flux,iostat=ierr)
1288  if( ierr < 0 ) then !--- missing
1289  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1290  elseif( ierr > 0 ) then !--- fatal error
1291  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_FLUX. Check!'
1292  call prc_mpistop
1293  endif
1294  if( io_lnml ) write(io_fid_log,nml=param_mkinit_flux)
1295 
1296  do j = js, je
1297  do i = is, ie
1298  sflx_rain(i,j) = flx_rain
1299  sflx_snow(i,j) = flx_snow
1300 
1301  sflx_lw_up(i,j) = 0.0_rp
1302  sflx_lw_dn(i,j) = flx_lw_dn
1303  sflx_sw_up(i,j) = 0.0_rp
1304  sflx_sw_dn(i,j) = flx_sw_dn
1305 
1306  toaflx_lw_up(i,j) = 0.0_rp
1307  toaflx_lw_dn(i,j) = 0.0_rp
1308  toaflx_sw_up(i,j) = 0.0_rp
1309  toaflx_sw_dn(i,j) = 0.0_rp
1310 
1311  sflx_rad_dn(i,j,1,1) = flx_sw_dn
1312  sflx_rad_dn(i,j,1,2) = 0.0_rp
1313  sflx_rad_dn(i,j,2,1) = flx_lw_dn
1314  sflx_rad_dn(i,j,2,2) = 0.0_rp
1315  enddo
1316  enddo
1317 
1318  return
subroutine, public prc_mpistop
Abort MPI.
module Atmosphere / Physics Cloud Microphysics
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_up
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_rain
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_sw_up
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_lw_dn
module Atmosphere / Physics Radiation
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_up
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_sflx_sw_dn
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_lw_dn
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_rd_sflx_downall
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_snow
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_sw_dn
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_lw_up
Here is the call graph for this function:
Here is the caller graph for this function:

◆ land_setup()

subroutine mod_mkinit::land_setup ( )

Land setup.

Definition at line 1324 of file mod_mkinit.f90.

References scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, mod_land_vars::land_sfc_albedo, mod_land_vars::land_sfc_temp, mod_land_vars::land_temp, mod_land_vars::land_water, and scale_process::prc_mpistop().

Referenced by interporation_fact().

1324  use mod_land_vars, only: &
1325  land_temp, &
1326  land_water, &
1327  land_sfc_temp, &
1329  implicit none
1330  ! Land state
1331  real(RP) :: lnd_temp ! soil temperature [K]
1332  real(RP) :: lnd_water = 0.15_rp ! soil moisture [m3/m3]
1333  real(RP) :: sfc_temp ! land skin temperature [K]
1334  real(RP) :: sfc_albedo_lw = 0.01_rp ! land surface albedo for LW [0-1]
1335  real(RP) :: sfc_albedo_sw = 0.20_rp ! land surface albedo for SW [0-1]
1336 
1337  integer :: i, j
1338  integer :: ierr
1339 
1340  namelist /param_mkinit_land/ &
1341  lnd_temp, &
1342  lnd_water, &
1343  sfc_temp, &
1344  sfc_albedo_lw, &
1345  sfc_albedo_sw
1346 
1347  lnd_temp = thetastd
1348  sfc_temp = thetastd
1349 
1350  !--- read namelist
1351  rewind(io_fid_conf)
1352  read(io_fid_conf,nml=param_mkinit_land,iostat=ierr)
1353  if( ierr < 0 ) then !--- missing
1354  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1355  elseif( ierr > 0 ) then !--- fatal error
1356  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_LAND. Check!'
1357  call prc_mpistop
1358  endif
1359  if( io_lnml ) write(io_fid_log,nml=param_mkinit_land)
1360 
1361  do j = js, je
1362  do i = is, ie
1363  land_temp(:,i,j) = lnd_temp
1364  land_water(:,i,j) = lnd_water
1365  land_sfc_temp(i,j) = sfc_temp
1366  land_sfc_albedo(i,j,i_lw) = sfc_albedo_lw
1367  land_sfc_albedo(i,j,i_sw) = sfc_albedo_sw
1368  end do
1369  end do
1370 
1371  return
subroutine, public prc_mpistop
Abort MPI.
integer, parameter, public i_lw
real(rp), dimension(:,:), allocatable, public land_sfc_temp
land surface skin temperature [K]
integer, parameter, public i_sw
real(rp), dimension(:,:,:), allocatable, public land_temp
temperature of each soil layer [K]
real(rp), dimension(:,:,:), allocatable, public land_sfc_albedo
land surface albedo [0-1]
module LAND Variables
real(rp), dimension(:,:,:), allocatable, public land_water
moisture of each soil layer [m3/m3]
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ocean_setup()

subroutine mod_mkinit::ocean_setup ( )

Ocean setup.

Definition at line 1377 of file mod_mkinit.f90.

References scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, mod_ocean_vars::ocean_sfc_albedo, mod_ocean_vars::ocean_sfc_temp, mod_ocean_vars::ocean_sfc_z0e, mod_ocean_vars::ocean_sfc_z0h, mod_ocean_vars::ocean_sfc_z0m, mod_ocean_vars::ocean_temp, and scale_process::prc_mpistop().

Referenced by interporation_fact(), and read_sounding().

1377  use mod_ocean_vars, only: &
1378  ocean_temp, &
1379  ocean_sfc_temp, &
1380  ocean_sfc_albedo, &
1381  ocean_sfc_z0m, &
1382  ocean_sfc_z0h, &
1384  implicit none
1385  ! Ocean state
1386  real(RP) :: ocn_temp ! ocean temperature [K]
1387  real(RP) :: sfc_temp ! ocean skin temperature [K]
1388  real(RP) :: sfc_albedo_lw = 0.04_rp ! ocean surface albedo for LW [0-1]
1389  real(RP) :: sfc_albedo_sw = 0.05_rp ! ocean surface albedo for SW [0-1]
1390  real(RP) :: sfc_z0m = 1.0e-4_rp ! ocean surface roughness length (momentum) [m]
1391  real(RP) :: sfc_z0h = 1.0e-4_rp ! ocean surface roughness length (heat) [m]
1392  real(RP) :: sfc_z0e = 1.0e-4_rp ! ocean surface roughness length (vapor) [m]
1393 
1394  integer :: i, j
1395  integer :: ierr
1396 
1397  namelist /param_mkinit_ocean/ &
1398  ocn_temp, &
1399  sfc_temp, &
1400  sfc_albedo_lw, &
1401  sfc_albedo_sw, &
1402  sfc_z0m, &
1403  sfc_z0h, &
1404  sfc_z0e
1405 
1406  ocn_temp = thetastd
1407  sfc_temp = thetastd
1408 
1409  !--- read namelist
1410  rewind(io_fid_conf)
1411  read(io_fid_conf,nml=param_mkinit_ocean,iostat=ierr)
1412  if( ierr < 0 ) then !--- missing
1413  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1414  elseif( ierr > 0 ) then !--- fatal error
1415  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_OCEAN. Check!'
1416  call prc_mpistop
1417  endif
1418  if( io_lnml ) write(io_fid_log,nml=param_mkinit_ocean)
1419 
1420 
1421  do j = js, je
1422  do i = is, ie
1423  ocean_temp(i,j) = ocn_temp
1424  ocean_sfc_temp(i,j) = sfc_temp
1425  ocean_sfc_albedo(i,j,i_lw) = sfc_albedo_lw
1426  ocean_sfc_albedo(i,j,i_sw) = sfc_albedo_sw
1427  ocean_sfc_z0m(i,j) = sfc_z0m
1428  ocean_sfc_z0h(i,j) = sfc_z0h
1429  ocean_sfc_z0e(i,j) = sfc_z0e
1430  end do
1431  end do
1432 
1433  return
real(rp), dimension(:,:,:), allocatable, public ocean_sfc_albedo
ocean surface albedo [0-1]
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0e
ocean surface roughness length for vapor [m]
real(rp), dimension(:,:), allocatable, public ocean_temp
temperature at uppermost ocean layer [K]
integer, parameter, public i_lw
integer, parameter, public i_sw
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0m
ocean surface roughness length for momentum [m]
real(rp), dimension(:,:), allocatable, public ocean_sfc_temp
ocean surface skin temperature [K]
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0h
ocean surface roughness length for heat [m]
module OCEAN Variables
Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_setup()

subroutine mod_mkinit::urban_setup ( )

Urban setup.

Definition at line 1439 of file mod_mkinit.f90.

References scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_process::prc_mpistop(), mod_urban_vars::urban_qc, mod_urban_vars::urban_rainb, mod_urban_vars::urban_raing, mod_urban_vars::urban_rainr, mod_urban_vars::urban_roff, mod_urban_vars::urban_sfc_albedo, mod_urban_vars::urban_sfc_temp, mod_urban_vars::urban_tb, mod_urban_vars::urban_tbl, mod_urban_vars::urban_tc, mod_urban_vars::urban_tg, mod_urban_vars::urban_tgl, mod_urban_vars::urban_tr, mod_urban_vars::urban_trl, and mod_urban_vars::urban_uc.

Referenced by interporation_fact().

1439  use mod_urban_vars, only: &
1440  urban_tr, &
1441  urban_tb, &
1442  urban_tg, &
1443  urban_tc, &
1444  urban_qc, &
1445  urban_uc, &
1446  urban_trl, &
1447  urban_tbl, &
1448  urban_tgl, &
1449  urban_rainr, &
1450  urban_rainb, &
1451  urban_raing, &
1452  urban_roff, &
1453  urban_sfc_temp, &
1455  implicit none
1456  ! urban state
1457  real(RP) :: urb_roof_temp ! Surface temperature of roof [K]
1458  real(RP) :: urb_bldg_temp ! Surface temperature of building [K
1459  real(RP) :: urb_grnd_temp ! Surface temperature of ground [K]
1460  real(RP) :: urb_cnpy_temp ! Diagnostic canopy air temperature
1461  real(RP) :: urb_cnpy_hmdt = 0.0_rp ! Diagnostic canopy humidity [-]
1462  real(RP) :: urb_cnpy_wind = 0.0_rp ! Diagnostic canopy wind [m/s]
1463  real(RP) :: urb_roof_layer_temp ! temperature in layer of roof [K]
1464  real(RP) :: urb_bldg_layer_temp ! temperature in layer of building [
1465  real(RP) :: urb_grnd_layer_temp ! temperature in layer of ground [K]
1466  real(RP) :: urb_roof_rain = 0.0_rp ! temperature in layer of roof [K]
1467  real(RP) :: urb_bldg_rain = 0.0_rp ! temperature in layer of building [
1468  real(RP) :: urb_grnd_rain = 0.0_rp ! temperature in layer of ground [K]
1469  real(RP) :: urb_runoff = 0.0_rp ! temperature in layer of ground [K]
1470  real(RP) :: urb_sfc_temp ! Grid average of surface temperature [K]
1471  real(RP) :: urb_alb_lw = 0.0_rp ! Grid average of surface albedo for LW [0-1]
1472  real(RP) :: urb_alb_sw = 0.0_rp ! Grid average of surface albedo for SW [0-1]
1473 
1474  integer :: i, j
1475  integer :: ierr
1476 
1477  namelist /param_mkinit_urban/ &
1478  urb_roof_temp, &
1479  urb_bldg_temp, &
1480  urb_grnd_temp, &
1481  urb_cnpy_temp, &
1482  urb_cnpy_hmdt, &
1483  urb_cnpy_wind, &
1484  urb_roof_layer_temp, &
1485  urb_bldg_layer_temp, &
1486  urb_grnd_layer_temp, &
1487  urb_roof_rain, &
1488  urb_bldg_rain, &
1489  urb_grnd_rain, &
1490  urb_runoff, &
1491  urb_sfc_temp, &
1492  urb_alb_lw, &
1493  urb_alb_sw
1494 
1495  urb_roof_temp = thetastd
1496  urb_bldg_temp = thetastd
1497  urb_grnd_temp = thetastd
1498  urb_cnpy_temp = thetastd
1499  urb_roof_layer_temp = thetastd
1500  urb_bldg_layer_temp = thetastd
1501  urb_grnd_layer_temp = thetastd
1502 
1503  urb_sfc_temp = thetastd
1504 
1505  !--- read namelist
1506  rewind(io_fid_conf)
1507  read(io_fid_conf,nml=param_mkinit_urban,iostat=ierr)
1508  if( ierr < 0 ) then !--- missing
1509  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1510  elseif( ierr > 0 ) then !--- fatal error
1511  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_URBAN. Check!'
1512  call prc_mpistop
1513  endif
1514  if( io_lnml ) write(io_fid_log,nml=param_mkinit_urban)
1515 
1516 
1517  do j = js, je
1518  do i = is, ie
1519  urban_tr(i,j) = urb_roof_temp
1520  urban_tb(i,j) = urb_bldg_temp
1521  urban_tg(i,j) = urb_grnd_temp
1522  urban_tc(i,j) = urb_cnpy_temp
1523  urban_qc(i,j) = urb_cnpy_hmdt
1524  urban_uc(i,j) = urb_cnpy_wind
1525  urban_trl(:,i,j) = urb_roof_layer_temp
1526  urban_tbl(:,i,j) = urb_bldg_layer_temp
1527  urban_tgl(:,i,j) = urb_grnd_layer_temp
1528  urban_rainr(i,j) = urb_roof_rain
1529  urban_rainb(i,j) = urb_bldg_rain
1530  urban_raing(i,j) = urb_grnd_rain
1531  urban_roff(i,j) = urb_runoff
1532  urban_sfc_temp(i,j) = urb_sfc_temp
1533  urban_sfc_albedo(i,j,i_lw) = urb_alb_lw
1534  urban_sfc_albedo(i,j,i_sw) = urb_alb_sw
1535  end do
1536  end do
1537 
1538  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:), allocatable, public urban_qc
real(rp), dimension(:,:), allocatable, public urban_tb
module URBAN Variables
real(rp), dimension(:,:), allocatable, public urban_raing
real(rp), dimension(:,:), allocatable, public urban_uc
integer, parameter, public i_lw
integer, parameter, public i_sw
real(rp), dimension(:,:), allocatable, public urban_tr
real(rp), dimension(:,:,:), allocatable, public urban_tgl
real(rp), dimension(:,:), allocatable, public urban_roff
real(rp), dimension(:,:,:), allocatable, public urban_sfc_albedo
real(rp), dimension(:,:), allocatable, public urban_tc
real(rp), dimension(:,:), allocatable, public urban_rainr
real(rp), dimension(:,:), allocatable, public urban_tg
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
real(rp), dimension(:,:,:), allocatable, public urban_trl
real(rp), dimension(:,:,:), allocatable, public urban_tbl
real(rp), dimension(:,:), allocatable, public urban_rainb
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_sounding()

subroutine mod_mkinit::read_sounding ( real(rp), dimension(ka), intent(out)  DENS,
real(rp), dimension(ka), intent(out)  VELX,
real(rp), dimension(ka), intent(out)  VELY,
real(rp), dimension(ka), intent(out)  POTT,
real(rp), dimension (ka), intent(out)  QV 
)

Read sounding data from file.

Definition at line 1545 of file mod_mkinit.f90.

References scale_tracer::aetracer_type, scale_tracer::aq_name, mod_atmos_vars::dens, gtool_file::filegetshape(), flux_setup(), scale_grid::grid_cx, scale_grid::grid_cy, scale_grid::grid_cz, scale_grid::grid_fx, scale_grid::grid_fy, scale_grid::grid_fz, scale_tracer::i_nc, scale_tracer::i_qc, scale_tracer::i_qv, scale_grid_index::ia, scale_grid_index::ie, interporation_fact(), scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, dc_log::log(), mod_atmos_vars::momx, mod_atmos_vars::momy, mod_atmos_vars::momz, ocean_setup(), scale_process::prc_mpistop(), scale_tracer::qa, scale_tracer::qqa, mod_atmos_vars::qtrc, scale_random::random_get(), scale_grid_real::real_cz, scale_grid_real::real_fz, rect_setup(), mod_atmos_vars::rhot, and scale_tracer::tracer_type.

Referenced by interporation_fact().

1545  implicit none
1546  real(RP), intent(out) :: dens(ka)
1547  real(RP), intent(out) :: velx(ka)
1548  real(RP), intent(out) :: vely(ka)
1549  real(RP), intent(out) :: pott(ka)
1550  real(RP), intent(out) :: qv (ka)
1551 
1552  real(RP) :: temp(ka)
1553  real(RP) :: pres(ka)
1554  real(RP) :: qc (ka)
1555 
1556  character(len=H_LONG) :: env_in_sounding_file = ''
1557 
1558  integer, parameter :: exp_klim = 100
1559  integer :: exp_kmax
1560 
1561  real(RP) :: sfc_theta ! surface potential temperature [K]
1562  real(RP) :: sfc_pres ! surface pressure [hPa]
1563  real(RP) :: sfc_qv ! surface watervapor [g/kg]
1564 
1565  real(RP) :: exp_z (exp_klim+1) ! height [m]
1566  real(RP) :: exp_pott(exp_klim+1) ! potential temperature [K]
1567  real(RP) :: exp_qv (exp_klim+1) ! water vapor [g/kg]
1568  real(RP) :: exp_u (exp_klim+1) ! velocity u [m/s]
1569  real(RP) :: exp_v (exp_klim+1) ! velocity v [m/s]
1570 
1571  real(RP) :: fact1, fact2
1572  integer :: k, kref
1573  integer :: fid
1574  integer :: ierr
1575 
1576  namelist /param_mkinit_sounding/ &
1577  env_in_sounding_file
1578 
1579  !--- read namelist
1580  rewind(io_fid_conf)
1581  read(io_fid_conf,nml=param_mkinit_sounding,iostat=ierr)
1582 
1583  if( ierr < 0 ) then !--- missing
1584  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1585  elseif( ierr > 0 ) then !--- fatal error
1586  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_SOUNDING. Check!'
1587  call prc_mpistop
1588  endif
1589  if( io_lnml ) write(io_fid_log,nml=param_mkinit_sounding)
1590 
1591  !--- prepare sounding profile
1592  if( io_l ) write(io_fid_log,*) '+++ Input sounding file:', trim(env_in_sounding_file)
1593  fid = io_get_available_fid()
1594  open( fid, &
1595  file = trim(env_in_sounding_file), &
1596  form = 'formatted', &
1597  status = 'old', &
1598  iostat = ierr )
1599 
1600  if ( ierr /= 0 ) then
1601  if( io_l ) write(*,*) 'xxx Input file not found!'
1602  endif
1603 
1604  !--- read sounding file till end
1605  read(fid,*) sfc_pres, sfc_theta, sfc_qv
1606 
1607  if( io_l ) write(io_fid_log,*) '+++ Surface pressure [hPa]', sfc_pres
1608  if( io_l ) write(io_fid_log,*) '+++ Surface pot. temp [K]', sfc_theta
1609  if( io_l ) write(io_fid_log,*) '+++ Surface water vapor [g/kg]', sfc_qv
1610 
1611  do k = 2, exp_klim
1612  read(fid,*,iostat=ierr) exp_z(k), exp_pott(k), exp_qv(k), exp_u(k), exp_v(k)
1613  if ( ierr /= 0 ) exit
1614  enddo
1615 
1616  exp_kmax = k - 1
1617  close(fid)
1618 
1619  ! Boundary
1620  exp_z(1) = 0.0_rp
1621  exp_pott(1) = sfc_theta
1622  exp_qv(1) = sfc_qv
1623  exp_u(1) = exp_u(2)
1624  exp_v(1) = exp_v(2)
1625  exp_z(exp_kmax+1) = 100.e3_rp
1626  exp_pott(exp_kmax+1) = exp_pott(exp_kmax)
1627  exp_qv(exp_kmax+1) = exp_qv(exp_kmax)
1628  exp_u(exp_kmax+1) = exp_u(exp_kmax)
1629  exp_v(exp_kmax+1) = exp_v(exp_kmax)
1630 
1631  do k = 1, exp_kmax+1
1632  exp_qv(k) = exp_qv(k) * 1.e-3_rp ! [g/kg]->[kg/kg]
1633  enddo
1634 
1635  ! calc in dry condition
1636  pres_sfc = sfc_pres * 1.e2_rp ! [hPa]->[Pa]
1637  pott_sfc = sfc_theta
1638  qv_sfc = sfc_qv * 1.e-3_rp ! [g/kg]->[kg/kg]
1639  qc_sfc = 0.0_rp
1640 
1641  !--- linear interpolate to model grid
1642  do k = ks, ke
1643  qc(k) = 0.0_rp
1644 
1645  do kref = 2, exp_kmax+1
1646  if ( grid_cz(k) > exp_z(kref-1) &
1647  .AND. grid_cz(k) <= exp_z(kref ) ) then
1648 
1649  fact1 = ( exp_z(kref) - grid_cz(k) ) / ( exp_z(kref)-exp_z(kref-1) )
1650  fact2 = ( grid_cz(k) - exp_z(kref-1) ) / ( exp_z(kref)-exp_z(kref-1) )
1651 
1652  pott(k) = exp_pott(kref-1) * fact1 &
1653  + exp_pott(kref ) * fact2
1654  qv(k) = exp_qv(kref-1) * fact1 &
1655  + exp_qv(kref ) * fact2
1656  velx(k) = exp_u(kref-1) * fact1 &
1657  + exp_u(kref ) * fact2
1658  vely(k) = exp_v(kref-1) * fact1 &
1659  + exp_v(kref ) * fact2
1660  endif
1661  enddo
1662  enddo
1663 
1664  ! make density & pressure profile in moist condition
1665  call hydrostatic_buildrho( dens(:), & ! [OUT]
1666  temp(:), & ! [OUT]
1667  pres(:), & ! [OUT]
1668  pott(:), & ! [IN]
1669  qv(:), & ! [IN]
1670  qc(:), & ! [IN]
1671  temp_sfc(1,1,1), & ! [OUT]
1672  pres_sfc(1,1,1), & ! [IN]
1673  pott_sfc(1,1,1), & ! [IN]
1674  qv_sfc(1,1,1), & ! [IN]
1675  qc_sfc(1,1,1) ) ! [IN]
1676 
1677  return
subroutine, public prc_mpistop
Abort MPI.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ interporation_fact()

subroutine mod_mkinit::interporation_fact ( real(rp), intent(out)  fact0,
real(rp), intent(out)  fact1,
integer, intent(out)  idx0,
integer, intent(out)  idx1,
real(rp), intent(in)  x,
real(rp), dimension(nx), intent(in)  x_org,
integer, intent(in)  nx,
logical, intent(in)  loop 
)

Definition at line 4586 of file mod_mkinit.f90.

References scale_tracer::aetracer_type, mod_atmos_vars::dens, flux_setup(), scale_grid::grid_cx, scale_grid::grid_cxg, scale_grid::grid_cz, scale_tracer::i_qv, scale_grid_index::ia, scale_grid_index::ie, scale_grid_index::imax, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, land_setup(), scale_landuse::landuse_calc_fact(), scale_landuse::landuse_frac_land, scale_landuse::landuse_frac_urban, mod_atmos_vars::momx, mod_atmos_vars::momy, mod_atmos_vars::momz, ocean_setup(), scale_process::prc_mpistop(), scale_rm_process::prc_num_x, mod_atmos_vars::qtrc, scale_random::random_get(), read_sounding(), mod_realinput::realinput_atmos(), mod_realinput::realinput_surface(), mod_atmos_vars::rhot, and urban_setup().

Referenced by read_sounding().

4586  implicit none
4587 
4588  real(RP), intent(out) :: fact0
4589  real(RP), intent(out) :: fact1
4590  integer, intent(out) :: idx0
4591  integer, intent(out) :: idx1
4592  real(RP), intent(in) :: x
4593  integer, intent(in) :: nx
4594  real(RP), intent(in) :: x_org(nx)
4595  logical, intent(in) :: loop
4596 
4597  real(RP) :: xwork
4598  integer :: i
4599 
4600  if ( x < x_org(1) ) then
4601  if ( loop ) then
4602  xwork = x_org(1) - ( x_org(2) - x_org(1) )**2 / ( x_org(3) - x_org(2) )
4603  fact0 = ( x_org(1) - x ) / ( x_org(1) - xwork )
4604  fact1 = ( x - xwork ) / ( x_org(1) - xwork )
4605  idx0 = nx
4606  idx1 = 1
4607  else
4608  fact0 = ( x_org(2) - x ) / ( x_org(2) - x_org(1) )
4609  fact1 = ( x - x_org(1) ) / ( x_org(2) - x_org(1) )
4610  idx0 = 1
4611  idx1 = 2
4612  end if
4613  else if ( x > x_org(nx) ) then
4614  if ( loop ) then
4615  xwork = x_org(nx) + ( x_org(nx) - x_org(nx-1) )**2 / ( x_org(nx-1) - x_org(nx-2) )
4616  fact0 = ( xwork - x ) / ( xwork - x_org(nx) )
4617  fact1 = ( x - x_org(nx) ) / ( xwork - x_org(nx) )
4618  idx0 = nx
4619  idx1 = 1
4620  else
4621  fact0 = ( x_org(nx) - x ) / ( x_org(nx) - x_org(nx-1) )
4622  fact1 = ( x - x_org(nx-1) ) / ( x_org(nx) - x_org(nx-1) )
4623  idx0 = nx-1
4624  idx1 = nx
4625  end if
4626  else
4627  do i = 2, nx
4628  if ( x <= x_org(i) ) then
4629  fact0 = ( x_org(i) - x ) / ( x_org(i) - x_org(i-1) )
4630  fact1 = ( x - x_org(i-1) ) / ( x_org(i) - x_org(i-1) )
4631  idx0 = i-1
4632  idx1 = i
4633  exit
4634  end if
4635  end do
4636  end if
4637 
4638  return
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ mkinit_type

integer, public mod_mkinit::mkinit_type = -1

Definition at line 102 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

102  integer, public :: mkinit_type = -1

◆ i_ignore

integer, parameter, public mod_mkinit::i_ignore = 0

Definition at line 103 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

103  integer, public, parameter :: i_ignore = 0

◆ i_planestate

integer, parameter, public mod_mkinit::i_planestate = 1

Definition at line 105 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

105  integer, public, parameter :: i_planestate = 1

◆ i_tracerbubble

integer, parameter, public mod_mkinit::i_tracerbubble = 2

Definition at line 106 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

106  integer, public, parameter :: i_tracerbubble = 2

◆ i_coldbubble

integer, parameter, public mod_mkinit::i_coldbubble = 3

Definition at line 107 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

107  integer, public, parameter :: i_coldbubble = 3

◆ i_lambwave

integer, parameter, public mod_mkinit::i_lambwave = 4

Definition at line 109 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

109  integer, public, parameter :: i_lambwave = 4

◆ i_gravitywave

integer, parameter, public mod_mkinit::i_gravitywave = 5

Definition at line 110 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

110  integer, public, parameter :: i_gravitywave = 5

◆ i_khwave

integer, parameter, public mod_mkinit::i_khwave = 6

Definition at line 111 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

111  integer, public, parameter :: i_khwave = 6

◆ i_turbulence

integer, parameter, public mod_mkinit::i_turbulence = 7

Definition at line 112 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

112  integer, public, parameter :: i_turbulence = 7

◆ i_mountainwave

integer, parameter, public mod_mkinit::i_mountainwave = 8

Definition at line 113 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

113  integer, public, parameter :: i_mountainwave = 8

◆ i_warmbubble

integer, parameter, public mod_mkinit::i_warmbubble = 9

Definition at line 115 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

115  integer, public, parameter :: i_warmbubble = 9

◆ i_supercell

integer, parameter, public mod_mkinit::i_supercell = 10

Definition at line 116 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

116  integer, public, parameter :: i_supercell = 10

◆ i_squallline

integer, parameter, public mod_mkinit::i_squallline = 11

Definition at line 117 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

117  integer, public, parameter :: i_squallline = 11

◆ i_wk1982

integer, parameter, public mod_mkinit::i_wk1982 = 12

Definition at line 118 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

118  integer, public, parameter :: i_wk1982 = 12

◆ i_dycoms2_rf01

integer, parameter, public mod_mkinit::i_dycoms2_rf01 = 13

Definition at line 119 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

119  integer, public, parameter :: i_dycoms2_rf01 = 13

◆ i_dycoms2_rf02

integer, parameter, public mod_mkinit::i_dycoms2_rf02 = 14

Definition at line 120 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

120  integer, public, parameter :: i_dycoms2_rf02 = 14

◆ i_rico

integer, parameter, public mod_mkinit::i_rico = 15

Definition at line 121 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

121  integer, public, parameter :: i_rico = 15

◆ i_interporation

integer, parameter, public mod_mkinit::i_interporation = 16

Definition at line 123 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

123  integer, public, parameter :: i_interporation = 16

◆ i_landcouple

integer, parameter, public mod_mkinit::i_landcouple = 17

Definition at line 125 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

125  integer, public, parameter :: i_landcouple = 17

◆ i_oceancouple

integer, parameter, public mod_mkinit::i_oceancouple = 18

Definition at line 126 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

126  integer, public, parameter :: i_oceancouple = 18

◆ i_urbancouple

integer, parameter, public mod_mkinit::i_urbancouple = 19

Definition at line 127 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

127  integer, public, parameter :: i_urbancouple = 19

◆ i_triplecouple

integer, parameter, public mod_mkinit::i_triplecouple = 20

Definition at line 128 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

128  integer, public, parameter :: i_triplecouple = 20

◆ i_bubblecouple

integer, parameter, public mod_mkinit::i_bubblecouple = 21

Definition at line 129 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

129  integer, public, parameter :: i_bubblecouple = 21

◆ i_seabreeze

integer, parameter, public mod_mkinit::i_seabreeze = 22

Definition at line 131 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

131  integer, public, parameter :: i_seabreeze = 22

◆ i_heatisland

integer, parameter, public mod_mkinit::i_heatisland = 23

Definition at line 132 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

132  integer, public, parameter :: i_heatisland = 23

◆ i_dycoms2_rf02_dns

integer, parameter, public mod_mkinit::i_dycoms2_rf02_dns = 24

Definition at line 134 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

134  integer, public, parameter :: i_dycoms2_rf02_dns = 24

◆ i_real

integer, parameter, public mod_mkinit::i_real = 25

Definition at line 136 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

136  integer, public, parameter :: i_real = 25

◆ i_grayzone

integer, parameter, public mod_mkinit::i_grayzone = 26

Definition at line 138 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

138  integer, public, parameter :: i_grayzone = 26

◆ i_boxaero

integer, parameter, public mod_mkinit::i_boxaero = 27

Definition at line 139 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

139  integer, public, parameter :: i_boxaero = 27

◆ i_warmbubbleaero

integer, parameter, public mod_mkinit::i_warmbubbleaero = 28

Definition at line 140 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

140  integer, public, parameter :: i_warmbubbleaero = 28

◆ i_cavityflow

integer, parameter, public mod_mkinit::i_cavityflow = 29

Definition at line 142 of file mod_mkinit.f90.

Referenced by mkinit(), and mkinit_setup().

142  integer, public, parameter :: i_cavityflow = 29