Skip to content

Commit

Permalink
Merge pull request #269 from SpeedyWeather/mk/heldsuarez
Browse files Browse the repository at this point in the history
Boundary layer scheme from Held and Suarez
  • Loading branch information
milankl committed Mar 27, 2023
2 parents 212d19b + f71cc86 commit 12fd32e
Show file tree
Hide file tree
Showing 10 changed files with 238 additions and 180 deletions.
284 changes: 146 additions & 138 deletions src/SpeedyWeather.jl
Original file line number Diff line number Diff line change
@@ -1,143 +1,151 @@
module SpeedyWeather

# STRUCTURE
import Parameters: @unpack

# NUMERICS
import Random
import FastGaussQuadrature
import LinearAlgebra: LinearAlgebra, Diagonal

# GPU, PARALLEL
import Base.Threads: Threads, @threads
import FLoops: FLoops, @floop
import KernelAbstractions
import CUDA
import CUDAKernels
import Adapt: Adapt, adapt, adapt_structure

# INPUT OUTPUT
import TOML
import Dates: Dates, DateTime
import Printf: @sprintf
import NetCDF: NetCDF, NcFile, NcDim, NcVar
import JLD2: jldopen
import CodecZlib
import BitInformation: round, round!
import UnicodePlots
import ProgressMeter

# EXPORT MAIN INTERFACE TO SPEEDY
export run_speedy,
run_speedy!,
initialize_speedy

# EXPORT MODELS
export Barotropic,
BarotropicModel,
ShallowWater,
ShallowWaterModel,
PrimitiveEquation,
PrimitiveDryCore,
PrimitiveWetCore,
PrimitiveDryCoreModel,
PrimitiveWetCoreModel

# EXPORT GRIDS
export LowerTriangularMatrix,
FullClenshawGrid,
FullGaussianGrid,
FullHEALPixGrid,
FullOctaHEALPixGrid,
OctahedralGaussianGrid,
OctahedralClenshawGrid,
HEALPixGrid,
OctaHEALPixGrid

# EXPORT OROGRAPHIES
export NoOrography,
EarthOrography,
ZonalRidge

# EXPORT INITIAL CONDITIONS
export StartFromFile,
StartFromRest,
ZonalJet,
ZonalWind,
StartWithVorticity

# EXPORT STRUCTS
export Parameters,
DynamicsConstants,
ParameterizationConstants,
Geometry,
SpectralTransform,
Boundaries,
PrognosticVariables,
DiagnosticVariables,
ColumnVariables

# EXPORT SPECTRAL FUNCTIONS
export spectral,
gridded,
spectral_truncation

include("utility_functions.jl")
include("abstract_types.jl")

# LowerTriangularMatrices for spherical harmonics
export LowerTriangularMatrices
include("LowerTriangularMatrices/LowerTriangularMatrices.jl")
using .LowerTriangularMatrices

# RingGrids
export RingGrids
include("RingGrids/RingGrids.jl")
using .RingGrids

# SpeedyTransforms
export SpeedyTransforms
include("SpeedyTransforms/SpeedyTransforms.jl")
using .SpeedyTransforms
# STRUCTURE
import Parameters: @unpack

# NUMERICS
import Random
import FastGaussQuadrature
import LinearAlgebra: LinearAlgebra, Diagonal

# GPU, PARALLEL
import Base.Threads: Threads, @threads
import FLoops: FLoops, @floop
import KernelAbstractions
import CUDA
import CUDAKernels
import Adapt: Adapt, adapt, adapt_structure

# INPUT OUTPUT
import TOML
import Dates: Dates, DateTime
import Printf: @sprintf
import NetCDF: NetCDF, NcFile, NcDim, NcVar
import JLD2: jldopen
import CodecZlib
import BitInformation: round, round!
import UnicodePlots
import ProgressMeter

# EXPORT MAIN INTERFACE TO SPEEDY
export run_speedy,
run_speedy!,
initialize_speedy

# EXPORT MODELS
export Barotropic,
BarotropicModel,
ShallowWater,
ShallowWaterModel,
PrimitiveEquation,
PrimitiveDryCore,
PrimitiveWetCore,
PrimitiveDryCoreModel,
PrimitiveWetCoreModel

# EXPORT GRIDS
export LowerTriangularMatrix,
FullClenshawGrid,
FullGaussianGrid,
FullHEALPixGrid,
FullOctaHEALPixGrid,
OctahedralGaussianGrid,
OctahedralClenshawGrid,
HEALPixGrid,
OctaHEALPixGrid

# EXPORT OROGRAPHIES
export NoOrography,
EarthOrography,
ZonalRidge

# EXPORT INITIAL CONDITIONS
export StartFromFile,
StartFromRest,
ZonalJet,
ZonalJetCoefs,
ZonalWind,
ZonalWindCoefs,
StartWithVorticity

