SCALE-RM
scalelib
src
ocean
scale_ocean_sfc.F90
Go to the documentation of this file.
1
!-------------------------------------------------------------------------------
10
!-------------------------------------------------------------------------------
11
module
scale_ocean_sfc
12
!-----------------------------------------------------------------------------
13
!
14
!++ used modules
15
!
16
use
scale_precision
17
use
scale_stdio
18
use
scale_grid_index
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
62
use
scale_grid_index
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/m2/s]
67
real(RP)
,
intent(out)
:: xmflx(
ia
,
ja
)
! x-momentum flux at the surface [kg/m2/s]
68
real(RP)
,
intent(out)
:: ymflx(
ia
,
ja
)
! y-momentum flux at the surface [kg/m2/s]
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
100
subroutine
alb
( &
101
SFC_albedo_t, &
102
SFC_albedo, &
103
cosSZA, &
104
dt )
105
use
scale_precision
106
use
scale_grid_index
107
implicit none
108
109
real(RP)
,
intent(out)
:: sfc_albedo_t(
ia
,
ja
,2)
! tendency of sea surface albedo [0-1]
110
real(RP)
,
intent(in)
:: sfc_albedo (
ia
,
ja
,2)
! sea surface [0-1]
111
real(RP)
,
intent(in)
:: cossza (
ia
,
ja
)
! cos(solar zenith angle) [0-1]
112
real(DP)
,
intent(in)
:: dt
! delta time
113
end subroutine
alb
114
end interface
115
procedure
(ocnsfc),
pointer
::
ocean_sfc
=> null()
116
procedure
(
alb
),
pointer
::
ocean_sfc_simplealbedo
=> null()
117
public
::
ocean_sfc
118
public
::
ocean_sfc_simplealbedo
119
120
!-----------------------------------------------------------------------------
121
!
122
!++ Public parameters & variables
123
!
124
!-----------------------------------------------------------------------------
125
!
126
!++ Private procedure
127
!
128
!-----------------------------------------------------------------------------
129
!
130
!++ Private parameters & variables
131
!
132
!-----------------------------------------------------------------------------
133
contains
134
135
subroutine
ocean_sfc_setup
( OCEAN_TYPE )
136
use
scale_process
, only
: &
137
prc_mpistop
138
use
scale_ocean_sfc_slab
, only
: &
139
ocean_sfc_slab_setup
, &
140
ocean_sfc_slab
, &
141
ocean_sfc_slab_simplealbedo
142
implicit none
143
144
character(len=*)
,
intent(in)
:: OCEAN_TYPE
145
!---------------------------------------------------------------------------
146
147
select case
( ocean_type )
148
case
(
'CONST'
,
'FILE'
)
149
call
ocean_sfc_slab_setup
( ocean_type )
150
ocean_sfc
=>
ocean_sfc_slab
151
ocean_sfc_simplealbedo
=>
ocean_sfc_slab_simplealbedo
152
case
(
'SLAB'
)
153
call
ocean_sfc_slab_setup
( ocean_type )
154
ocean_sfc
=>
ocean_sfc_slab
155
ocean_sfc_simplealbedo
=>
ocean_sfc_slab_simplealbedo
156
end select
157
158
return
159
end subroutine
ocean_sfc_setup
160
161
end module
scale_ocean_sfc
scale_ocean_sfc::ocean_sfc_simplealbedo
procedure(alb), pointer, public ocean_sfc_simplealbedo
Definition:
scale_ocean_sfc.F90:116
scale_ocean_sfc_slab::ocean_sfc_slab
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)
Definition:
scale_ocean_sfc_slab.F90:124
scale_process::prc_mpistop
subroutine, public prc_mpistop
Abort MPI.
Definition:
scale_process.F90:234
scale_ocean_sfc::ocean_sfc_setup
subroutine, public ocean_sfc_setup(OCEAN_TYPE)
Definition:
scale_ocean_sfc.F90:136
scale_ocean_sfc::ocean_sfc
procedure(ocnsfc), pointer, public ocean_sfc
Definition:
scale_ocean_sfc.F90:115
scale_stdio
module STDIO
Definition:
scale_stdio.F90:12
scale_ocean_sfc_slab::ocean_sfc_slab_simplealbedo
subroutine, public ocean_sfc_slab_simplealbedo(SFC_albedo_t, SFC_albedo, cosSZA, dt)
Definition:
scale_ocean_sfc_slab.F90:286
scale_ocean_sfc_slab
module OCEAN / Surface flux with slab ocean model
Definition:
scale_ocean_sfc_slab.F90:11
scale_grid_index
module grid index
Definition:
scale_grid_index.F90:14
scale_grid_index::ia
integer, public ia
of x whole cells (local, with HALO)
Definition:
scale_grid_index.F90:55
scale_process
module PROCESS
Definition:
scale_process.F90:14
scale_ocean_sfc::alb
Definition:
scale_ocean_sfc.F90:100
scale_ocean_sfc
module OCEAN / Surface fluxes
Definition:
scale_ocean_sfc.F90:11
scale_ocean_sfc_slab::ocean_sfc_slab_setup
subroutine, public ocean_sfc_slab_setup(OCEAN_TYPE)
Setup.
Definition:
scale_ocean_sfc_slab.F90:51
scale_precision
module PRECISION
Definition:
scale_precision.F90:13
scale_grid_index::ja
integer, public ja
of y whole cells (local, with HALO)
Definition:
scale_grid_index.F90:56
Generated by
1.8.13