SCALE-RM
Functions/Subroutines
scale_ocean_sfc_slab Module Reference

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

Functions/Subroutines

subroutine, public ocean_sfc_slab_setup (OCEAN_TYPE)
 Setup. More...
 
subroutine, public ocean_sfc_slab (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)
 
subroutine, public ocean_sfc_slab_simplealbedo (SFC_albedo_t, SFC_albedo, cosSZA, dt)
 

Detailed Description

module OCEAN / Surface flux with slab ocean model

Description
Surface flux with slab ocean model
Author
Team SCALE

Function/Subroutine Documentation

◆ ocean_sfc_slab_setup()

subroutine, public scale_ocean_sfc_slab::ocean_sfc_slab_setup ( character(len=*), intent(in)  OCEAN_TYPE)

Setup.

Definition at line 51 of file scale_ocean_sfc_slab.F90.

References scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::ja, 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().

51  use scale_process, only: &
53  use scale_landuse, only: &
55  implicit none
56 
57  character(len=*), intent(in) :: ocean_type
58 
59  integer :: i, j
60  !---------------------------------------------------------------------------
61 
62  if( io_l ) write(io_fid_log,*)
63  if( io_l ) write(io_fid_log,*) '++++++ Module[SLAB] / Categ[OCEAN SFC] / Origin[SCALElib]'
64 
65  select case( ocean_type )
66  case ( 'CONST' )
67  sst_update = .false.
68  case ( 'SLAB', 'FILE' )
69  sst_update = .true.
70  case default
71  write(*,*) 'xxx wrong OCEAN_TYPE. Check!'
72  call prc_mpistop
73  end select
74 
75  ! judge to run slab ocean model
76  allocate( is_ocn(ia,ja) )
77 
78  do j = js, je
79  do i = is, ie
80  if( landuse_fact_ocean(i,j) > 0.0_rp ) then
81  is_ocn(i,j) = .true.
82  else
83  is_ocn(i,j) = .false.
84  end if
85  end do
86  end do
87 
88  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.
integer, public ia
of x whole cells (local, with HALO)
module LANDUSE
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
integer, public js
start point of inner domain: y, local
module PROCESS
integer, public ie
end point of inner domain: x, local
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:

◆ ocean_sfc_slab()