# EXPORT BOUNDARY LAYER SCHEMES
export BoundaryLayer,
NoBoundaryLayer,
LinearDrag

# EXPORT STRUCTS
export Parameters,
DynamicsConstants,
ParameterizationConstants,
Geometry,
SpectralTransform,
Boundaries,
PrognosticVariables,
DiagnosticVariables,
ColumnVariables

# EXPORT SPECTRAL FUNCTIONS
export spectral,
gridded,
spectral_truncation

include("gpu.jl") # defines utility for GPU / KernelAbstractions
include("default_parameters.jl") # defines Parameters

# DYNAMICS
include("dynamics/constants.jl") # defines DynamicsConstants
include("dynamics/geometry.jl") # defines Geometry
include("dynamics/boundaries.jl") # defines Boundaries
include("dynamics/define_diffusion.jl") # defines HorizontalDiffusion
include("dynamics/define_implicit.jl") # defines ImplicitShallowWater, ImplicitPrimitiveEq
include("dynamics/parameter_structs.jl") # defines GenLogisticCoefs
include("dynamics/models.jl") # defines ModelSetups
include("dynamics/prognostic_variables.jl") # defines PrognosticVariables
include("dynamics/diagnostic_variables.jl") # defines DiagnosticVariables
include("dynamics/initial_conditions.jl")
include("dynamics/scaling.jl")
include("dynamics/geopotential.jl")
include("dynamics/tendencies_dynamics.jl")
include("dynamics/tendencies.jl")
include("dynamics/implicit.jl")
include("dynamics/diffusion.jl")
include("dynamics/time_integration.jl")

# PHYSICS
include("physics/constants.jl") # defines Parameterization Constants
include("physics/parameter_structs.jl") # defines MagnusCoefs, RadiationCoefs
include("physics/column_variables.jl") # defines ColumnVariables
include("physics/thermodynamics.jl")
include("physics/tendencies_parametrizations.jl")
include("physics/convection.jl")
include("physics/large_scale_condensation.jl")
include("physics/longwave_radiation.jl")
include("physics/shortwave_radiation.jl")

# OUTPUT
include("output/output.jl") # defines Output
include("output/feedback.jl") # defines Feedback
include("output/pretty_printing.jl")
include("utility_functions.jl")
include("abstract_types.jl")

# LowerTriangularMatrices for spherical harmonics
export LowerTriangularMatrices
include("LowerTriangularMatrices/LowerTriangularMatrices.jl")
using .LowerTriangularMatrices

# RingGrids
export RingGrids
include("RingGrids/RingGrids.jl")
using .RingGrids

# SpeedyTransforms
export SpeedyTransforms
include("SpeedyTransforms/SpeedyTransforms.jl")
using .SpeedyTransforms

# INTERFACE
include("run_speedy.jl")
include("gpu.jl") # defines utility for GPU / KernelAbstractions
include("default_parameters.jl") # defines Parameters

# DYNAMICS
include("dynamics/constants.jl") # defines DynamicsConstants
include("dynamics/geometry.jl") # defines Geometry
include("dynamics/boundaries.jl") # defines Boundaries
include("dynamics/define_diffusion.jl") # defines HorizontalDiffusion
include("dynamics/define_implicit.jl") # defines ImplicitShallowWater, ImplicitPrimitiveEq
include("dynamics/parameter_structs.jl") # defines GenLogisticCoefs
include("dynamics/models.jl") # defines ModelSetups
include("dynamics/prognostic_variables.jl") # defines PrognosticVariables
include("dynamics/diagnostic_variables.jl") # defines DiagnosticVariables
include("dynamics/initial_conditions.jl")
include("dynamics/scaling.jl")
include("dynamics/geopotential.jl")
include("dynamics/tendencies_dynamics.jl")
include("dynamics/tendencies.jl")
include("dynamics/implicit.jl")
include("dynamics/diffusion.jl")
include("dynamics/time_integration.jl")

# PHYSICS
include("physics/constants.jl") # defines Parameterization Constants
include("physics/parameter_structs.jl") # defines MagnusCoefs, RadiationCoefs
include("physics/column_variables.jl") # defines ColumnVariables
include("physics/thermodynamics.jl")
include("physics/tendencies.jl")
include("physics/convection.jl")
include("physics/large_scale_condensation.jl")
include("physics/longwave_radiation.jl")
include("physics/shortwave_radiation.jl")
include("physics/boundary_layer.jl")

# OUTPUT
include("output/output.jl") # defines Output
include("output/feedback.jl") # defines Feedback
include("output/pretty_printing.jl")

# INTERFACE
include("run_speedy.jl")
end
3 changes: 3 additions & 0 deletions src/abstract_types.jl
Original file line number Diff line number Diff line change
Expand Up @@ -25,6 +25,9 @@ abstract type AbstractBoundaries{NF} end
# ATMOSPHERIC COLUMN FOR PARAMETERIZATIONS
abstract type AbstractColumnVariables{NF} end

# PARAMETERIZATIONS
abstract type BoundaryLayer end

# INPUT/OUTPUT
abstract type AbstractFeedback end
abstract type AbstractOutput end
Expand Down
3 changes: 3 additions & 0 deletions src/default_parameters.jl
Original file line number Diff line number Diff line change
Expand Up @@ -95,6 +95,9 @@ Base.@kwdef struct Parameters{Model<:ModelSetup} <: AbstractParameters{Model}
# Radiation
radiation_coefs::Coefficients = RadiationCoefs{NF}()

# BOUNDARY LAYER
boundary_layer::BoundaryLayer = LinearDrag{NF}()

# TIME STEPPING
startdate::DateTime = DateTime(2000,1,1) # time at which the integration starts
n_days::Float64 = 10 # number of days to integrate for
Expand Down
35 changes: 15 additions & 20 deletions src/dynamics/tendencies.jl
Original file line number Diff line number Diff line change
@@ -1,23 +1,24 @@
"""
get_tendencies!(diagn,model)
dynamics_tendencies!(diagn,model)
Calculate all tendencies for the BarotropicModel."""
function get_tendencies!( diagn::DiagnosticVariablesLayer,
model::BarotropicModel)
function dynamics_tendencies!( diagn::DiagnosticVariablesLayer,
model::BarotropicModel)

# only (absolute) vorticity advection for the barotropic model
vorticity_flux_divcurl!(diagn,model,curl=false) # = -∇⋅(u(ζ+f),v(ζ+f))
end

"""
get_tendencies!(diagn,surface,pres,time,model)
dynamics_tendencies!(diagn,surface,pres,time,model)
Calculate all tendencies for the ShallowWaterModel."""
function get_tendencies!( diagn::DiagnosticVariablesLayer,
surface::SurfaceVariables,
pres::LowerTriangularMatrix, # spectral pressure/η for geopotential
time::DateTime, # time to evaluate the tendencies at
model::ShallowWaterModel) # struct containing all constants
function dynamics_tendencies!( diagn::DiagnosticVariablesLayer,
surface::SurfaceVariables,
pres::LowerTriangularMatrix, # spectral pressure/η for geopotential
time::DateTime, # time to evaluate the tendencies at
model::ShallowWaterModel) # struct containing all constants

S,C = model.spectral_transform, model.constants

# for compatibility with other ModelSetups pressure pres = interface displacement η here
Expand All @@ -33,23 +34,17 @@ function get_tendencies!( diagn::DiagnosticVariablesLayer,
end

"""
get_tendencies!(diagn,surface,pres,time,model)
dynamics_tendencies!(diagn,surface,pres,time,model)
Calculate all tendencies for the primitive equation model (wet or dry)."""
function get_tendencies!( diagn::DiagnosticVariables,
progn::PrognosticVariables,
time::DateTime,
model::PrimitiveEquation,
lf::Int=2 # leapfrog index for tendencies
)
function dynamics_tendencies!( diagn::DiagnosticVariables,
progn::PrognosticVariables,
model::PrimitiveEquation,
lf::Int=2) # leapfrog index for tendencies

B, G, S = model.boundaries, model.geometry, model.spectral_transform
@unpack surface = diagn

# PARAMETERIZATIONS
# parameterization_tendencies!(diagn,time,model)

# DYNAMICS
# for semi-implicit corrections (α >= 0.5) linear gravity-wave related tendencies are
# evaluated at previous timestep i-1 (i.e. lf=1 leapfrog time step)
# nonlinear terms and parameterizations are always evaluated at lf
Expand Down
4 changes: 3 additions & 1 deletion src/dynamics/tendencies_dynamics.jl
Original file line number Diff line number Diff line change
Expand Up @@ -297,7 +297,9 @@ function _vertical_advection!( ξ_tend::Grid, # tendency of quantity
Δσₖ::NF # layer thickness on σ levels
) where {NF<:AbstractFloat,Grid<:AbstractGrid{NF}}
Δσₖ2⁻¹ = -1/2Δσₖ # precompute
@. ξ_tend = Δσₖ2⁻¹ * (σ_tend_below*(ξ_below - ξ) + σ_tend_above*- ξ_above))

# += as the tendencies already contain the parameterizations
@. ξ_tend += Δσₖ2⁻¹ * (σ_tend_below*(ξ_below - ξ) + σ_tend_above*- ξ_above))
end

function vordiv_tendencies!(diagn::DiagnosticVariablesLayer,
Expand Down
Loading

0 comments on commit 12fd32e

Please sign in to comment.