SCALE-RM
Functions/Subroutines
scale_ocean_sfc_const Module Reference

module OCEAN / Surface flux with constant ocean model More...

Functions/Subroutines

subroutine, public ocean_sfc_const_setup (OCEAN_TYPE)
 Setup. More...
 
subroutine, public ocean_sfc_const (SST_t, ZMFLX, XMFLX, YMFLX, SHFLX, LHFLX, WHFLX, U10, V10, T2, Q2, TMPA, PRSA, WA, UA, VA, RHOA, QVA, Z1, PBL, PRSS, LWD, SWD, TW, SST, ALB_LW, ALB_SW, Z0M, Z0H, Z0E, dt)
 

Detailed Description

module OCEAN / Surface flux with constant ocean model

Description
Surface flux with constant ocean model
Author
Team SCALE

Function/Subroutine Documentation

◆ ocean_sfc_const_setup()

subroutine, public scale_ocean_sfc_const::ocean_sfc_const_setup ( character(len=*), intent(in)  OCEAN_TYPE)

Setup.

Definition at line 46 of file scale_ocean_sfc_const.F90.

References scale_stdio::io_fid_log, and scale_stdio::io_l.

Referenced by scale_ocean_sfc::ocean_sfc_setup().

46  implicit none
47 
48  character(len=*), intent(in) :: OCEAN_TYPE
49  !---------------------------------------------------------------------------
50 
51  if( io_l ) write(io_fid_log,*)
52  if( io_l ) write(io_fid_log,*) '++++++ Module[CONST] / Categ[OCEAN SFC] / Origin[SCALElib]'
53 
54  return
Here is the caller graph for this function:

◆ ocean_sfc_const()

subroutine, public scale_ocean_sfc_const::ocean_sfc_const ( real(rp), dimension(ia,ja), intent(out)  SST_t,
real(rp), dimension(ia,ja), intent(out)  ZMFLX,
real(rp), dimension(ia,ja), intent(out)  XMFLX,
real(rp), dimension(ia,ja), intent(out)  YMFLX,
real(rp), dimension(ia,ja), intent(out)  SHFLX,
real(rp), dimension(ia,ja), intent(out)  LHFLX,
real(rp), dimension(ia,ja), intent(out)  WHFLX,
real(rp), dimension (ia,ja), intent(out)  U10,
real(rp), dimension (ia,ja), intent(out)  V10,
real(rp), dimension (ia,ja), intent(out)  T2,
real(rp), dimension (ia,ja), intent(out)  Q2,
real(rp), dimension(ia,ja), intent(in)  TMPA,
real(rp), dimension(ia,ja), intent(in)  PRSA,
real(rp), dimension (ia,ja), intent(in)  WA,
real(rp), dimension (ia,ja), intent(in)  UA,
real(rp), dimension (ia,ja), intent(in)  VA,
real(rp), dimension(ia,ja), intent(in)  RHOA,
real(rp), dimension (ia,ja), intent(in)  QVA,
real(rp), dimension (ia,ja), intent(in)  Z1,
real(rp), dimension (ia,ja), intent(in)  PBL,
real(rp), dimension(ia,ja), intent(in)  PRSS,
real(rp), dimension (ia,ja), intent(in)  LWD,
real(rp), dimension (ia,ja), intent(in)  SWD,
real(rp), dimension (ia,ja), intent(in)  TW,
real(rp), dimension (ia,ja), intent(in)  SST,
real(rp), dimension(ia,ja), intent(in)  ALB_LW,
real(rp), dimension(ia,ja), intent(in)  ALB_SW,
real(rp), dimension (ia,ja), intent(in)  Z0M,
real(rp), dimension (ia,ja), intent(in)  Z0H,
real(rp), dimension (ia,ja), intent(in)  Z0E,
real(dp), intent(in)  dt 
)

Definition at line 90 of file scale_ocean_sfc_const.F90.

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_stb, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_landuse::landuse_fact_ocean, and scale_process::prc_mpistop().

Referenced by scale_ocean_sfc::ocean_sfc_setup().

90  use scale_process, only: &
92  use scale_const, only: &
93  pre00 => const_pre00, &
94  rdry => const_rdry, &
95  cpdry => const_cpdry, &
96  stb => const_stb
97  use scale_landuse, only: &
98  landuse_fact_ocean
99  use scale_bulkflux, only: &
100  bulkflux
101  use scale_atmos_hydrometeor, only: &
102  hydrometeor_lhv => atmos_hydrometeor_lhv
103  use scale_atmos_saturation, only: &
104  qsat => atmos_saturation_pres2qsat_all
105  implicit none
106 
107  ! arguments
108  real(RP), intent(out) :: SST_t(IA,JA) ! tendency of sea surface temperature
109  real(RP), intent(out) :: ZMFLX(IA,JA) ! z-momentum flux at the surface [kg/m/s2]
110  real(RP), intent(out) :: XMFLX(IA,JA) ! x-momentum flux at the surface [kg/m/s2]
111  real(RP), intent(out) :: YMFLX(IA,JA) ! y-momentum flux at the surface [kg/m/s2]
112  real(RP), intent(out) :: SHFLX(IA,JA) ! sensible heat flux at the surface [W/m2]
113  real(RP), intent(out) :: LHFLX(IA,JA) ! latent heat flux at the surface [W/m2]
114  real(RP), intent(out) :: WHFLX(IA,JA) ! water heat flux at the surface [W/m2]
115  real(RP), intent(out) :: U10 (IA,JA) ! velocity u at 10m [m/s]
116  real(RP), intent(out) :: V10 (IA,JA) ! velocity v at 10m [m/s]
117  real(RP), intent(out) :: T2 (IA,JA) ! temperature at 2m [K]
118  real(RP), intent(out) :: Q2 (IA,JA) ! water vapor at 2m [kg/kg]
119 
120  real(RP), intent(in) :: TMPA(IA,JA) ! temperature at the lowest atmospheric layer [K]
121  real(RP), intent(in) :: PRSA(IA,JA) ! pressure at the lowest atmospheric layer [Pa]
122  real(RP), intent(in) :: WA (IA,JA) ! velocity w at the lowest atmospheric layer [m/s]
123  real(RP), intent(in) :: UA (IA,JA) ! velocity u at the lowest atmospheric layer [m/s]
124  real(RP), intent(in) :: VA (IA,JA) ! velocity v at the lowest atmospheric layer [m/s]
125  real(RP), intent(in) :: RHOA(IA,JA) ! density at the lowest atmospheric layer [kg/m3]
126  real(RP), intent(in) :: QVA (IA,JA) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
127  real(RP), intent(in) :: Z1 (IA,JA) ! cell center height at the lowest atmospheric layer [m]
128  real(RP), intent(in) :: PBL (IA,JA) ! the top of atmospheric mixing layer [m]
129  real(RP), intent(in) :: PRSS(IA,JA) ! pressure at the surface [Pa]
130  real(RP), intent(in) :: LWD (IA,JA) ! downward long-wave radiation flux at the surface [W/m2]
131  real(RP), intent(in) :: SWD (IA,JA) ! downward short-wave radiation flux at the surface [W/m2]
132 
133  real(RP), intent(in) :: TW (IA,JA) ! water temperature [K]
134  real(RP), intent(in) :: SST (IA,JA) ! sea surface temperature [K]
135  real(RP), intent(in) :: ALB_LW(IA,JA) ! surface albedo for LW (0-1)
136  real(RP), intent(in) :: ALB_SW(IA,JA) ! surface albedo for SW (0-1)
137  real(RP), intent(in) :: Z0M (IA,JA) ! roughness length for momentum [m]
138  real(RP), intent(in) :: Z0H (IA,JA) ! roughness length for heat [m]
139  real(RP), intent(in) :: Z0E (IA,JA) ! roughness length for vapor [m]
140  real(DP), intent(in) :: dt ! delta time
141 
142  ! works
143  real(RP) :: Ustar ! friction velocity [m]
144  real(RP) :: Tstar ! friction temperature [K]
145  real(RP) :: Qstar ! friction mixing rate [kg/kg]
146  real(RP) :: Uabs ! modified absolute velocity [m/s]
147 
148  real(RP) :: QVsat ! saturation water vapor mixing ratio at surface [kg/kg]
149  real(RP) :: LHV(IA,JA) ! latent heat of vaporization [J/kg]
150 
151  real(RP) :: FracU10 ! calculation parameter for U10 [-]
152  real(RP) :: FracT2 ! calculation parameter for T2 [-]
153  real(RP) :: FracQ2 ! calculation parameter for Q2 [-]
154 
155  integer :: i, j
156  !---------------------------------------------------------------------------
157 
158  if( io_l ) write(io_fid_log,*) '*** Ocean surface step: Const'
159 
160  call hydrometeor_lhv( lhv(:,:), tmpa(:,:) )
161 
162  ! calculate tendency
163  do j = js, je
164  do i = is, ie
165  sst_t(i,j) = 0.0_rp
166  enddo
167  enddo
168 
169  ! calculate surface flux
170  do j = js, je
171  do i = is, ie
172 
173  if( landuse_fact_ocean(i,j) > 0.0_rp ) then
174 
175  ! saturation at the surface
176  call qsat( qvsat, & ! [OUT]
177  sst(i,j), & ! [IN]
178  prss(i,j) ) ! [IN]
179 
180  call bulkflux( &
181  ustar, & ! [OUT]
182  tstar, & ! [OUT]
183  qstar, & ! [OUT]
184  uabs, & ! [OUT]
185  fracu10, & ! [OUT]
186  fract2, & ! [OUT]
187  fracq2, & ! [OUT]
188  tmpa(i,j), & ! [IN]
189  sst(i,j), & ! [IN]
190  prsa(i,j), & ! [IN]
191  prss(i,j), & ! [IN]
192  qva(i,j), & ! [IN]
193  qvsat, & ! [IN]
194  ua(i,j), & ! [IN]
195  va(i,j), & ! [IN]
196  z1(i,j), & ! [IN]
197  pbl(i,j), & ! [IN]
198  z0m(i,j), & ! [IN]
199  z0h(i,j), & ! [IN]
200  z0e(i,j) ) ! [IN]
201 
202  zmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * wa(i,j)
203  xmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * ua(i,j)
204  ymflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * va(i,j)
205 
206  shflx(i,j) = - rhoa(i,j) * ustar * tstar &
207  * cpdry * ( prss(i,j) / pre00 )**( rdry/cpdry )
208  lhflx(i,j) = - rhoa(i,j) * ustar * qstar * lhv(i,j)
209 
210  ! calculation for residual
211  whflx(i,j) = ( 1.0_rp - alb_sw(i,j) ) * swd(i,j) * (-1.0_rp) &
212  - ( 1.0_rp - alb_lw(i,j) ) * ( lwd(i,j) - stb * sst(i,j)**4 ) &
213  + shflx(i,j) + lhflx(i,j)
214 
215  ! diagnostic variables considering unstable/stable state
216  !U10(i,j) = FracU10 * UA(i,j)
217  !V10(i,j) = FracU10 * VA(i,j)
218  !T2 (i,j) = ( 1.0_RP - FracT2 ) * SST(i,j) + FracT2 * TMPA(i,j)
219  !Q2 (i,j) = ( 1.0_RP - FracQ2 ) * QVsat + FracQ2 * QVA (i,j)
220 
221  ! diagnostic variables for neutral state
222  u10(i,j) = ua(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
223  v10(i,j) = va(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
224  t2(i,j) = sst(i,j) + ( tmpa(i,j) - sst(i,j) ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0h(i,j) ) ) &
225  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0h(i,j) ) )
226  q2(i,j) = qvsat + ( qva(i,j) - qvsat ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0e(i,j) ) ) &
227  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0e(i,j) ) )
228 
229  else
230 
231  ! not calculate surface flux
232  zmflx(i,j) = 0.0_rp
233  xmflx(i,j) = 0.0_rp
234  ymflx(i,j) = 0.0_rp
235  shflx(i,j) = 0.0_rp
236  lhflx(i,j) = 0.0_rp
237  whflx(i,j) = 0.0_rp
238  u10(i,j) = 0.0_rp
239  v10(i,j) = 0.0_rp
240  t2(i,j) = 0.0_rp
241  q2(i,j) = 0.0_rp
242 
243  end if
244 
245  enddo
246  enddo
247 
248  return
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
module LANDUSE
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module Surface bulk flux
Here is the call graph for this function:
Here is the caller graph for this function: