ECCOtour
Base.replace!
Base.replace!
Base.write
ECCOtour.allstats
ECCOtour.annualcontrib!
ECCOtour.apply_hemisphere_mask!
ECCOtour.apply_regional_mask!
ECCOtour.basin_mask
ECCOtour.calc_UV_conv3D!
ECCOtour.centerlon!
ECCOtour.columnscale!
ECCOtour.columnscale!
ECCOtour.depthlevels
ECCOtour.exch_UV_cs3D
ECCOtour.extract_sst34
ECCOtour.extract_timeseries
ECCOtour.extract_timeseries
ECCOtour.extract_timeseries
ECCOtour.factors4regularpoles
ECCOtour.faststats
ECCOtour.get_filtermatrix
ECCOtour.hanncoeffs
ECCOtour.hannfilter
ECCOtour.hannfilter
ECCOtour.hannsum
ECCOtour.hannsum!
ECCOtour.hannsum!
ECCOtour.inrectangle
ECCOtour.land2nan!
ECCOtour.land2zero!
ECCOtour.landmask
ECCOtour.latlonC
ECCOtour.latlonG
ECCOtour.matrixfilter
ECCOtour.matrixsaveinterpolation
ECCOtour.matrixspray
ECCOtour.mdsio2dict
ECCOtour.mdsio2sigma
ECCOtour.mdsio2sigma1
ECCOtour.mdsio2sigma2
ECCOtour.nancount
ECCOtour.ncwritefromtemplate
ECCOtour.netcdf2dict
ECCOtour.netcdf2dict
ECCOtour.netcdf2sigma1
ECCOtour.patchmean
ECCOtour.position_label
ECCOtour.prereginterp
ECCOtour.pressurelevels
ECCOtour.read_netcdf
ECCOtour.readarea
ECCOtour.reginterp
ECCOtour.regional_mask
ECCOtour.regularpoles
ECCOtour.remove_climatology
ECCOtour.remove_seasonal
ECCOtour.remove_seasonal
ECCOtour.rotate_uv
ECCOtour.searchdir
ECCOtour.seasonal_matrices
ECCOtour.sigma
ECCOtour.sigma2grid
ECCOtour.smooth
ECCOtour.time_label
ECCOtour.timestamp_monthly_v4r4
ECCOtour.trend_matrices
ECCOtour.trend_theta
ECCOtour.trend_theta!
ECCOtour.vars2sigma
ECCOtour.velocity2center
ECCOtour.velocity2center!
ECCOtour.wrapdist
ECCOtour.write_vars
ECCOtour.writeregularpoles
ECCOtour.writeregularpoles
ECCOtour.writeregularpoles
ECCOtour.θbudg2sigma
ECCOtour.θbudg2sigma2
IsopycnalSurfaces.vars2sigma1
Statistics.mean
Statistics.mean
Statistics.mean
Statistics.std
Base.replace!
— Methodfunction replace!(f,a::MeshArrays.gcmarray{T,N,Array{T,2}}}
Arguments
f::function
= replace a with f(a)a::MeshArrays.gcmarray{T,N,Array{T,2}}
: gcmarray with variable type and time-dimension
Base.replace!
— Methodfunction replace!(a::MeshArrays.gcmarray{T,N,Array{T,2}}},b::Pair) where T
Arguments
a::MeshArrays.gcmarray{T,N,Array{T,2}}
: gcmarray with variable type and time-dimensionb::Pair
: replacea
elements of pair 1 with pair 2
Base.write
— Methodfunction write(vars::Dict{String,Array{Float32,3}}, params::RegularpolesParameters, γ::MeshArrays.gcmgrid, pathout, filesuffix, filelog, gridatts)
ECCOtour.allstats
— Methodfunction allstats(x)
Compute all statistics of gcmgrid type using function calls, same as faststats
Warning: assumes land is zero
Input
x::MeshArrays.gcmarray{T,N,Matrix{T}}
: input of gcmarray type
Output
xbar::T
: mean value after discarding dryvalxmax::T
: maximum valuexmin::T
: minimum valueσx::T
: standard deviationabsxbar::T
: mean of absolute value
ECCOtour.annualcontrib!
— Methodfunction matrixfilter_annualcontrib!(θout,F,θin)
Arguments
θout
: filtered valuesF
: filter function as a matrixθin
: input values at all times t (e.g., θ(t))
ECCOtour.apply_hemisphere_mask!
— Methodfunction apply_hemisphere_mask!(mask,hemisphere,γ)
overlay a hemispheric mask on a previous `mask`
in-place function
hemisphere options are `:north`,`:south`, and `:both`
Note: both has not been tested with this version
ECCOtour.apply_regional_mask!
— Methodfunction regional_mask(field,mask)
Arguments
field
: input field that is mutated/modifiedmask
: multiplicative factor, same spatial grid as field
ECCOtour.basin_mask
— Methodfunction basin_mask(basin_name,γ;hemisphere=nothing,southlat=nothing,northlat=nothing,Lsmooth=nothing)
Make a mask based on Gael Forget's definitions of ocean basins and sub-basins.
Note: This mask contains float values. It could be more efficient with a BitArray.
Arguments
basin_name::Union{String,Vector{String}}
: string options are Arctic, Atlantic, Baffin Bay, Barents Sea, Bering Sea,
East China Sea, GIN Seas, Gulf, Gulf of Mexico, Hudson Bay, indian, Japan Sea, Java Sea, Mediterranean Sea, North Sea, Okhotsk Sea, Pacific, Red Sea, South China Sea, Timor Sea. -hemisphere::Symbol
: optional argument with values :north
, :south
, and :both
-'southlat::Number': optional argument specifying southerly latitude of mask -'northlat::Number': optional argument specifying northerly latitude of mask -Lsmooth::Number
: smoothing lengthscale in grid points -'southlat::Number': optional argument specifying southerly latitude of mask -'northlat::Number': optional argument specifying northerly latitude of mask
Output
- 'mask': space and time field of surface forcing, value of zero inside
designated lat/lon rectangle and fading to 1 outside sponge zone on each edge. This is because this field ends up being SUBTRACTED from the total forcing
ECCOtour.calc_UV_conv3D!
— Methodfunction calc_UV_conv3D!(uFLD::MeshArrays.gcmarray{T, 2, Matrix{T}},
vFLD::MeshArrays.gcmarray{T, 2, Matrix{T}}, CONV::MeshArrays.gcmarray{T, 2, Matrix{T}}) where T<:Real
tmpU, tmpV = exch_UV_cs3D(uFLD,vFLD)
By Anthony Meza
ECCOtour.centerlon!
— Methodfunction centerlon!
Make a longitudinal coordinate system
centered at lonmid.
Makes handling wraparound easy later on.
ECCOtour.columnscale!
— Methodfunction columnscale!(product,M,flux)
product::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
M::Array{Float32,1}
flux::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
ECCOtour.columnscale!
— Methodfunction columnscale!(product,M,flux)
product::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
M::Array{Float32,1}
flux::MeshArrays.gcmarray{Float32,1,Array{Float32,2}}
ECCOtour.depthlevels
— Methodfunction depthlevels(γ)
Depths of ECCO v4r4 grid (positive vals)
ECCOtour.exch_UV_cs3D
— Methodfunction exch_UV_cs3D(fldU::MeshArrays.gcmarray{T, 2, Matrix{T}},
fldV::MeshArrays.gcmarray{T, 2, Matrix{T}}) where T<:Real
By Anthony Meza
ECCOtour.extract_sst34
— Methodfunction extract_sst34
extract by reading multiple files
ECCOtour.extract_timeseries
— Methodextract_timeseries(froot,years,γ,xval,yval,fval)
Arguments
froot::String
: filename rootyears::StepRange
: iterator for multiple filesγ::gcmarray
: GCM grid from MeshArrays.jlxval::Integer
: x grid indexyval::Integer
: y grid indexfval::Integer
: face index
Output
tseries
: timeseriesnseries
: nseries = number of timeseries elements in each year
```
ECCOtour.extract_timeseries
— Methodextract_timeseries(flux,xval,yval,fval)
Arguments
flux::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
: input array of arraysyears::StepRange
: iterator for multiple filesxval::Integer
: x grid indexyval::Integer
: y grid indexfval::Integer
: face index
Output
series
: timeseries
ECCOtour.extract_timeseries
— Methodextract_timeseries(flux,xval,yval,fval)
Arguments
flux::MeshArrays.gcmarray{Float64,2,Array{Float64,2}}
: input array of arraysyears::StepRange
: iterator for multiple filesxval::Integer
: x grid indexyval::Integer
: y grid indexfval::Integer
: face index
Output
series
: timeseries
ECCOtour.factors4regularpoles
— Methodfunction factors4regularpoles(γ)
Get interpolation factors for regularpoles grid in one place
ECCOtour.faststats
— Methodfunction faststats(x)
Compute fast statistics of gcmgrid type using function calls, eliminate redundancy
Input
x::MeshArrays.gcmarray{T,1,Matrix{T}}
: input of gcmarray type
Output
xbar::T
: mean value after discarding dryvalxmax::T
: maximum valuexmin::T
: minimum valueσx::T
: standard deviationabsxbar::T
: mean of absolute value
ECCOtour.get_filtermatrix
— Methodget_filtermatrix(tin, tout; ival = 3)
This function generates a filter matrix based on the input time arrays tin
and tout
.
Arguments
tin
: Input time array.tout
: Output time array.ival
: : Optional argument with a default value of 3. It is used to determine the starting index for the filter.
e.g ival = 3 # first two 14-day values of Eout2in are set to zero ival = 1 # all days of Eout2in are allowed to be non-zero
Output
Eout2in
: matrix that maps from tout to tinFin2out
: matrix that maps from tin to tout
ECCOtour.hanncoeffs
— Methodhanncoeffs(tin,tout,L)
Arguments
tin
: input dataset time (vector valued)tout
: time where filtered output will be createdL
: filter length
Output
h
: normalized coefficients of Hann(ing) filter
ECCOtour.hannfilter
— Methodfunction hannfilter(θin,tin,tout,L,γ)
Arguments
θin::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
tin::StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}
tout::StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}
L::Float64
γ::gcmgrid
Output
θout
ECCOtour.hannfilter
— Methodfunction hannfilter(θin,tin,tout,L)
Arguments
θin::Array{Float32,1}
: inputtin
tout
L
Output
θout
: output
ECCOtour.hannsum!
— Methodhannsum!(θout,θin,tin,tout,L)
use this when a timeseries is available
Arguments
θout
: output value at time th (e.g., θh(th))θin
: input value at all times t (e.g., θ(t))tin
: input dataset time (vector valued)tout
: time where filtered output will be createdL
: filter length
ECCOtour.hannsum!
— Methodfunction hannsum!(θout,θin,tin,tout,toutindex,L::Float64)
Call this version for mesharrays
Arguments
θout::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
: output value at time tout (e.g., θout(tout))θin::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
: input value at all times t (e.g., θ(t))tin::StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}
: input dataset time (vector valued)tout::StepRangeLen{Float64,Base.TwicePrecision{Float64},Base.TwicePrecision{Float64}}
: time where filtered output will be createdtoutindex::Integer
: output times?L::Float64
: filter length
ECCOtour.hannsum
— Methodhannsum(θin,tin,tout,L)
use this when a timeseries is available
Arguments
θ
: input value at all times t (e.g., θ(t))tout
: time where filtered output will be createdtin
: input dataset time (vector valued)L
: filter length
Output
θout
: output value at time th (e.g., θh(th))
ECCOtour.inrectangle
— Methodfunction inrectangle(lonin,latin,lons,lats)
find all gridcells within Cartesian rectangle
Arguments
lons
: tuple of longitude rangelats
: tuple of latitude range
Output
rectangle
: boolean of type gcmarray
ECCOtour.land2nan!
— Methodfunction land2nan!(msk,γ)
Replace surface land points with NaN
Deprecate this function: slower than MeshArrays.mask
ECCOtour.land2zero!
— Methodfunction land2zero!(msk,γ)
see tests for different masking from MeshArrays
that could replace this function
ECCOtour.landmask
— Methodfunction landmask(γ;level=1)
The land mask at a given depth
Default is the surface
ECCOtour.latlonC
— Methodfunction latlonC(γ)
Latitude-longitude of ECCOv4r4 "C" (tracer) grid
ECCOtour.latlonG
— Methodfunction latlongG(γ)
Latitude-longitude of ECCOv4r4 "G" (velocity) grid
ECCOtour.matrixfilter
— Methodmatrixfilter(F,froot,years,γ)
writing it in a funny way to save computation
issue with timeseries being read in different files
Arguments
F
: filter in matrix formfroot
: filename rootyears
: iterator for multiple filesγ
: GCM grid (meshArray type)
ECCOtour.matrixsaveinterpolation
— Methodfunction matrixsaveinterpolation(E,savefield,frootin,frootout,years,γ)
writing it in a funny way to save computation
issue with timeseries being read in different files.
Interpolates "savefield" using the E matrix onto a monthly grid.
Then saves the interpolated field "savefield".
By Anthony Meza
Arguments
E
: interpolating spray operator in matrix formsavefield
: field to be saved to multiple filesfrootin
: filename root of input fieldfrootout
: filename root of output fieldyears
: iterator for multiple filesγ
: GCM grid (meshArray type)
ECCOtour.matrixspray
— Methodmatrixspray(F,rmfield,frootin,frootout,years,γ)
writing it in a funny way to save computation
issue with timeseries being read in different files
Arguments
F
: spray operator in matrix formrmfield
: field to be removed from filesfrootin
: filename root of input fieldfrootout
: filename root of output fieldyears
: iterator for multiple filesγ
: GCM grid (meshArray type)
ECCOtour.mdsio2dict
— Methodfunction mdsio2dict(pathin,filein,γ)
Read native (i.e., mdsio) MITgcm files and return a Dictionary.
ECCOtour.mdsio2sigma
— Methodfunction mdsio2sigma(pathin::String,pathout::String,fileroots::Vector{String},γ,pstdz,siggrid;splorder=3,linearinterp=false,eos="TEOS10")
Take variables in a filelist, read, map to sigma1, write to file.
Arguments
pathin::String
: path of inputpathout::String
: path of outputfileroots::String
: beginning of file names in pathinγ
: grid description needed for preallocationsplorder::Integer
: optional keyword argument, default = 3, order of splinelinearinterp::Logical
: optional keyword argument, do linear interpolation?, default = falseeos::String
: optional key argument for equation of state, default = "JMD95"
ECCOtour.mdsio2sigma1
— Methodfunction mdsio2sigma1(pathin::String,pathout::String,fileroots::Vector{String},γ,pstdz,sig1grid;splorder=3,linearinterp=false,eos="JMD95")
Take variables in a filelist, read, map to sigma1, write to file.
Arguments
pathin::String
: path of inputpathout::String
: path of outputfileroots::String
: beginning of file names in pathinγ
: grid description needed for preallocationsplorder::Integer
: optional keyword argument, default = 3, order of splinelinearinterp::Logical
: optional keyword argument, do linear interpolation?, default = falseeos::String
: optional key argument for equation of state, default = "JMD95"
ECCOtour.mdsio2sigma2
— Methodmdsio2sigma2(pathin::String,pathout::String,fileroots::Vector{String},
γ::gcmgrid,pstdz::Vector{Float64},siggrid::Vector{Float64};
splorder=3,linearinterp=false,eos="JMD95") = mdsio2sigma(pathin,pathout,fileroots,γ,
pstdz,siggrid, 2000.0, "sigma2";
splorder=splorder,linearinterp=linearinterp,eos=eos)
By Anthony Meza
ECCOtour.nancount
— Methodfunction nancount
Count number of NaN's in gcmgrid type
Input
field::MeshArrays.gcmarray{T,N,Matrix{T}}
: input of gcmarray type, can have multiple records
Output
nannum::Int
: number of NaNs
ECCOtour.ncwritefromtemplate
— Methodfunction ncwritefromtemplate This version writes NetCDF output on the regularpoles grid
ECCOtour.netcdf2dict
— Methodnetcdf2dict array output (?)
ECCOtour.netcdf2dict
— Methodnetcdf2dict gcmarray output
ECCOtour.netcdf2sigma1
— Methodfunction netcdf2sigma1
Take variables in a NetCDF filelist, read, map to sigma1, write to file.
ECCOtour.patchmean
— Method function patchmean(x,area,ϕ,λ,ispatch,iswet)
get weight for rectangle region.
Arguments
x
: variable of interestarea
: weightingϕ
: latλ
: lonispatch
: true in patch of interestiswet
= function that is true if in ocean
Output
xbar
: weighted filtered average of x
ECCOtour.position_label
— Methodfunction position_label(lon,lat)
Arguments
lon
: longitudelat
: latitude
Output
lbl
: LaTex string used for figure labels
ECCOtour.prereginterp
— Methodfunction prereginterp(latgrid,longrid,γ)
prepare for regular interpolation
regular = onto rectangular 2d map
Arguments
latgrid
: 1d array of latitudelongrid
: 1d array of longitudeγ
: GCM grid
Output
f,i,j,w
: interpolation factors
ECCOtour.pressurelevels
— Methodfunction pressurelevels(z)
Standard pressures of ECCO v4r4 grid
ECCOtour.read_netcdf
— Methodread_netcdf(fileName,fldName,mygrid)
Read model output from NetCDF files that are global and convert to MeshArray instance. ```
ECCOtour.readarea
— Methodfunction readarea(γ)
area of ECCO v4r4 grid
ECCOtour.reginterp
— Methodfunction reginterp(fldin,nx,ny,f,i,j,w)
regular interpolation with precomputed factors
regular = onto rectangular 2d map
Arguments
fldin
: gcmarray field of inputnx,ny
: x,y dimension of output regular fieldλ=(f,i,j,w)
: precomputed interpolation factors
Output
fldout
: interpolated regular 2d field
ECCOtour.regional_mask
— Methodfunction regional_mask(latpt,lonpt,latrect,lonrect,dlat,dlon)
Arguments
latpt
: latitude gridlonpt
: longitude gridlatrect
: tuple of latitude range of central rectanglelonrect
: tuple of longitude range of central rectangledlon
: width in degrees longitude of East/West fade zones- 'dlat': width in degrees latitude of North/South fade zones
Output
- 'mask': space and time field of surface forcing, value of zero inside
designated lat/lon rectangle and fading to 1 outside sponge zone on each edge. This is because this field ends up being SUBTRACTED from the total forcing
ECCOtour.regularpoles
— Methodfunction mdsio2regularpoles(pathin,filein,γ,nx,ny,nyarc,λarc,nyantarc,λantarc)
ECCOtour.remove_climatology
— Methodfunction remove_climatology(x,xbar)
`x` = long monthly timeseries, starting in Jan
`xbar` = 12-month climatology, starting in Jan
remove climatology
ECCOtour.remove_seasonal
— Methodfunction remove_seasonal(x,Ecycle,Fcycle,γ) natively
`x` = long monthly timeseries, starting in Jan
remove seasonal cycle of this timeseries
ECCOtour.remove_seasonal
— Methodfunction remove_seasonal(x,Ecycle,Fcycle,γ)
remove seasonal cycle of this timeseries on the native gcmgrid
Arguments
x::gcmarray{T,N,Array{T,2}}
= gcmarray of a long timeseries where time is dimension NEcycle::Matrix
: matrix that operates on seasonal cycle parameters and return a seasonal cycle timeseriesFcycle::Matrix
: precomputed pseudo-inverse ofEcycle
γ::gcmgrid
: MITgcm grid
Output
x′::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
: potential density referenced to p₀ minus 1000 kg/m³
ECCOtour.rotate_uv
— Methodfunction rotate_uv(uvel,vvel,G) From Gael Forget, JuliaClimateNotebooks/Transport 3. Convert to Eastward/Northward
transport 4. Display Subdomain Arrays (optional)
ECCOtour.searchdir
— Methodfunction searchdir(path,key)
a useful one-line function
Arguments
path
: directory to search for filekey
: expression to search for in path
ECCOtour.seasonal_matrices
— Functionfunction seasonal_matrices(fcycle,t,overtones=1)
Arguments
fcycle
: frequency of seasonal cyclet
: timeovertones=1
: optional argument for number of overtones
Output
E
: matrix that solves E*parameters= seasonal cycleF=E†
: generalized inverse ofE
ECCOtour.sigma
— Methodfunction sigma(θ,S,pz,p₀)
sigma values from SeaWaterDensity for gcmarrays
Arguments
θ::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
: potential temperatureS::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
: practical salinitypz::Array{Float64,1}
: vertical profile of standard pressuresp₀::Int64
: reference pressure for sigma value
Output
σ::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
: potential density referenced to p₀ minus 1000 kg/m³
ECCOtour.sigma2grid
— Methodfunction sigma2grid() Standard chosen from the time-averaged Pacific Ocean Density Configuration. From this analysis, density extrema were found to be between σ2 ∈ (29, 37) kg per m^3
Arguments
z
: value
Output
σ₂ grid
: list (vector) of σ₁ values
ECCOtour.smooth
— Methodfunction smooth(msk::MeshArrays.gcmarray,lengthscale)
Smooth a gcmarray with a lengthscale of `X` points
Based off Gael Forget, MeshArrays.jl
ECCOtour.time_label
— Methodfunction time_label(index)
Time label useful in plots. For ECCOv4r4.
Input
index::Integer
: index number for monthly average field, number of months since Jan. 1992
Output
tlbl
: label with month name (abbrev. with 3 letters) and calendar year CE
ECCOtour.timestamp_monthly_v4r4
— Methodfunction timestamp_monthly_v4r4(t)
print year and month given time index
assuming using ECCOv4r4
ECCOtour.trend_matrices
— Methodfunction trend_matrices(t)
Arguments
t
: time
Output
E
: matrix that solves E*parameters= timeseriesF=E†
: generalized inverse ofE
ECCOtour.trend_theta!
— Methodfunction trend_theta!(β,diagpathexp,path_out,tecco,γ,F)
get linear trend of potential temperature
from monthly-average ECCOv4r4 gcmarray fields
Arguments
β
: gcmarray field of trendsdiagpathexp
: where to find filestecco
: time stamps of monthly outputγ
: GCM grid informationF
: linear estimator of trend
Output
β
: updated trends
ECCOtour.trend_theta
— Methodfunction trend_theta(diagpathexp,path_out,tecco,γ,F)
get linear trend of potential temperature
from monthly-average ECCOv4r4 gcmarray fields
Arguments
diagpathexp
: where to find filestecco
: time stamps of monthly outputγ
: GCM grid informationF
: linear estimator of trend
Output
β
: updated trends
ECCOtour.vars2sigma
— Methodfunction vars2sigma(vars,p,siggrid,γ,spline_order)
map variables onto sigma1 surfaces for gcmarrays
Arguments
vars::Dict{String,MeshArrays.gcmarray{T,N,Array{T,2}}}
: dict of gcmarraysp::Array{Float64,1}
: vertical profile of standard pressuressiggrid
: σ₁ surface valuesγ
: grid description needed for preallocationsplorder
: 1-5, order of splinelinearinterp
: optional logical
Output
varsσ::Dict{String,MeshArrays.gcmarray{T,N,Array{T,2}}}
: dict of gcmarrays of variables on sigma1 surfaces
ECCOtour.velocity2center!
— Methodvelocity2center!(uC,vC,u,v,G)
From Gael Forget, JuliaClimateNotebooks/Transport
#1. Convert to `Sv` units and mask out land
2. Interpolate `x/y` transport to grid cell center
ECCOtour.velocity2center
— Methodvelocity2center(u,v,G)
From Gael Forget, JuliaClimateNotebooks/Transport
#1. Convert to `Sv` units and mask out land
2. Interpolate `x/y` transport to grid cell center
ECCOtour.wrapdist
— Methodfunction wrapdist
Longitudinal distance in degrees
Passing the date line may be the shortest distance.
ECCOtour.write_vars
— Methodfunction write_vars This version writes mdsio output on the native grid Missing: write the accompanying meta file
ECCOtour.writeregularpoles
— Methodfunction writeregularpoles(vars,γ,pathout,filesuffix,filelog,λC,lonatts,ϕC,latatts,z,depthatts)
ECCOtour.writeregularpoles
— Methodfunction writeregularpoles(vars::Dict{String,Array{Float32,2}},γ,pathout,filesuffix,filelog,λC,lonatts,ϕC,latatts,z,depthatts)
ECCOtour.writeregularpoles
— Methodfunction writeregularpoles(vars,pathout,filesuffix,filelog,λC,lonatts,ϕC,latatts)
ECCOtour.θbudg2sigma
— Methodfunction θbudg2sigma(pathin::String,pathout::String,fileroots::Vector{String},γ,pstdz,siggrid;splorder=3,linearinterp=false,eos="TEOS10")
Take variables in a filelist, read, map to sigma1, write to file.
Arguments
inv_cell_volumes::MeshArrays.gcmarray{T, 2, Matrix{T}}
: Inverted cell volumes, land should be 0, wet pts should be 1/Vpathin::String
: path of inputpathout::String
: path of outputfileroots::String
: beginning of file names in pathinγ
: grid description needed for preallocationsplorder::Integer
: optional keyword argument, default = 3, order of splinelinearinterp::Logical
: optional keyword argument, do linear interpolation?, default = falseeos::String
: optional key argument for equation of state, default = "JMD95"
ECCOtour.θbudg2sigma2
— Methodθbudg2sigma2(invRAC::MeshArrays.gcmarray{T, 1, Matrix{T}}, invRAU::MeshArrays.gcmarray{T, 2, Matrix{T}}, invRAV::MeshArrays.gcmarray{T, 2, Matrix{T}}, pathin::String,pathout::String,fileroots::Vector{String}, γ::gcmgrid,pstdz::Vector{Float64},siggrid::Vector{Float64}; splorder=3,linearinterp=false,eos="JMD95") where T<:Real = θbudg2sigma(invRAC, invRAU, invRAV, pathin,pathout,fileroots,γ, pstdz,siggrid, 2000.0, "sigma2"; splorder=splorder,linearinterp=linearinterp,eos=eos)
By Anthony Meza
IsopycnalSurfaces.vars2sigma1
— Methodfunction vars2sigma1(vars,p,sig1grid,γ,spline_order)
map variables onto sigma1 surfaces for gcmarrays
Arguments
vars::Dict{String,MeshArrays.gcmarray{T,N,Array{T,2}}}
: dict of gcmarraysp::Array{Float64,1}
: vertical profile of standard pressuressig1grid
: σ₁ surface valuesγ
: grid description needed for preallocationsplorder
: 1-5, order of splinelinearinterp
: optional logical
Output
varsσ::Dict{String,MeshArrays.gcmarray{T,N,Array{T,2}}}
: dict of gcmarrays of variables on sigma1 surfaces
Statistics.mean
— Methodfunction mean
Compute area-weighted mean of gcmgrid type using filtered with function isgood
Area weighting = area
Eliminate all values = isgood(x) -> true
Input
x::MeshArrays.gcmarray{Float32,1,Array{Float32,2}}
: input of gcmarray typeweight::MeshArrays.gcmarray{Float32,1,Array{Float32,2}}
: weighting variable of gcmarray typeisgood::Function
: returns true is a value to be used in the mean
Output
xbar::Float32
: mean value (weighted and filtered)
Statistics.mean
— Methodfunction mean
Compute area-weighted mean of gcmgrid type using function calls, eliminate redundancy
Area weighting = area
Eliminate all values = dryval
Input
x::MeshArrays.gcmarray{Float32,1,Array{Float32,2}}
: input of gcmarray typeweight::MeshArrays.gcmarray{Float32,1,Array{Float32,2}}
: weighting variable of gcmarray typedryval::Float32
: land value (doesn't work for NaN32)
Output
xbar::Float32
: mean value (unweighted)
Statistics.mean
— Methodfunction mean
Compute mean of gcmgrid type using function calls, eliminate redundancy
Eliminate all values = dryval
Input
x::MeshArrays.gcmarray{T,N,Array{T,2}}
: input of gcmarray typedryval::T
: land value (doesn't work for NaN32)
Output
xbar::T
: mean value (unweighted)
Statistics.std
— Methodfunction std
Compute standard deviation of gcmgrid type using function calls, eliminate redundancy
Eliminate all values = dryval
Input
x::MeshArrays.gcmarray{T,N,Array{T,2}}
: input of gcmarray typexbar::T
: mean valuedryval::T
: land value
Output
σx::T
: standard deviation