SCALE-RM
scale_ocean_sfc.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_setup
27 
28  abstract interface
29  subroutine ocnsfc( &
30  SST_t, &
31  ZMFLX, &
32  XMFLX, &
33  YMFLX, &
34  SHFLX, &
35  LHFLX, &
36  WHFLX, &
37  U10, &
38  V10, &
39  T2, &
40  Q2, &
41  TMPA, &
42  PRSA, &
43  WA, &
44  UA, &
45  VA, &
46  RHOA, &
47  QVA, &
48  Z1, &
49  PBL, &
50  PRSS, &
51  LWD, &
52  SWD, &
53  TW, &
54  SST, &
55  ALB_LW, &
56  ALB_SW, &
57  Z0M, &
58  Z0H, &
59  Z0E, &
60  dt )
61  use scale_precision
63  implicit none
64 
65  real(RP), intent(out) :: SST_t(IA,JA) ! tendency of sea surface temperature
66  real(RP), intent(out) :: ZMFLX(IA,JA) ! z-momentum flux at the surface [kg/m/s2]
67  real(RP), intent(out) :: XMFLX(IA,JA) ! x-momentum flux at the surface [kg/m/s2]
68  real(RP), intent(out) :: YMFLX(IA,JA) ! y-momentum flux at the surface [kg/m/s2]
69  real(RP), intent(out) :: SHFLX(IA,JA) ! sensible heat flux at the surface [W/m2]
70  real(RP), intent(out) :: LHFLX(IA,JA) ! latent heat flux at the surface [W/m2]
71  real(RP), intent(out) :: WHFLX(IA,JA) ! water heat flux at the surface [W/m2]
72  real(RP), intent(out) :: U10 (IA,JA) ! velocity u at 10m [m/s]
73  real(RP), intent(out) :: V10 (IA,JA) ! velocity v at 10m [m/s]
74  real(RP), intent(out) :: T2 (IA,JA) ! temperature at 2m [K]
75  real(RP), intent(out) :: Q2 (IA,JA) ! water vapor at 2m [kg/kg]
76 
77  real(RP), intent(in) :: TMPA(IA,JA) ! temperature at the lowest atmospheric layer [K]
78  real(RP), intent(in) :: PRSA(IA,JA) ! pressure at the lowest atmospheric layer [Pa]
79  real(RP), intent(in) :: WA (IA,JA) ! velocity w at the lowest atmospheric layer [m/s]
80  real(RP), intent(in) :: UA (IA,JA) ! velocity u at the lowest atmospheric layer [m/s]
81  real(RP), intent(in) :: VA (IA,JA) ! velocity v at the lowest atmospheric layer [m/s]
82  real(RP), intent(in) :: RHOA(IA,JA) ! density at the lowest atmospheric layer [kg/m3]
83  real(RP), intent(in) :: QVA (IA,JA) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
84  real(RP), intent(in) :: Z1 (IA,JA) ! cell center height at the lowest atmospheric layer [m]
85  real(RP), intent(in) :: PBL (IA,JA) ! the top of atmospheric mixing layer [m]
86  real(RP), intent(in) :: PRSS(IA,JA) ! pressure at the surface [Pa]
87  real(RP), intent(in) :: LWD (IA,JA) ! downward long-wave radiation flux at the surface (upward positive) [W/m2]
88  real(RP), intent(in) :: SWD (IA,JA) ! downward short-wave radiation flux at the surface (upward positive) [W/m2]
89 
90  real(RP), intent(in) :: TW (IA,JA) ! water temperature [K]
91  real(RP), intent(in) :: SST (IA,JA) ! sea surface temperature [K]
92  real(RP), intent(in) :: ALB_LW(IA,JA) ! surface albedo for LW (0-1)
93  real(RP), intent(in) :: ALB_SW(IA,JA) ! surface albedo for SW (0-1)
94  real(RP), intent(in) :: Z0M (IA,JA) ! roughness length for momentum [m]
95  real(RP), intent(in) :: Z0H (IA,JA) ! roughness length for heat [m]
96  real(RP), intent(in) :: Z0E (IA,JA) ! roughness length for vapor [m]
97  real(DP), intent(in) :: dt ! delta time
98  end subroutine ocnsfc
99  end interface
100 
101  procedure(ocnsfc), pointer :: ocean_sfc => null()
102 
103  public :: ocean_sfc
104  public :: ocean_sfc_simplealbedo
105 
106  !-----------------------------------------------------------------------------
107  !
108  !++ Public parameters & variables
109  !
110  !-----------------------------------------------------------------------------
111  !
112  !++ Private procedure
113  !
114  !-----------------------------------------------------------------------------
115  !
116  !++ Private parameters & variables
117  !
118  !-----------------------------------------------------------------------------
119 contains
120 
121  subroutine ocean_sfc_setup( OCEAN_TYPE )
122  use scale_process, only: &
124  use scale_ocean_sfc_const, only: &
127  use scale_ocean_sfc_slab, only: &
130  implicit none
131 
132  character(len=*), intent(in) :: ocean_type
133  !---------------------------------------------------------------------------
134 
135  select case( ocean_type )
136  case( 'CONST' )
137  call ocean_sfc_const_setup( ocean_type )
139  case( 'SLAB', 'FILE' )
140  call ocean_sfc_slab_setup( ocean_type )
142  end select
143 
144  return
145  end subroutine ocean_sfc_setup
146 
147  !-----------------------------------------------------------------------------
148  subroutine ocean_sfc_simplealbedo( &
149  SFC_albedo_t, &
150  SFC_albedo, &
151  cosSZA, &
152  dt )
154  use scale_const, only: &
155  i_sw => const_i_sw, &
156  i_lw => const_i_lw
157  implicit none
158 
159  ! arguments
160  real(RP), intent(out) :: sfc_albedo_t(ia,ja,2) ! tendency of sea surface albedo (0-1)
161  real(RP), intent(in) :: sfc_albedo (ia,ja,2) ! sea surface (0-1)
162  real(RP), intent(in) :: cossza (ia,ja) ! cos(solar zenith angle) (0-1)
163  real(DP), intent(in) :: dt ! delta time
164 
165  ! works
166  real(RP) :: sfc_albedo1(ia,ja,2)
167  real(RP) :: am1
168 
169  real(RP), parameter :: c_ocean_albedo(3) = (/ -0.747900_rp, &
170  -4.677039_rp, &
171  1.583171_rp /)
172 
173  integer :: i, j
174  !---------------------------------------------------------------------------
175 
176  do j = js, je
177  do i = is, ie
178  am1 = max( min( cossza(i,j), 0.961_rp ), 0.0349_rp )
179 
180  ! SFC_albedo1(i,j,I_LW) = 0.5_RP do nothing
181  sfc_albedo1(i,j,i_sw) = exp( ( c_ocean_albedo(3)*am1 + c_ocean_albedo(2) )*am1 + c_ocean_albedo(1) )
182  enddo
183  enddo
184 
185  ! calculate tendency
186  do j = js, je
187  do i = is, ie
188  sfc_albedo_t(i,j,i_lw) = 0.0_rp
189  sfc_albedo_t(i,j,i_sw) = ( sfc_albedo1(i,j,i_sw) - sfc_albedo(i,j,i_sw) ) / dt
190  enddo
191  enddo
192 
193  return
194  end subroutine ocean_sfc_simplealbedo
195 
196 end module scale_ocean_sfc
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:95
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 prc_mpistop
Abort MPI.
subroutine, public ocean_sfc_setup(OCEAN_TYPE)
subroutine, public ocean_sfc_simplealbedo(SFC_albedo_t, SFC_albedo, cosSZA, dt)
procedure(ocnsfc), pointer, public ocean_sfc
module STDIO
Definition: scale_stdio.F90:12
module OCEAN / Surface flux with slab ocean model
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
integer, public js
start point of inner domain: y, local
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module OCEAN / Surface fluxes
subroutine, public ocean_sfc_slab_setup(OCEAN_TYPE)
Setup.
integer, public ie
end point of inner domain: x, local
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:96
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 ja
of whole cells: y, local, with HALO