Library
Documenting the user interface.
Wall stress model
Walrus.WallStressModel
— ModuleWallStressModel
A wall stress model for LES simulation with default parameters similar to that proposed in Schumann (1975), Hartel (1996), Piomelli et al. (1989), and Taylor and Sarkar (2007).
Walrus.WallStressModel.WallStress
— MethodWallStress(; von_Karman_constant = 0.4,
kinematic_viscosity = 1e-6,
B = 5.2,
precomputed_friction_velocities = false,
precompute_speeds = [0:25/100000:25;],
grid = nothing,
arch = isnothing(grid) ? CPU() : architecture(grid))
Returns a wall stress model for LES simulation with default parameters similar to that proposed in Schumann (1975), Hartel (1996), Piomelli et al. (1989), and Taylor and Sarkar (2007).
Friction velocities will be precomputed at precompute_speeds
if precomputed_friction_velocities
is true and grid
is provided.
Keyword Arguments
von_Karman_constant
: the von Karman wall stress constantkinematic_viscosity
: kinematic viscosity of the water above the wallB
: wall stress constantprecomputed_friction_velocities
: precompute friction velocities?precompute_speeds
: bottom water speeds to precompute friction velocities for, this should encompas the range of speeds possible in your simulationgrid
: the grid to precompute the friction velocities forarch
: architecture to adapt precomputed friction velocities for
Example
julia> using Walrus: WallStress
julia> using Oceananigans
julia> wall_stress = WallStress()
(::WallStress{Float64, Nothing}) (generic function with 1 method)
julia> boundary_conditions = (u = FieldBoundaryConditions(bottom = FluxBoundaryCondition(wall_stress, discrete_form = true, parameters = Val(:x))),
v = FieldBoundaryConditions(bottom = FluxBoundaryCondition(wall_stress, discrete_form = true, parameters = Val(:y))))
(u = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: FluxBoundaryCondition: DiscreteBoundaryFunction (::WallStress{Float64, Nothing}) with parameters Val{:x}
├── top: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing), v = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: FluxBoundaryCondition: DiscreteBoundaryFunction (::WallStress{Float64, Nothing}) with parameters Val{:y}
├── top: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing))
Walrus.WallStressModel.WallStressBoundaryConditions
— MethodWallStressBoundaryConditions(; von_Karman_constant = 0.4,
kinematic_viscosity = 1e-6,
B = 5.2,
precompute_speeds = [0:25/100000:25;],
grid = nothing)
Convenience constructor to setup WallStress
boundary conditions.
Keyword Arguments
von_Karman_constant
: the von Karman wall stress constantkinematic_viscosity
: kinematic viscosity of the water above the wallB
: wall stress constantprecomputed_friction_velocities
: precompute friction velocities?precompute_speeds
: bottom water speeds to precompute friction velocities for, this should encompas the range of speeds possible in your simulationgrid
: the grid to precompute the friction velocities for
Example
julia> using Walrus: WallStressBoundaryConditions
julia> using Oceananigans
julia> stress_boundary_conditions = WallStressBoundaryConditions()
(u = FluxBoundaryCondition: DiscreteBoundaryFunction (::WallStress{Float64, Nothing}) with parameters Val{:x}, v = FluxBoundaryCondition: DiscreteBoundaryFunction (::WallStress{Float64, Nothing}) with parameters Val{:y})
julia> boundary_conditions = (u = FieldBoundaryConditions(bottom = stress_boundary_conditions.u),
v = FieldBoundaryConditions(bottom = stress_boundary_conditions.v))
(u = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: FluxBoundaryCondition: DiscreteBoundaryFunction (::WallStress{Float64, Nothing}) with parameters Val{:x}
├── top: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing), v = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: FluxBoundaryCondition: DiscreteBoundaryFunction (::WallStress{Float64, Nothing}) with parameters Val{:y}
├── top: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing))
Radiative transfer models
Walrus.RadiativeTransfer
— ModuleRadiativeTransfer
Includes models for ratiative transfer through water which can induce body heating
Walrus.RadiativeTransfer.HomogeneousBodyHeating
— TypeHomogeneousBodyHeating
A model for single band light attenuation which heats the water in the form:
\[I(x, y, z) = I_0(x, y) * \exp\left(-\alpha z\right),\]
where $I$ is the radiation intensity and $\alpha$ is the attenuation coefficient. This heats the water like $\frac{\partial T(x, y, z)}{\partial t} = \frac{I(x, y, z)A}{c^p\rho}$ where $A$ is the area of the cell, $c^p$ is the specific heat capacity and $\rho$ is the water density.
Walrus.RadiativeTransfer.HomogeneousBodyHeating
— MethodHomogeneousBodyHeating(; surface_flux,
water_attenuation_coefficient = 1.8,
water_heat_capacity = 3991.0, # J K⁻¹ kg⁻¹
water_density = 1026.0) # kg m⁻³
Creates a model in which a surface_flux
(W / m²) is attenuated by and heats the water. This interacts with Oceananigans as a body forcing.
Keyword Arguments
surface_flux
(required): a function returning the surface radiaiton flux in the formsurface_flux(x, y, t)
or single valuewater_attenuation_coefficient
: the radiation attenuation coefficient of the waterwater_heat_capacity
: the specific heat capacity of the waterwater_density
: density of the water
Example
julia> using Walrus: HomogeneousBodyHeating
julia> using Oceananigans
julia> grid = RectilinearGrid(size = (128, 128, 128), extent = (1000, 1000, 1000))
128×128×128 RectilinearGrid{Float64, Oceananigans.Grids.Periodic, Oceananigans.Grids.Periodic, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 3×3×3 halo
├── Periodic x ∈ [0.0, 1000.0) regularly spaced with Δx=7.8125
├── Periodic y ∈ [0.0, 1000.0) regularly spaced with Δy=7.8125
└── Bounded z ∈ [-1000.0, 0.0] regularly spaced with Δz=7.8125
julia> body_heating = HomogeneousBodyHeating(; surface_flux = (x, y, t) -> 100)
(::HomogeneousBodyHeating{Float64, var"#1#2"}) (generic function with 1 method)
julia> model = NonhydrostaticModel(; grid, forcing = (; T = Forcing(body_heating, discrete_form=true)), tracers = :T)
NonhydrostaticModel{CPU, RectilinearGrid}(time = 0 seconds, iteration = 0)
├── grid: 128×128×128 RectilinearGrid{Float64, Oceananigans.Grids.Periodic, Oceananigans.Grids.Periodic, Oceananigans.Grids.Bounded} on Oceananigans.Architectures.CPU with 3×3×3 halo
├── timestepper: QuasiAdamsBashforth2TimeStepper
├── advection scheme: Centered reconstruction order 2
├── tracers: T
├── closure: Nothing
├── buoyancy: Nothing
└── coriolis: Nothing
Tidal forcing models
Walrus.TidalForcings
— ModuleTidalForcing
Provides quick setup of tidal forcing
Walrus.TidalForcings.Tide
— TypeTide(; x_amplitude,
y_amplitude,
period = 12.3782216453hours,
nodal_time = 0.,
x_lag = 0.,
y_lag = 0.,
coriolis = nothing)
Sets up a model of tidal forcing with default parameters of an M2 tide.
Keyword Arguments
x_amplitude
: the tidal amplitude in the x directiony_amplitude
: the tidal amplitude in the x directionperiod
: the tidal period (defaults to that of an M2 tide)nodal_time
: the time at which peak flow occursx_lag
: the phase lag for the tidal component in the x directiony_lag
: the phase lag for the tidal component in the y directioncoriolis
: a model for the coriolis parameter
Example
julia> using Walrus: Tide
julia> using Oceananigans
julia> tide = Tide(x_amplitude = 0.1, y_amplitude = 0.)
(::Tide{Float64, Nothing}) (generic function with 2 methods)
julia> forcing = (u = Forcing(tide, parameters = Val(:x), discrete_form = true),
v = Forcing(tide, parameters = Val(:y), discrete_form = true))
(u = DiscreteForcing{Val{:x}}
├── func: (::Tide{Float64, Nothing}) (generic function with 2 methods)
└── parameters: Val{:x}(), v = DiscreteForcing{Val{:y}}
├── func: (::Tide{Float64, Nothing}) (generic function with 2 methods)
└── parameters: Val{:y}())
Walrus.TidalForcings.TidalForcing
— MethodTidalForcing(; x_amplitude,
y_amplitude,
period = 12.3782216453hours,
nodal_time = 0.,
x_lag = 0.,
y_lag = 0.,
coriolis = nothing)
A convenience constructor for Tide
which returns the forcings pre wrapped.
Keyword Arguments
x_amplitude
: the tidal amplitude in the x directiony_amplitude
: the tidal amplitude in the x directionperiod
: the tidal period (defaults to that of an M2 tide)nodal_time
: the time at which peak flow occursx_lag
: the phase lag for the tidal component in the x directiony_lag
: the phase lag for the tidal component in the y directioncoriolis
: a model for the coriolis parameter
Example
julia> using Walrus: TidalForcing
julia> tidal_forcing = TidalForcing(x_amplitude = 0.1, y_amplitude = 0.)
(u = DiscreteForcing{Val{:x}}
├── func: (::Tide{Float64, Nothing}) (generic function with 2 methods)
└── parameters: Val{:x}(), v = DiscreteForcing{Val{:y}}
├── func: (::Tide{Float64, Nothing}) (generic function with 2 methods)
└── parameters: Val{:y}())
Wind stress model
Walrus.WindStressModel.LogarithmicNeutralWind
— MethodLogarithmicNeutralWind(; monin_obukhov_stability_length = 0.4
charnock_coefficient = 0.014
air_kinematic_viscosity = 1.488e-5
gravity_wave_coefficient = 0.11
gravity = g_Earth,
precomputed_roughness_length = false,
precompute_wind_speeds = [0:25/100000:25;],
arch = CPU())
Returns a LogarithmicNeutralWind
parameterisation for the surface drag coefficient
$C_d$ is parameterised as,
\[C_d = \left(\frac{\kappa}{\log{\frac{10}{z_0}}}\right)^2,\]
where $\kappa$ is the Monin‐Obukhov stability length and $z_0$ is the velocity roughness length. This is the roughness length scale which logarithmically brings the relative velocity to zero at the surface, i.e.
\[U=\frac{u\star}{\kappa}\log\frac{z}{z_0},\]
where $u\star$ is the friction velocity. Additionally $z_0$ is given as,
\[z_0=b\frac{\nu}{u\star} + \frac{a_c}{g}u\star^2,\]
where $\nu$ is the kinematic viscosity of air and g is the acceleration of gravity.
This model itterativly solves these equations to find $z_0$. Alternativly, if the flag precomputed_roughness_length
is set to they are pre computed at precompute_wind_speeds
between which $z_0$ is then interpolated during run time. Precomputed velocities are converted to appropriate types for arch
(i.e. CPU()
or GPU()
)
This parameterisaion is described in Smith (1988)
Walrus.WindStressModel.WindStress
— MethodWindStress(; reference_wind_speed,
reference_wind_direction,
drag_coefficient = LogarithmicNeutralWind(),
air_density = 1.225,
water_density = 1026.)
Returns a wind stress model where the stress is given by,
\[\frac{\tau}{\rho_o} = \rho_aC_dSU_{x/y},\]
where $\rho_o$ is the water density, $\rho_a$ is the air density, $C_d$ is the drag coefficient, $U_{x/y}$ are the x and y components of relative wind speed, and $S=\sqrt{U_x^2+U_y^2}$.
$C_d$ is calculated from a parameterisation, by default this is a "log neutral" wind parameterisation with velocity roughness length parameterisaion like Smith (1988).
In the default configuration this is the same as described in Fairall et al. (2011).
Keyword Arguments
reference_wind_speed
(required): a function returning the (10m neutral) wind speed in the formreference_wind_speed(x, y, t)
or single valuereference_wind_direction
(required): a function returning the (10m neutral) wind direction in the formreference_wind_direction(x, y, t)
or single valuedrag_coefficient
: the drag coefficient parameterisationair_density
: air density in kg/m³water_density
: water density in kg/m³
Example
julia> using Walrus: WindStress
julia> using Oceananigans
julia> reference_wind_speed = 0.1
0.1
julia> reference_wind_direction = 0.
0.0
julia> wind_stress = WindStress(; reference_wind_speed, reference_wind_direction)
(::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) (generic function with 2 methods)
julia> boundary_conditions = (u = FieldBoundaryConditions(top = FluxBoundaryCondition(wind_stress, parameters = Val(:x))),
v = FieldBoundaryConditions(top = FluxBoundaryCondition(wind_stress, parameters = Val(:y))))
(u = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── top: FluxBoundaryCondition: ContinuousBoundaryFunction (::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) at (Nothing, Nothing, Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing), v = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── top: FluxBoundaryCondition: ContinuousBoundaryFunction (::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) at (Nothing, Nothing, Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing))
Walrus.WindStressModel.WindStressBoundaryConditions
— MethodWindStressBoundaryConditions(; reference_wind_speed,
reference_wind_direction,
drag_coefficient = LogarithmicNeutralWind(),
air_density = 1.225,
water_density = 1026.)
Convenience constructor to setup WindStress
boundary conditions.
Keyword Arguments
reference_wind_speed
(required): a function returning the (10m neutral) wind speed in the formreference_wind_speed(x, y, t)
or single valuereference_wind_direction
(required): a function returning the (10m neutral) wind direction in the formreference_wind_direction(x, y, t)
or single valuedrag_coefficient
: the drag coefficient parameterisationair_density
: air density in kg/m³water_density
: water density in kg/m³
Example
julia> using Walrus: WindStressBoundaryConditions
julia> using Oceananigans
julia> wind_stress_boundary_conditions = WindStressBoundaryConditions(; reference_wind_speed = 0.1, reference_wind_direction = 90.)
(u = FluxBoundaryCondition: DiscreteBoundaryFunction (::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) with parameters Val{:x}, v = FluxBoundaryCondition: DiscreteBoundaryFunction (::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) with parameters Val{:y})
julia> boundary_conditions = (u = FieldBoundaryConditions(top = wind_stress_boundary_conditions.u),
v = FieldBoundaryConditions(top = wind_stress_boundary_conditions.v))
(u = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── top: FluxBoundaryCondition: DiscreteBoundaryFunction (::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) with parameters Val{:x}
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing), v = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── top: FluxBoundaryCondition: DiscreteBoundaryFunction (::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) with parameters Val{:y}
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing))
Walrus.WindStressModel.find_velocity_roughness_length
— Methodfind_velocity_roughness_length(wind_speed, reference_height, params)
A function that finds the velocity roughness length for the LogarithmicNeutralWind
drag coefficient model.
This will sometimes fail as the function is not well behaved at either low reference heights (it has been tuned for 10m wind), or high (⪆ 20 m/s).
Surface heat exchange model
Walrus.SurfaceHeatingModel.SurfaceHeatExchange
— MethodSurfaceHeatExchange(; wind_stress,
air_temperature = 18, # °C
latent_heat_vaporisation = EmpiricalLatentHeatVaporisation(),
vapour_pressure = AugustRocheMagnusVapourPressure(),
water_specific_heat_capacity = 3991., # J / K / kg
water_density = 1026., # kg / m³
air_specific_heat_capacity = 1003.5, # J / K / kg
air_density = 1.204, # kg
air_water_mixing_ratio = 0.001, # kg / kg
stephan_boltzman_constant = 5.670374419e-8) # W / K⁴
Specifies surface heat exhange in the form:
\[Q = Qᵢᵣ + Qₛ + Qₗ,\]
where $Qᵢᵣ$ is the heat flux due to long wave (infra red) radiation, $Qₛ$ is the sensible heat flux, and $Qₗ$ is the latent heat flux. Notably the short wave radiation flux is neglegted here as it is assumed to penatrate far enough into the water that it is treated separately by HomogeneousBodyHeating
(we therefore are also assuming that the short wave penitration is suitably short that it is neglected).
The short wave term is given by the Stephan-Boltzman equation so:
\[Qᵢᵣ = σ(T⁴ - Tₐ⁴),\]
where $σ$ is the Stephan-Boltzman constant, $T$ is the ocean temperature, and $Tₐ$ is the 2m air temperature.
The sensible and latent heat flux's are given by the bluk parameterisations described in Fairall et al. (2011) and are given by:
\[Qₛ = ρₐcₚₐCₕS(T - Tₐ),\]
and $Qₗ = ρₐLₑCₕS(q(T) - qₐ)$, where $ρₐ$ is the density of the air, $cₚₐ$ is the specific heat capacity of air, $Cₕ$ is the heat transfer coefficient, $S$ is the wind speed, $Lₑ$ is the latent heat of vaporizaion parameterised by latent_heat_vaporisation
, $q(T)$ is the saturation vapour pressure of water and is parameterised by vapour_pressure
.
The heat transfer coefficients are given by the same parameterisation as the WindStress
. For example for the LogarithmicNeutralWind
the transfer coefficient is given as: $C_h = \frac{\kappa}{\log{\frac{2}{z_0}}}\frac{\kappa}{\log{\frac{2}{z_{ot}}}}$, where $z_{ot}$ is the sclar roughness parameter given by: $z_{ot} = \min\left(1.15\cdot10^{-4}, 5.5\cdot10^{-5}R_r^{-0.6}\right)$, where $R_r$ is the roughness reynolds number given as $R_r = \frac{u\star z_0}{\nu}$ where $\nu$ is the kinematic viscosity of air.
The heat flux is then given by:
\[F = \frac{Q}{\rho_oc_{po}},\]
where $\rho_o$ and $c_{po}$ are the density and specific heat capacity of water.
(Note: we will retain the Oceananigans convention that negative heat flux at a top boundary increases temeprature).
Keyword Arguments
wind_stress
: wind stress modelair_temperature
: the air temperature in °C as a function with signature(x, y, t)
or a constantlatent_heat_vaporisation
: latent heat of vaporisation in J / kgvapour_pressure
: parameterisation for saturation vapour pressure in waterwater_specific_heat_capacity
: the specific heat capacity of water in J / K / kgwater_density
: water density in kg / m³air_specific_heat_capacit
: the specific heat capacity of air in J / K / kgair_density
: air density in kg / m³air_water_mixing_ratio
: water content of air in kg / kgstephan_boltzman_constant
: the Stephan-Boltzman constant in W / K⁴
Example
julia> using Walrus
julia> using Oceananigans
julia> wind_stress = WindStress(; reference_wind_speed = 0., reference_wind_direction = 90.)
(::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) (generic function with 2 methods)
julia> surface_heat_exchange = SurfaceHeatExchange(; wind_stress)
(::SurfaceHeatExchange{WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}, Int64, Walrus.SurfaceHeatingModel.EmpiricalLatentHeatVaporisation{Float64}, Walrus.SurfaceHeatingModel.AugustRocheMagnusVapourPressure{Float64}, Float64}) (generic function with 1 method)
julia> boundary_conditions = (; T = FieldBoundaryConditions(top = FluxBoundaryCondition(surface_heat_exchange, field_dependencies = (:T, :u, :v))))
(T = Oceananigans.FieldBoundaryConditions, with boundary conditions
├── west: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── east: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── south: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── north: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── bottom: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing)
├── top: FluxBoundaryCondition: ContinuousBoundaryFunction (::SurfaceHeatExchange{WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}, Int64, Walrus.SurfaceHeatingModel.EmpiricalLatentHeatVaporisation{Float64}, Walrus.SurfaceHeatingModel.AugustRocheMagnusVapourPressure{Float64}, Float64}) at (Nothing, Nothing, Nothing)
└── immersed: DefaultBoundaryCondition (FluxBoundaryCondition: Nothing),)
Walrus.SurfaceHeatingModel.SurfaceHeatExchangeBoundaryCondition
— MethodSurfaceHeatExchange(; wind_stress,
air_temperature = 18, # °C
latent_heat_vaporisation = EmpiricalLatentHeatVaporisation(),
vapour_pressure = AugustRocheMagnusVapourPressure(),
water_specific_heat_capacity = 3991., # J / K / kg
water_density = 1026., # kg / m³
air_specific_heat_capacity = 1003.5, # J / K / kg
air_density = 1.204, # kg
air_water_mixing_ratio = 0.001, # kg / kg
stephan_boltzman_constant = 5.670374419e-8) # W / K⁴
A convenience constructor returning SurfaceHeatExchange
as a boundary condition
Keyword Arguments
wind_stress
: wind stress modelair_temperature
: the air temperature in °C as a function with signature(x, y, t)
or a constantlatent_heat_vaporisation
: latent heat of vaporisation in J / kgvapour_pressure
: parameterisation for saturation vapour pressure in waterwater_specific_heat_capacity
: the specific heat capacity of water in J / K / kgwater_density
: water density in kg / m³air_specific_heat_capacit
: the specific heat capacity of air in J / K / kgair_density
: air density in kg / m³air_water_mixing_ratio
: water content of air in kg / kgstephan_boltzman_constant
: the Stephan-Boltzman constant in W / K⁴
Example
julia> using Walrus
julia> wind_stress = WindStress(; reference_wind_speed = 0., reference_wind_direction = 90.)
(::WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}) (generic function with 2 methods)
julia> surface_heat_exchange = SurfaceHeatExchangeBoundaryCondition(; wind_stress)
FluxBoundaryCondition: DiscreteBoundaryFunction with (::SurfaceHeatExchange{WindStress{Float64, Float64, LogarithmicNeutralWind{Float64, Nothing}, Float64}, Int64, Walrus.SurfaceHeatingModel.EmpiricalLatentHeatVaporisation{Float64}, Walrus.SurfaceHeatingModel.AugustRocheMagnusVapourPressure{Float64}, Float64})