subroutine, public scale_ocean_sfc_slab::ocean_sfc_slab ( 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 124 of file scale_ocean_sfc_slab.F90.

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_stb, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, dc_log::log(), and scale_process::prc_mpistop().

Referenced by scale_ocean_sfc::ocean_sfc_setup().

124  use scale_grid_index
125  use scale_process, only: &
127  use scale_const, only: &
128  cpdry => const_cpdry, &
129  stb => const_stb
130  use scale_bulkflux, only: &
131  bulkflux
132  use scale_atmos_thermodyn, only: &
133  atmos_thermodyn_templhv
134  use scale_atmos_saturation, only: &
135  qsat => atmos_saturation_pres2qsat_all
136  implicit none
137 
138  ! arguments
139  real(RP), intent(out) :: sst_t(ia,ja) ! tendency of sea surface temperature
140  real(RP), intent(out) :: zmflx(ia,ja) ! z-momentum flux at the surface [kg/m2/s]
141  real(RP), intent(out) :: xmflx(ia,ja) ! x-momentum flux at the surface [kg/m2/s]
142  real(RP), intent(out) :: ymflx(ia,ja) ! y-momentum flux at the surface [kg/m2/s]
143  real(RP), intent(out) :: shflx(ia,ja) ! sensible heat flux at the surface [W/m2]
144  real(RP), intent(out) :: lhflx(ia,ja) ! latent heat flux at the surface [W/m2]
145  real(RP), intent(out) :: whflx(ia,ja) ! water heat flux at the surface [W/m2]
146  real(RP), intent(out) :: u10 (ia,ja) ! velocity u at 10m [m/s]
147  real(RP), intent(out) :: v10 (ia,ja) ! velocity v at 10m [m/s]
148  real(RP), intent(out) :: t2 (ia,ja) ! temperature at 2m [K]
149  real(RP), intent(out) :: q2 (ia,ja) ! water vapor at 2m [kg/kg]
150 
151  real(RP), intent(in) :: tmpa(ia,ja) ! temperature at the lowest atmospheric layer [K]
152  real(RP), intent(in) :: prsa(ia,ja) ! pressure at the lowest atmospheric layer [Pa]
153  real(RP), intent(in) :: wa (ia,ja) ! velocity w at the lowest atmospheric layer [m/s]
154  real(RP), intent(in) :: ua (ia,ja) ! velocity u at the lowest atmospheric layer [m/s]
155  real(RP), intent(in) :: va (ia,ja) ! velocity v at the lowest atmospheric layer [m/s]
156  real(RP), intent(in) :: rhoa(ia,ja) ! density at the lowest atmospheric layer [kg/m3]
157  real(RP), intent(in) :: qva (ia,ja) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
158  real(RP), intent(in) :: z1 (ia,ja) ! cell center height at the lowest atmospheric layer [m]
159  real(RP), intent(in) :: pbl (ia,ja) ! the top of atmospheric mixing layer [m]
160  real(RP), intent(in) :: prss(ia,ja) ! pressure at the surface [Pa]
161  real(RP), intent(in) :: lwd (ia,ja) ! downward long-wave radiation flux at the surface [W/m2]
162  real(RP), intent(in) :: swd (ia,ja) ! downward short-wave radiation flux at the surface [W/m2]
163 
164  real(RP), intent(in) :: tw (ia,ja) ! water temperature [K]
165  real(RP), intent(in) :: sst (ia,ja) ! sea surface temperature [K]
166  real(RP), intent(in) :: alb_lw(ia,ja) ! surface albedo for LW [0-1]
167  real(RP), intent(in) :: alb_sw(ia,ja) ! surface albedo for SW [0-1]
168  real(RP), intent(in) :: z0m (ia,ja) ! roughness length for momentum [m]
169  real(RP), intent(in) :: z0h (ia,ja) ! roughness length for heat [m]
170  real(RP), intent(in) :: z0e (ia,ja) ! roughness length for vapor [m]
171  real(DP), intent(in) :: dt ! delta time
172 
173  ! works
174  real(RP) :: sst1(ia,ja)
175 
176  real(RP) :: ustar ! friction velocity [m]
177  real(RP) :: tstar ! friction temperature [K]
178  real(RP) :: qstar ! friction mixing rate [kg/kg]
179  real(RP) :: uabs ! modified absolute velocity [m/s]
180 
181  real(RP) :: sqv ! saturation water vapor mixing ratio at surface [kg/kg]
182  real(RP) :: lhv(ia,ja) ! latent heat for vaporization depending on temperature [J/kg]
183 
184  integer :: i, j
185  !---------------------------------------------------------------------------
186 
187  ! update surface temperature
188  do j = js, je
189  do i = is, ie
190  sst1(i,j) = tw(i,j) ! assumed well-mixed condition
191  end do
192  end do
193 
194  call atmos_thermodyn_templhv( lhv, tmpa )
195 
196  do j = js, je
197  do i = is, ie
198 
199  if( is_ocn(i,j) ) then
200 
201  ! saturation at the surface
202  call qsat( sqv, & ! [OUT]
203  sst1(i,j), & ! [IN]
204  prss(i,j) ) ! [IN]
205 
206  call bulkflux( &
207  ustar, & ! [OUT]
208  tstar, & ! [OUT]
209  qstar, & ! [OUT]
210  uabs, & ! [OUT]
211  tmpa(i,j), & ! [IN]
212  sst1(i,j), & ! [IN]
213  prsa(i,j), & ! [IN]
214  prss(i,j), & ! [IN]
215  qva(i,j), & ! [IN]
216  sqv, & ! [IN]
217  ua(i,j), & ! [IN]
218  va(i,j), & ! [IN]
219  z1(i,j), & ! [IN]
220  pbl(i,j), & ! [IN]
221  z0m(i,j), & ! [IN]
222  z0h(i,j), & ! [IN]
223  z0e(i,j) ) ! [IN]
224 
225  zmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * wa(i,j)
226  xmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * ua(i,j)
227  ymflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * va(i,j)
228 
229  shflx(i,j) = -cpdry * rhoa(i,j) * ustar * tstar
230  lhflx(i,j) = -lhv(i,j) * rhoa(i,j) * ustar * qstar
231 
232  ! calculation for residual
233  whflx(i,j) = ( 1.0_rp - alb_sw(i,j) ) * swd(i,j) * (-1.0_rp) &
234  - ( 1.0_rp - alb_lw(i,j) ) * ( lwd(i,j) - stb * sst1(i,j)**4 ) &
235  + shflx(i,j) + lhflx(i,j)
236 
237  ! diagnositc variables
238  u10(i,j) = ua(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
239  v10(i,j) = va(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
240  t2(i,j) = sst1(i,j) + ( tmpa(i,j) - sst1(i,j) ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0h(i,j) ) ) &
241  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0h(i,j) ) )
242  q2(i,j) = sqv + ( qva(i,j) - sqv ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0e(i,j) ) ) &
243  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0e(i,j) ) )
244 
245  else
246  ! not calculate surface flux
247  zmflx(i,j) = 0.0_rp
248  xmflx(i,j) = 0.0_rp
249  ymflx(i,j) = 0.0_rp
250  shflx(i,j) = 0.0_rp
251  lhflx(i,j) = 0.0_rp
252  whflx(i,j) = 0.0_rp
253  u10(i,j) = 0.0_rp
254  v10(i,j) = 0.0_rp
255  t2(i,j) = 0.0_rp
256  q2(i,j) = 0.0_rp
257  end if
258 
259  enddo
260  enddo
261 
262  ! calculate tendency
263  if( sst_update ) then
264  do j = js, je
265  do i = is, ie
266  sst_t(i,j) = ( sst1(i,j) - sst(i,j) ) / dt
267  enddo
268  enddo
269  else
270  do j = js, je
271  do i = is, ie
272  sst_t(i,j) = 0.0_rp
273  enddo
274  enddo
275  end if
276 
277  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
integer, public va
Definition: scale_index.F90:38
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
Definition: scale_const.F90:51
module grid index
integer, public ia
of x whole cells (local, with HALO)
procedure(bc), pointer, public bulkflux
integer, public js
start point of inner domain: y, local
module PROCESS
subroutine, public log(type, message)
Definition: dc_log.f90:133
module CONSTANT
Definition: scale_const.F90:14
integer, public ie
end point of inner domain: x, local
module Surface bulk flux
module ATMOSPHERE / Thermodynamics
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:

