SCALE-RM
scale_ocean_sfc_const.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
19  !-----------------------------------------------------------------------------
20  implicit none
21  private
22  !-----------------------------------------------------------------------------
23  !
24  !++ Public procedure
25  !
26  public :: ocean_sfc_const_setup
27  public :: ocean_sfc_const
28 
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public parameters & variables
32  !
33  !-----------------------------------------------------------------------------
34  !
35  !++ Private procedure
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private parameters & variables
40  !
41  !-----------------------------------------------------------------------------
42 contains
43  !-----------------------------------------------------------------------------
45  subroutine ocean_sfc_const_setup( OCEAN_TYPE )
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
55  end subroutine ocean_sfc_const_setup
56 
57  !-----------------------------------------------------------------------------
58  subroutine ocean_sfc_const( &
59  SST_t, & ! [OUT]
60  ZMFLX, & ! [OUT]
61  XMFLX, & ! [OUT]
62  YMFLX, & ! [OUT]
63  SHFLX, & ! [OUT]
64  LHFLX, & ! [OUT]
65  WHFLX, & ! [OUT]
66  U10, & ! [OUT]
67  V10, & ! [OUT]
68  T2, & ! [OUT]
69  Q2, & ! [OUT]
70  TMPA, & ! [IN]
71  PRSA, & ! [IN]
72  WA, & ! [IN]
73  UA, & ! [IN]
74  VA, & ! [IN]
75  RHOA, & ! [IN]
76  QVA, & ! [IN]
77  Z1, & ! [IN]
78  PBL, & ! [IN]
79  PRSS, & ! [IN]
80  LWD, & ! [IN]
81  SWD, & ! [IN]
82  TW, & ! [IN]
83  SST, & ! [IN]
84  ALB_LW, & ! [IN]
85  ALB_SW, & ! [IN]
86  Z0M, & ! [IN]
87  Z0H, & ! [IN]
88  Z0E, & ! [IN]
89  dt ) ! [IN]
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: &
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
249  end subroutine ocean_sfc_const
250 
251 end module scale_ocean_sfc_const
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.
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
Definition: scale_const.F90:51
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module STDIO
Definition: scale_stdio.F90:12
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module OCEAN / Surface flux with constant ocean model
subroutine, public ocean_sfc_const_setup(OCEAN_TYPE)
Setup.
module grid index
integer, public ia
of whole cells: x, local, with HALO
module LANDUSE
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
procedure(bc), pointer, public bulkflux
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
integer, public js
start point of inner domain: y, local
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public ie
end point of inner domain: x, local
module Surface bulk flux
module PRECISION
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)
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of whole cells: y, local, with HALO