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.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)
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_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, 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  implicit none
372 
373  integer :: iq
374  !---------------------------------------------------------------------------
375 
376  if ( mkinit_type == i_ignore ) then
377  if( io_l ) write(io_fid_log,*)
378  if( io_l ) write(io_fid_log,*) '++++++ SKIP MAKING INITIAL DATA ++++++'
379  else
380  if( io_l ) write(io_fid_log,*)
381  if( io_l ) write(io_fid_log,*) '++++++ START MAKING INITIAL DATA ++++++'
382 
383  !--- Initialize variables
384 #ifndef DRY
385  do iq = 2, qa
386  qtrc(:,:,:,iq) = 0.0_rp
387  enddo
388 #endif
389 
390  pres(:,:,:) = const_undef8
391  temp(:,:,:) = const_undef8
392  pott(:,:,:) = const_undef8
393  qsat(:,:,:) = const_undef8
394  qv(:,:,:) = const_undef8
395  qc(:,:,:) = const_undef8
396  velx(:,:,:) = const_undef8
397  vely(:,:,:) = const_undef8
398 
399  rndm(:,:,:) = const_undef8
400 
401  pres_sfc(:,:,:) = const_undef8
402  temp_sfc(:,:,:) = const_undef8
403  pott_sfc(:,:,:) = const_undef8
404  qsat_sfc(:,:,:) = const_undef8
405  qv_sfc(:,:,:) = const_undef8
406  qc_sfc(:,:,:) = const_undef8
407 
408  if ( tracer_type == 'SUZUKI10' ) then
409  if( io_l ) write(io_fid_log,*) '*** Aerosols for SBM are included ***'
410  call sbmaero_setup
411  endif
412 
413  select case(mkinit_type)
414  case(i_planestate)
415  call mkinit_planestate
416  case(i_tracerbubble)
417  call mkinit_tracerbubble
418  case(i_coldbubble)
419  call mkinit_coldbubble
420  case(i_lambwave)
421  call mkinit_lambwave
422  case(i_gravitywave)
423  call mkinit_gravitywave
424  case(i_khwave)
425  call mkinit_khwave
426  case(i_turbulence)
427  call mkinit_turbulence
428  case(i_mountainwave)
429  call mkinit_mountainwave
430  case(i_warmbubble)
431  call mkinit_warmbubble
432  case(i_supercell)
433  call mkinit_supercell
434  case(i_squallline)
435  call mkinit_squallline
436  case(i_wk1982)
437  call mkinit_wk1982
438  case(i_dycoms2_rf01)
439  call mkinit_dycoms2_rf01
440  case(i_dycoms2_rf02)
441  call mkinit_dycoms2_rf02
442  case(i_rico)
443  call mkinit_rico
444  case(i_interporation)
445  call mkinit_interporation
446  case(i_oceancouple)
447  call mkinit_planestate
448  call mkinit_oceancouple
449  case(i_landcouple)
450  call mkinit_planestate
451  call mkinit_landcouple
452  case(i_urbancouple)
453  call mkinit_planestate
454  call mkinit_landcouple ! tentative
455  call mkinit_urbancouple
456  case(i_triplecouple)
457  call mkinit_planestate
458  call mkinit_oceancouple
459  call mkinit_landcouple
460  call mkinit_urbancouple
461  case(i_bubblecouple)
462  call mkinit_planestate
463  call mkinit_warmbubble
464  call mkinit_oceancouple
465  call mkinit_landcouple
466  call mkinit_urbancouple
467  case(i_seabreeze)
468  call mkinit_planestate
469  call mkinit_seabreeze
470  case(i_heatisland)
471  call mkinit_planestate
472  call mkinit_heatisland
473  case(i_dycoms2_rf02_dns)
474  call mkinit_dycoms2_rf02_dns
475  case(i_real)
476  call mkinit_real
477  case(i_grayzone)
478  call mkinit_grayzone
479  case(i_boxaero)
480  call mkinit_boxaero
481  case(i_warmbubbleaero)
482  call mkinit_warmbubbleaero
483  case(i_cavityflow)
484  call mkinit_cavityflow
485  case default
486  write(*,*) ' xxx Unsupported TYPE:', mkinit_type
487  call prc_mpistop
488  endselect
489 
490  if( io_l ) write(io_fid_log,*) '++++++ END MAKING INITIAL DATA ++++++'
491 
492  ! setup surface condition
493  call ocean_surface_set( countup = .false. )
494  call land_surface_set ( countup = .false. )
495  call urban_surface_set( countup = .false. )
496  call atmos_surface_get
497 
498  ! output boundary file
499  call landuse_write
500 
501  ! output restart file
502  if( atmos_sw_restart ) call atmos_vars_restart_write
503  if( ocean_sw_restart ) call ocean_vars_restart_write
504  if( land_sw_restart ) call land_vars_restart_write
505  if( urban_sw_restart ) call urban_vars_restart_write
506  endif
507 
508  return
subroutine, public prc_mpistop
Abort MPI.
logical, public land_restart_output
output restart file?
module URBAN driver
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine, public ocean_vars_restart_write
Write ocean restart.
module ATMOSPHERIC Variables
subroutine, public atmos_vars_restart_write
Write restart of atmospheric variables.
integer, public qa
module URBAN Variables
character(len=h_short), public tracer_type
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 CONSTANT
Definition: scale_const.F90:14
logical, public urban_restart_output
output restart file?
real(dp), parameter, public const_undef8
undefined value (REAL8)
Definition: scale_const.F90:42
logical, public ocean_restart_output
output restart file?
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
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 599 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().

599  use scale_const, only: &
601  implicit none
602 
603  ! Bubble
604  logical :: rct_eachnode = .false. ! Arrange rectangle at each node? [kg/kg]
605  real(RP) :: rct_cz = 2.e3_rp ! center location [m]: z
606  real(RP) :: rct_cx = 2.e3_rp ! center location [m]: x
607  real(RP) :: rct_cy = 2.e3_rp ! center location [m]: y
608  real(RP) :: rct_rz = 2.e3_rp ! rectangle z width [m]: z
609  real(RP) :: rct_rx = 2.e3_rp ! rectangle x width [m]: x
610  real(RP) :: rct_ry = 2.e3_rp ! rectangle y width [m]: y
611 
612  namelist / param_rect / &
613  rct_eachnode, &
614  rct_cz, &
615  rct_cx, &
616  rct_cy, &
617  rct_rz, &
618  rct_rx, &
619  rct_ry
620 
621  real(RP) :: cz_offset
622  real(RP) :: cx_offset
623  real(RP) :: cy_offset
624  real(RP) :: dist
625 
626  integer :: ierr
627  integer :: k, i, j
628  !---------------------------------------------------------------------------
629 
630  if( io_l ) write(io_fid_log,*)
631  if( io_l ) write(io_fid_log,*) '++++++ Module[mkinit rectangle] / Categ[preprocess] / Origin[SCALE-RM]'
632 
633  !--- read namelist
634  rewind(io_fid_conf)
635  read(io_fid_conf,nml=param_rect,iostat=ierr)
636 
637  if( ierr < 0 ) then !--- missing
638  if( io_l ) write(io_fid_log,*) 'xxx Not found namelist. Check!'
639  call prc_mpistop
640  elseif( ierr > 0 ) then !--- fatal error
641  write(*,*) 'xxx Not appropriate names in namelist PARAM_RECT. Check!'
642  call prc_mpistop
643  endif
644  if( io_lnml ) write(io_fid_log,nml=param_rect)
645 
646  rect(:,:,:) = const_undef8
647 
648  if ( rct_eachnode ) then
649  cz_offset = grid_cz(ks)
650  cx_offset = grid_cx(is)
651  cy_offset = grid_cy(js)
652  else
653  cz_offset = 0.0_rp
654  cx_offset = 0.0_rp
655  cy_offset = 0.0_rp
656  endif
657 
658  do j = 1, ja
659  do i = 1, ia
660  do k = ks, ke
661 
662  ! make tracer rectangle
663  dist = 2.0_rp * max( &
664  abs(grid_cz(k) - cz_offset - rct_cz)/rct_rz, &
665  abs(grid_cx(i) - cx_offset - rct_cx)/rct_rx, &
666  abs(grid_cy(j) - cy_offset - rct_cy)/rct_ry &
667  & )
668  if ( dist <= 1.0_rp ) then
669  rect(k,i,j) = 1.0_rp
670  else
671  rect(k,i,j) = 0.0_rp
672  end if
673  enddo
674  enddo
675  enddo
676 
677  return
integer, public is
start point of inner domain: x, local
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public ke
end point of inner domain: z, local
integer, public ia
of x whole cells (local, with HALO)
integer, public js
start point of inner domain: y, local
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
real(dp), parameter, public const_undef8
undefined value (REAL8)
Definition: scale_const.F90:42
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)
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 969 of file mod_mkinit.f90.

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

Referenced by rect_setup().

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

1059  implicit none
1060  real(RP) :: m0,m2,m3,dg,sg,m3_bar,m2_bar
1061  real(RP) :: m2_new,m2_old,dm2
1062  real(RP), parameter :: sgmax=2.5_rp
1063  real(RP), parameter :: rk1=2._rp
1064  real(RP), parameter :: rk2=3._rp
1065  real(RP), parameter :: ratio =rk1/rk2
1066  real(RP), parameter :: rk1_hat=1._rp/(ratio*(rk2-rk1))
1067  real(RP), parameter :: rk2_hat=ratio/(rk1-rk2)
1068 
1069  dm2=0._rp
1070 
1071  if (m0 <= 0._rp .or. m2 <= 0._rp .or. m3 <= 0._rp) then
1072  m0=0._rp
1073  m2=0._rp
1074  m3=0._rp
1075  dg=-1._rp
1076  sg=-1._rp
1077  return
1078  endif
1079 
1080  m2_old = m2
1081  m3_bar = m3/m0
1082  m2_bar = m2/m0
1083  dg = m2_bar**rk1_hat*m3_bar**rk2_hat
1084 
1085  if (m2_bar/m3_bar**ratio < 1._rp) then !stdev > 1.
1086  sg = exp(sqrt(2._rp/(rk1*(rk1-rk2)) &
1087  * log(m2_bar/m3_bar**ratio) ))
1088  endif
1089 
1090  if (sg > sgmax) then
1091  print *,'sg=',sg
1092  sg = sgmax
1093  call diag_d2(m0,m3,sg, & !i
1094  m2,dg ) !o
1095  ! print *,'warning! modified as sg exceeded sgmax (diag_ds)'
1096  endif
1097 
1098  if (m2_bar/m3_bar**ratio >= 1._rp) then !stdev < 1.
1099  sg = 1._rp
1100  call diag_d2(m0,m3,sg, & !i
1101  m2,dg ) !o
1102  ! print *,'warning! modified as sg < 1. (diag_ds)'
1103  endif
1104 
1105  m2_new = m2
1106  dm2 = m2_old - m2_new !m2_pres - m2_diag
1107 
1108  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 1223 of file mod_mkinit.f90.

References scale_const::const_pi.

Referenced by diag_ds().

1223  use scale_const, only: &
1224  pi => const_pi
1225  implicit none
1226 
1227  real(RP), intent(in) :: x, f0, r0, alpha, rhoa
1228  real(RP) :: faero
1229  real(RP) :: rad
1230  !---------------------------------------------------------------------------
1231 
1232  rad = ( exp(x) * 3.0_rp / 4.0_rp / pi / rhoa )**(1.0_rp/3.0_rp)
1233 
1234  faero = f0 * (rad/r0)**(-alpha)
1235 
1236  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 1242 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().

1242  use mod_atmos_phy_mp_vars, only: &
1243  sflx_rain => atmos_phy_mp_sflx_rain, &
1244  sflx_snow => atmos_phy_mp_sflx_snow
1245  use mod_atmos_phy_rd_vars, only: &
1246  sflx_lw_up => atmos_phy_rd_sflx_lw_up, &
1247  sflx_lw_dn => atmos_phy_rd_sflx_lw_dn, &
1248  sflx_sw_up => atmos_phy_rd_sflx_sw_up, &
1249  sflx_sw_dn => atmos_phy_rd_sflx_sw_dn, &
1250  toaflx_lw_up => atmos_phy_rd_toaflx_lw_up, &
1251  toaflx_lw_dn => atmos_phy_rd_toaflx_lw_dn, &
1252  toaflx_sw_up => atmos_phy_rd_toaflx_sw_up, &
1253  toaflx_sw_dn => atmos_phy_rd_toaflx_sw_dn, &
1254  sflx_rad_dn => atmos_phy_rd_sflx_downall
1255  implicit none
1256  ! Flux from Atmosphere
1257  real(RP) :: flx_rain = 0.0_rp ! surface rain flux [kg/m2/s]
1258  real(RP) :: flx_snow = 0.0_rp ! surface snow flux [kg/m2/s]
1259  real(RP) :: flx_lw_dn = 0.0_rp ! surface downwad long-wave radiation flux [J/m2/s]
1260  real(RP) :: flx_sw_dn = 0.0_rp ! surface downwad short-wave radiation flux [J/m2/s]
1261 
1262  namelist / param_mkinit_flux / &
1263  flx_rain, &
1264  flx_snow, &
1265  flx_lw_dn, &
1266  flx_sw_dn
1267 
1268  integer :: i, j
1269  integer :: ierr
1270  !---------------------------------------------------------------------------
1271 
1272  !--- read namelist
1273  rewind(io_fid_conf)
1274  read(io_fid_conf,nml=param_mkinit_flux,iostat=ierr)
1275  if( ierr < 0 ) then !--- missing
1276  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1277  elseif( ierr > 0 ) then !--- fatal error
1278  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_FLUX. Check!'
1279  call prc_mpistop
1280  endif
1281  if( io_lnml ) write(io_fid_log,nml=param_mkinit_flux)
1282 
1283  do j = js, je
1284  do i = is, ie
1285  sflx_rain(i,j) = flx_rain
1286  sflx_snow(i,j) = flx_snow
1287 
1288  sflx_lw_up(i,j) = 0.0_rp
1289  sflx_lw_dn(i,j) = flx_lw_dn
1290  sflx_sw_up(i,j) = 0.0_rp
1291  sflx_sw_dn(i,j) = flx_sw_dn
1292 
1293  toaflx_lw_up(i,j) = 0.0_rp
1294  toaflx_lw_dn(i,j) = 0.0_rp
1295  toaflx_sw_up(i,j) = 0.0_rp
1296  toaflx_sw_dn(i,j) = 0.0_rp
1297 
1298  sflx_rad_dn(i,j,1,1) = flx_sw_dn
1299  sflx_rad_dn(i,j,1,2) = 0.0_rp
1300  sflx_rad_dn(i,j,2,1) = flx_lw_dn
1301  sflx_rad_dn(i,j,2,2) = 0.0_rp
1302  enddo
1303  enddo
1304 
1305  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
module Atmosphere / Physics Cloud Microphysics
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
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
integer, public js
start point of inner domain: y, local
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
integer, public ie
end point of inner domain: x, local
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(:,:), allocatable, public atmos_phy_rd_toaflx_sw_dn
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
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 1311 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().

1311  use mod_land_vars, only: &
1312  land_temp, &
1313  land_water, &
1314  land_sfc_temp, &
1316  implicit none
1317  ! Land state
1318  real(RP) :: lnd_temp ! soil temperature [K]
1319  real(RP) :: lnd_water = 0.15_rp ! soil moisture [m3/m3]
1320  real(RP) :: sfc_temp ! land skin temperature [K]
1321  real(RP) :: sfc_albedo_lw = 0.01_rp ! land surface albedo for LW [0-1]
1322  real(RP) :: sfc_albedo_sw = 0.20_rp ! land surface albedo for SW [0-1]
1323 
1324  integer :: i, j
1325  integer :: ierr
1326 
1327  namelist /param_mkinit_land/ &
1328  lnd_temp, &
1329  lnd_water, &
1330  sfc_temp, &
1331  sfc_albedo_lw, &
1332  sfc_albedo_sw
1333 
1334  lnd_temp = thetastd
1335  sfc_temp = thetastd
1336 
1337  !--- read namelist
1338  rewind(io_fid_conf)
1339  read(io_fid_conf,nml=param_mkinit_land,iostat=ierr)
1340  if( ierr < 0 ) then !--- missing
1341  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1342  elseif( ierr > 0 ) then !--- fatal error
1343  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_LAND. Check!'
1344  call prc_mpistop
1345  endif
1346  if( io_lnml ) write(io_fid_log,nml=param_mkinit_land)
1347 
1348  do j = js, je
1349  do i = is, ie
1350  land_temp(:,i,j) = lnd_temp
1351  land_water(:,i,j) = lnd_water
1352  land_sfc_temp(i,j) = sfc_temp
1353  land_sfc_albedo(i,j,i_lw) = sfc_albedo_lw
1354  land_sfc_albedo(i,j,i_sw) = sfc_albedo_sw
1355  end do
1356  end do
1357 
1358  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
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]
integer, public js
start point of inner domain: y, local
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]
integer, public ie
end point of inner domain: x, local
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
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 1364 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().

1364  use mod_ocean_vars, only: &
1365  ocean_temp, &
1366  ocean_sfc_temp, &
1367  ocean_sfc_albedo, &
1368  ocean_sfc_z0m, &
1369  ocean_sfc_z0h, &
1371  implicit none
1372  ! Ocean state
1373  real(RP) :: ocn_temp ! ocean temperature [K]
1374  real(RP) :: sfc_temp ! ocean skin temperature [K]
1375  real(RP) :: sfc_albedo_lw = 0.04_rp ! ocean surface albedo for LW [0-1]
1376  real(RP) :: sfc_albedo_sw = 0.05_rp ! ocean surface albedo for SW [0-1]
1377  real(RP) :: sfc_z0m = 1.0e-4_rp ! ocean surface roughness length (momentum) [m]
1378  real(RP) :: sfc_z0h = 1.0e-4_rp ! ocean surface roughness length (heat) [m]
1379  real(RP) :: sfc_z0e = 1.0e-4_rp ! ocean surface roughness length (vapor) [m]
1380 
1381  integer :: i, j
1382  integer :: ierr
1383 
1384  namelist /param_mkinit_ocean/ &
1385  ocn_temp, &
1386  sfc_temp, &
1387  sfc_albedo_lw, &
1388  sfc_albedo_sw, &
1389  sfc_z0m, &
1390  sfc_z0h, &
1391  sfc_z0e
1392 
1393  ocn_temp = thetastd
1394  sfc_temp = thetastd
1395 
1396  !--- read namelist
1397  rewind(io_fid_conf)
1398  read(io_fid_conf,nml=param_mkinit_ocean,iostat=ierr)
1399  if( ierr < 0 ) then !--- missing
1400  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1401  elseif( ierr > 0 ) then !--- fatal error
1402  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_OCEAN. Check!'
1403  call prc_mpistop
1404  endif
1405  if( io_lnml ) write(io_fid_log,nml=param_mkinit_ocean)
1406 
1407 
1408  do j = js, je
1409  do i = is, ie
1410  ocean_temp(i,j) = ocn_temp
1411  ocean_sfc_temp(i,j) = sfc_temp
1412  ocean_sfc_albedo(i,j,i_lw) = sfc_albedo_lw
1413  ocean_sfc_albedo(i,j,i_sw) = sfc_albedo_sw
1414  ocean_sfc_z0m(i,j) = sfc_z0m
1415  ocean_sfc_z0h(i,j) = sfc_z0h
1416  ocean_sfc_z0e(i,j) = sfc_z0e
1417  end do
1418  end do
1419 
1420  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), dimension(:,:,:), allocatable, public ocean_sfc_albedo
ocean surface albedo [0-1]
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
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
integer, public js
start point of inner domain: y, local
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:,:), allocatable, public ocean_sfc_z0m
ocean surface roughness length for momentum [m]
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(:,:), allocatable, public ocean_sfc_temp
ocean surface skin temperature [K]
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
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 1426 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().

1426  use mod_urban_vars, only: &
1427  urban_tr, &
1428  urban_tb, &
1429  urban_tg, &
1430  urban_tc, &
1431  urban_qc, &
1432  urban_uc, &
1433  urban_trl, &
1434  urban_tbl, &
1435  urban_tgl, &
1436  urban_rainr, &
1437  urban_rainb, &
1438  urban_raing, &
1439  urban_roff, &
1440  urban_sfc_temp, &
1442  implicit none
1443  ! urban state
1444  real(RP) :: urb_roof_temp ! Surface temperature of roof [K]
1445  real(RP) :: urb_bldg_temp ! Surface temperature of building [K
1446  real(RP) :: urb_grnd_temp ! Surface temperature of ground [K]
1447  real(RP) :: urb_cnpy_temp ! Diagnostic canopy air temperature
1448  real(RP) :: urb_cnpy_hmdt = 0.0_rp ! Diagnostic canopy humidity [-]
1449  real(RP) :: urb_cnpy_wind = 0.0_rp ! Diagnostic canopy wind [m/s]
1450  real(RP) :: urb_roof_layer_temp ! temperature in layer of roof [K]
1451  real(RP) :: urb_bldg_layer_temp ! temperature in layer of building [
1452  real(RP) :: urb_grnd_layer_temp ! temperature in layer of ground [K]
1453  real(RP) :: urb_roof_rain = 0.0_rp ! temperature in layer of roof [K]
1454  real(RP) :: urb_bldg_rain = 0.0_rp ! temperature in layer of building [
1455  real(RP) :: urb_grnd_rain = 0.0_rp ! temperature in layer of ground [K]
1456  real(RP) :: urb_runoff = 0.0_rp ! temperature in layer of ground [K]
1457  real(RP) :: urb_sfc_temp ! Grid average of surface temperature [K]
1458  real(RP) :: urb_alb_lw = 0.0_rp ! Grid average of surface albedo for LW [0-1]
1459  real(RP) :: urb_alb_sw = 0.0_rp ! Grid average of surface albedo for SW [0-1]
1460 
1461  integer :: i, j
1462  integer :: ierr
1463 
1464  namelist /param_mkinit_urban/ &
1465  urb_roof_temp, &
1466  urb_bldg_temp, &
1467  urb_grnd_temp, &
1468  urb_cnpy_temp, &
1469  urb_cnpy_hmdt, &
1470  urb_cnpy_wind, &
1471  urb_roof_layer_temp, &
1472  urb_bldg_layer_temp, &
1473  urb_grnd_layer_temp, &
1474  urb_roof_rain, &
1475  urb_bldg_rain, &
1476  urb_grnd_rain, &
1477  urb_runoff, &
1478  urb_sfc_temp, &
1479  urb_alb_lw, &
1480  urb_alb_sw
1481 
1482  urb_roof_temp = thetastd
1483  urb_bldg_temp = thetastd
1484  urb_grnd_temp = thetastd
1485  urb_cnpy_temp = thetastd
1486  urb_roof_layer_temp = thetastd
1487  urb_bldg_layer_temp = thetastd
1488  urb_grnd_layer_temp = thetastd
1489 
1490  urb_sfc_temp = thetastd
1491 
1492  !--- read namelist
1493  rewind(io_fid_conf)
1494  read(io_fid_conf,nml=param_mkinit_urban,iostat=ierr)
1495  if( ierr < 0 ) then !--- missing
1496  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
1497  elseif( ierr > 0 ) then !--- fatal error
1498  write(*,*) 'xxx Not appropriate names in namelist PARAM_MKINIT_URBAN. Check!'
1499  call prc_mpistop
1500  endif
1501  if( io_lnml ) write(io_fid_log,nml=param_mkinit_urban)
1502 
1503 
1504  do j = js, je
1505  do i = is, ie
1506  urban_tr(i,j) = urb_roof_temp
1507  urban_tb(i,j) = urb_bldg_temp
1508  urban_tg(i,j) = urb_grnd_temp
1509  urban_tc(i,j) = urb_cnpy_temp
1510  urban_qc(i,j) = urb_cnpy_hmdt
1511  urban_uc(i,j) = urb_cnpy_wind
1512  urban_trl(:,i,j) = urb_roof_layer_temp
1513  urban_tbl(:,i,j) = urb_bldg_layer_temp
1514  urban_tgl(:,i,j) = urb_grnd_layer_temp
1515  urban_rainr(i,j) = urb_roof_rain
1516  urban_rainb(i,j) = urb_bldg_rain
1517  urban_raing(i,j) = urb_grnd_rain
1518  urban_roff(i,j) = urb_runoff
1519  urban_sfc_temp(i,j) = urb_sfc_temp
1520  urban_sfc_albedo(i,j,i_lw) = urb_alb_lw
1521  urban_sfc_albedo(i,j,i_sw) = urb_alb_sw
1522  end do
1523  end do
1524 
1525  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:), allocatable, public urban_qc
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
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
integer, public js
start point of inner domain: y, local
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
integer, public ie
end point of inner domain: x, local
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(:,:), allocatable, public urban_tg
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
real(rp), dimension(:,:,:), allocatable, public urban_trl
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
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 1532 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().

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

4573  implicit none
4574 
4575  real(RP), intent(out) :: fact0
4576  real(RP), intent(out) :: fact1
4577  integer, intent(out) :: idx0
4578  integer, intent(out) :: idx1
4579  real(RP), intent(in) :: x
4580  integer, intent(in) :: nx
4581  real(RP), intent(in) :: x_org(nx)
4582  logical, intent(in) :: loop
4583 
4584  real(RP) :: xwork
4585  integer :: i
4586 
4587  if ( x < x_org(1) ) then
4588  if ( loop ) then
4589  xwork = x_org(1) - ( x_org(2) - x_org(1) )**2 / ( x_org(3) - x_org(2) )
4590  fact0 = ( x_org(1) - x ) / ( x_org(1) - xwork )
4591  fact1 = ( x - xwork ) / ( x_org(1) - xwork )
4592  idx0 = nx
4593  idx1 = 1
4594  else
4595  fact0 = ( x_org(2) - x ) / ( x_org(2) - x_org(1) )
4596  fact1 = ( x - x_org(1) ) / ( x_org(2) - x_org(1) )
4597  idx0 = 1
4598  idx1 = 2
4599  end if
4600  else if ( x > x_org(nx) ) then
4601  if ( loop ) then
4602  xwork = x_org(nx) + ( x_org(nx) - x_org(nx-1) )**2 / ( x_org(nx-1) - x_org(nx-2) )
4603  fact0 = ( xwork - x ) / ( xwork - x_org(nx) )
4604  fact1 = ( x - x_org(nx) ) / ( xwork - x_org(nx) )
4605  idx0 = nx
4606  idx1 = 1
4607  else
4608  fact0 = ( x_org(nx) - x ) / ( x_org(nx) - x_org(nx-1) )
4609  fact1 = ( x - x_org(nx-1) ) / ( x_org(nx) - x_org(nx-1) )
4610  idx0 = nx-1
4611  idx1 = nx
4612  end if
4613  else
4614  do i = 2, nx
4615  if ( x <= x_org(i) ) then
4616  fact0 = ( x_org(i) - x ) / ( x_org(i) - x_org(i-1) )
4617  fact1 = ( x - x_org(i-1) ) / ( x_org(i) - x_org(i-1) )
4618  idx0 = i-1
4619  idx1 = i
4620  exit
4621  end if
4622  end do
4623  end if
4624 
4625  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