◆ ocean_sfc_slab_simplealbedo()

subroutine, public scale_ocean_sfc_slab::ocean_sfc_slab_simplealbedo ( real(rp), dimension(ia,ja,2), intent(out)  SFC_albedo_t,
real(rp), dimension (ia,ja,2), intent(in)  SFC_albedo,
real(rp), dimension (ia,ja), intent(in)  cosSZA,
real(dp), intent(in)  dt 
)

Definition at line 286 of file scale_ocean_sfc_slab.F90.

References scale_const::const_i_lw, scale_const::const_i_sw, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, and scale_grid_index::js.

Referenced by scale_ocean_sfc::ocean_sfc_setup().

286  use scale_grid_index
287  use scale_const, only: &
288  i_sw => const_i_sw, &
289  i_lw => const_i_lw
290  implicit none
291 
292  ! arguments
293  real(RP), intent(out) :: sfc_albedo_t(ia,ja,2) ! tendency of sea surface albedo [0-1]
294  real(RP), intent(in) :: sfc_albedo (ia,ja,2) ! sea surface [0-1]
295  real(RP), intent(in) :: cossza (ia,ja) ! cos(solar zenith angle) [0-1]
296  real(DP), intent(in) :: dt ! delta time
297 
298  ! works
299  real(RP) :: sfc_albedo1(ia,ja,2)
300  real(RP) :: am1
301 
302  real(RP), parameter :: c_ocean_albedo(3) = (/ -0.747900_rp, &
303  -4.677039_rp, &
304  1.583171_rp /)
305 
306  integer :: i, j
307  !---------------------------------------------------------------------------
308 
309  do j = js, je
310  do i = is, ie
311  am1 = max( min( cossza(i,j), 0.961_rp ), 0.0349_rp )
312 
313  ! SFC_albedo1(i,j,I_LW) = 0.5_RP do nothing
314  sfc_albedo1(i,j,i_sw) = exp( ( c_ocean_albedo(3)*am1 + c_ocean_albedo(2) )*am1 + c_ocean_albedo(1) )
315  enddo
316  enddo
317 
318  ! calculate tendency
319  do j = js, je
320  do i = is, ie
321  sfc_albedo_t(i,j,i_lw) = 0.0_rp
322  sfc_albedo_t(i,j,i_sw) = ( sfc_albedo1(i,j,i_sw) - sfc_albedo(i,j,i_sw) ) / dt
323  enddo
324  enddo
325 
326  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:98
integer, parameter, public i_lw
module grid index
integer, parameter, public i_sw
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 ie
end point of inner domain: x, local
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function: