FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines
scale_atm_phy_sf_bulk_simple Module Reference

module atmosphere / physics / surface / simple More...

Functions/Subroutines

subroutine, public atmos_phy_sf_simple_setup
 Setup.
 
subroutine, public atmos_phy_sf_simple_flux (ia, is, ie, ja, js, je, atm_w, atm_u, atm_v, atm_temp, atm_pres, atm_qv, sfc_dens, sfc_temp, sfc_pres, atm_z1, sflx_mw, sflx_mu, sflx_mv, sflx_sh, sflx_lh, sflx_qv, u10, v10)
 Constant flux.
 

Detailed Description

module atmosphere / physics / surface / simple

Description
Flux from/to bottom wall of atmosphere (surface) Constant bulk coefficient
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_SF_SIMPLE
    nametypedefault valuecomment
    ATMOS_PHY_SF_U_MINM real(RP) 0.0_RP minimum limit of absolute velocity for momentum [m/s]
    ATMOS_PHY_SF_CONST_CM real(RP) 0.0011_RP constant bulk coefficient for momentum [NIL]
    ATMOS_PHY_SF_CONST_CH real(RP) 0.0044_RP constant bulk coefficient for heat [NIL]
    ATMOS_PHY_SF_CONST_CE real(RP) 0.0044_RP constant bulk coefficient for evaporation [NIL]
    ATMOS_PHY_SF_BULK_BETA real(RP) 1.0_RP evaporation efficiency (0-1)

History Output
No history output

Function/Subroutine Documentation

◆ atmos_phy_sf_simple_setup()

subroutine, public scale_atm_phy_sf_bulk_simple::atmos_phy_sf_simple_setup

Setup.

Definition at line 64 of file scale_atm_phy_sf_bulk_simple.F90.

65 use scale_prc, only: &
66 prc_abort
67 implicit none
68
69 namelist / param_atmos_phy_sf_simple / &
70 atmos_phy_sf_u_minm, &
71 atmos_phy_sf_const_cm, &
72 atmos_phy_sf_const_ch, &
73 atmos_phy_sf_const_ce, &
74 atmos_phy_sf_bulk_beta
75
76 integer :: ierr
77 !---------------------------------------------------------------------------
78
79 log_newline
80 log_info("ATMOS_PHY_SF_simple_setup",*) 'Setup'
81 log_info("ATMOS_PHY_SF_simple_setup",*) 'Simple flux'
82
83 !--- read namelist
84 rewind(io_fid_conf)
85 read(io_fid_conf,nml=param_atmos_phy_sf_simple,iostat=ierr)
86 if( ierr < 0 ) then !--- missing
87 log_info("ATMOS_PHY_SF_simple_setup",*) 'Not found namelist. Default used.'
88 elseif( ierr > 0 ) then !--- fatal error
89 log_error("ATMOS_PHY_SF_simple_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_SF_SIMPLE. Check!'
90 call prc_abort
91 endif
92 log_nml(param_atmos_phy_sf_simple)
93
94 return

◆ atmos_phy_sf_simple_flux()

subroutine, public scale_atm_phy_sf_bulk_simple::atmos_phy_sf_simple_flux ( integer, intent(in) ia,
integer, intent(in) is,
integer, intent(in) ie,
integer, intent(in) ja,
integer, intent(in) js,
integer, intent(in) je,
real(rp), dimension (ia,ja), intent(in) atm_w,
real(rp), dimension (ia,ja), intent(in) atm_u,
real(rp), dimension (ia,ja), intent(in) atm_v,
real(rp), dimension(ia,ja), intent(in) atm_temp,
real(rp), dimension(ia,ja), intent(in) atm_pres,
real(rp), dimension (ia,ja), intent(in) atm_qv,
real(rp), dimension(ia,ja), intent(in) sfc_dens,
real(rp), dimension(ia,ja), intent(in) sfc_temp,
real(rp), dimension(ia,ja), intent(in) sfc_pres,
real(rp), dimension (ia,ja), intent(in) atm_z1,
real(rp), dimension(ia,ja), intent(out) sflx_mw,
real(rp), dimension(ia,ja), intent(out) sflx_mu,
real(rp), dimension(ia,ja), intent(out) sflx_mv,
real(rp), dimension(ia,ja), intent(out) sflx_sh,
real(rp), dimension(ia,ja), intent(out) sflx_lh,
real(rp), dimension(ia,ja), intent(out) sflx_qv,
real(rp), dimension (ia,ja), intent(out) u10,
real(rp), dimension (ia,ja), intent(out) v10 )

Constant flux.

Definition at line 99 of file scale_atm_phy_sf_bulk_simple.F90.

107 use scale_const, only: &
108 pi => const_pi
109 use scale_atmos_hydrometeor, only: &
110 hydrometeor_lhv => atmos_hydrometeor_lhv
111 use scale_atmos_saturation, only: &
112 saturation_psat_all => atmos_saturation_psat_all
113 use scale_time, only: &
114 time_nowsec
115 implicit none
116 integer, intent(in) :: IA, IS, IE
117 integer, intent(in) :: JA, JS, JE
118
119 real(RP), intent(in) :: ATM_W (IA,JA) ! velocity w at the lowermost layer (cell center) [m/s]
120 real(RP), intent(in) :: ATM_U (IA,JA) ! velocity u at the lowermost layer (cell center) [m/s]
121 real(RP), intent(in) :: ATM_V (IA,JA) ! velocity v at the lowermost layer (cell center) [m/s]
122 real(RP), intent(in) :: ATM_TEMP(IA,JA) ! temperature at the lowermost layer (cell center) [K]
123 real(RP), intent(in) :: ATM_PRES(IA,JA) ! pressure at the lowermost layer (cell center) [Pa]
124 real(RP), intent(in) :: ATM_QV (IA,JA) ! qv at the lowermost layer (cell center) [kg/kg]
125 real(RP), intent(in) :: SFC_DENS(IA,JA) ! density at the surface atmosphere [kg/m3]
126 real(RP), intent(in) :: SFC_TEMP(IA,JA) ! tempertire at the surface atmosphere [K]
127 real(RP), intent(in) :: SFC_PRES(IA,JA) ! pressure at the surface atmosphere [Pa]
128 real(RP), intent(in) :: ATM_Z1 (IA,JA) ! height of the lowermost grid from surface (cell center) [m]
129
130 real(RP), intent(out) :: SFLX_MW(IA,JA) ! surface flux for z-momentum (area center) [m/s*kg/m2/s]
131 real(RP), intent(out) :: SFLX_MU(IA,JA) ! surface flux for x-momentum (area center) [m/s*kg/m2/s]
132 real(RP), intent(out) :: SFLX_MV(IA,JA) ! surface flux for y-momentum (area center) [m/s*kg/m2/s]
133 real(RP), intent(out) :: SFLX_SH(IA,JA) ! surface flux for sensible heat (area center) [J/m2/s]
134 real(RP), intent(out) :: SFLX_LH(IA,JA) ! surface flux for latent heat (area center) [J/m2/s]
135 real(RP), intent(out) :: SFLX_QV(IA,JA) ! surface flux for qv (area center) [kg/m2/s]
136 real(RP), intent(out) :: U10 (IA,JA) ! velocity u at 10m height
137 real(RP), intent(out) :: V10 (IA,JA) ! velocity v at 10m height
138
139 real(RP) :: ATM_Uabs(IA,JA) ! absolute velocity at z1 [m/s]
140 real(RP) :: R10
141
142 real(RP) :: SFC_PSAT (IA,JA) ! saturatad water vapor pressure [Pa]
143 real(RP) :: LHV(IA,JA)
144
145 real(RP) :: SFC_QSAT ! saturatad water vapor mixing ratio [kg/kg]
146 real(RP) :: SFC_QV(IA,JA) ! water vapor mixing ratio [kg/kg]
147
148 integer :: i, j
149 !---------------------------------------------------------------------------
150
151 log_progress(*) 'atmosphere / physics / surface flux / simple'
152
153 !$omp parallel do
154 do j = js, je
155 do i = is, ie
156 atm_uabs(i,j) = min( atmos_phy_sf_u_maxm, max( atmos_phy_sf_u_minm, &
157 sqrt( atm_w(i,j)**2 + atm_u(i,j)**2 + atm_v(i,j)**2 ) ) ) ! at cell center
158 enddo
159 enddo
160
161 !-----< momentum >-----
162
163 !$omp parallel do
164 do j = js, je
165 do i = is, ie
166 sflx_mw(i,j) = - atmos_phy_sf_const_cm * atm_uabs(i,j) * sfc_dens(i,j) * atm_w(i,j)
167 sflx_mu(i,j) = - atmos_phy_sf_const_cm * atm_uabs(i,j) * sfc_dens(i,j) * atm_u(i,j)
168 sflx_mv(i,j) = - atmos_phy_sf_const_cm * atm_uabs(i,j) * sfc_dens(i,j) * atm_v(i,j)
169 enddo
170 enddo
171
172 !-----< heat & mass flux >-----
173 call saturation_psat_all( ia, is, ie, ja, js, je, &
174 sfc_temp(:,:), & ! [IN]
175 sfc_psat(:,:) ) ! [OUT]
176
177 call hydrometeor_lhv( &
178 ia, is, ie, ja, js, je, &
179 sfc_temp(:,:), & ! [IN]
180 lhv(:,:) ) ! [OUT]
181
182 !$omp parallel do private( SFC_QSAT )
183 do j = js, je
184 do i = is, ie
185 sfc_qsat = epsvap * sfc_psat(i,j) / ( sfc_pres(i,j) - ( 1.0_rp-epsvap ) * sfc_psat(i,j) )
186 sfc_qv(i,j) = ( 1.0_rp - atmos_phy_sf_bulk_beta ) * atm_qv(i,j) + atmos_phy_sf_bulk_beta * sfc_qsat
187
188 sflx_sh(i,j) = atmos_phy_sf_const_ch * atm_uabs(i,j) * sfc_dens(i,j) * cpdry * ( sfc_temp(i,j) - atm_temp(i,j) )
189 sflx_lh(i,j) = atmos_phy_sf_const_ce * atm_uabs(i,j) * sfc_dens(i,j) * lhv(i,j) * ( sfc_qv(i,j) - atm_qv(i,j) )
190
191 sflx_qv(i,j) = sflx_lh(i,j) / lhv(i,j)
192 enddo
193 enddo
194
195 !-----< U10, V10 >-----
196
197 !$omp parallel do &
198 !$omp private(R10)
199 do j = js, je
200 do i = is, ie
201 r10 = 10.0_rp / atm_z1(i,j)
202
203 u10(i,j) = r10 * atm_u(i,j)
204 v10(i,j) = r10 * atm_v(i,j)
205 enddo
206 enddo
207
208 return