ECCOtour

Base.replace!Method

function 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
source
Base.replace!Method

function 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-dimension
  • b::Pair: replace a elements of pair 1 with pair 2
source
Base.writeMethod

function write(vars::Dict{String,Array{Float32,3}}, params::RegularpolesParameters, γ::MeshArrays.gcmgrid, pathout, filesuffix, filelog, gridatts)

source
ECCOtour.allstatsMethod
function 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 dryval
  • xmax::T: maximum value
  • xmin::T: minimum value
  • σx::T: standard deviation
  • absxbar::T: mean of absolute value
source
ECCOtour.annualcontrib!Method
function matrixfilter_annualcontrib!(θout,F,θin)

Arguments

  • θout: filtered values
  • F: filter function as a matrix
  • θin: input values at all times t (e.g., θ(t))
source
ECCOtour.apply_hemisphere_mask!Method
function 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
source
ECCOtour.apply_regional_mask!Method
function regional_mask(field,mask)

Arguments

  • field: input field that is mutated/modified
  • mask: multiplicative factor, same spatial grid as field
source
ECCOtour.basin_maskMethod
function 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

source
ECCOtour.calc_UV_conv3D!Method
function 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
source
ECCOtour.centerlon!Method
function centerlon!
Make a longitudinal coordinate system
centered at lonmid.
Makes handling wraparound easy later on.
source
ECCOtour.columnscale!Method
function columnscale!(product,M,flux)
  • product::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
  • M::Array{Float32,1}
  • flux::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
source
ECCOtour.columnscale!Method
function columnscale!(product,M,flux)
  • product::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}
  • M::Array{Float32,1}
  • flux::MeshArrays.gcmarray{Float32,1,Array{Float32,2}}
source
ECCOtour.exch_UV_cs3DMethod
function exch_UV_cs3D(fldU::MeshArrays.gcmarray{T, 2, Matrix{T}},
    fldV::MeshArrays.gcmarray{T, 2, Matrix{T}}) where T<:Real

By Anthony Meza
source
ECCOtour.extract_timeseriesMethod
extract_timeseries(froot,years,γ,xval,yval,fval)

Arguments

  • froot::String: filename root
  • years::StepRange: iterator for multiple files
  • γ::gcmarray: GCM grid from MeshArrays.jl
  • xval::Integer: x grid index
  • yval::Integer: y grid index
  • fval::Integer: face index

Output

  • tseries: timeseries
  • nseries: nseries = number of timeseries elements in each year

```

source
ECCOtour.extract_timeseriesMethod
extract_timeseries(flux,xval,yval,fval)

Arguments

  • flux::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}: input array of arrays
  • years::StepRange: iterator for multiple files
  • xval::Integer: x grid index
  • yval::Integer: y grid index
  • fval::Integer: face index

Output

  • series: timeseries
source
ECCOtour.extract_timeseriesMethod
extract_timeseries(flux,xval,yval,fval)

Arguments

  • flux::MeshArrays.gcmarray{Float64,2,Array{Float64,2}}: input array of arrays
  • years::StepRange: iterator for multiple files
  • xval::Integer: x grid index
  • yval::Integer: y grid index
  • fval::Integer: face index

Output

  • series: timeseries
source
ECCOtour.faststatsMethod
function 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 dryval
  • xmax::T: maximum value
  • xmin::T: minimum value
  • σx::T: standard deviation
  • absxbar::T: mean of absolute value
source
ECCOtour.get_filtermatrixMethod
get_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 tin
  • Fin2out: matrix that maps from tin to tout
source
ECCOtour.hanncoeffsMethod
hanncoeffs(tin,tout,L)

Arguments

  • tin: input dataset time (vector valued)
  • tout: time where filtered output will be created
  • L: filter length

Output

  • h: normalized coefficients of Hann(ing) filter
source
ECCOtour.hannfilterMethod
function 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
source
ECCOtour.hannfilterMethod
function hannfilter(θin,tin,tout,L)

Arguments

  • θin::Array{Float32,1}: input
  • tin
  • tout
  • L

Output

  • θout: output
source
ECCOtour.hannsum!Method
hannsum!(θ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 created
  • L: filter length
source
ECCOtour.hannsum!Method
function 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 created
  • toutindex::Integer: output times?
  • L::Float64: filter length
source
ECCOtour.hannsumMethod
hannsum(θ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 created
  • tin: input dataset time (vector valued)
  • L: filter length

Output

  • θout: output value at time th (e.g., θh(th))
source
ECCOtour.inrectangleMethod
function inrectangle(lonin,latin,lons,lats)
find all gridcells within Cartesian rectangle

Arguments

  • lons: tuple of longitude range
  • lats: tuple of latitude range

Output

  • rectangle: boolean of type gcmarray
source
ECCOtour.land2nan!Method
function land2nan!(msk,γ)

Replace surface land points with NaN

Deprecate this function: slower than MeshArrays.mask

source
ECCOtour.land2zero!Method
function land2zero!(msk,γ)

see tests for different masking from MeshArrays

that could replace this function

source
ECCOtour.landmaskMethod
function landmask(γ;level=1)

 The land mask at a given depth
Default is the surface
source
ECCOtour.matrixfilterMethod
matrixfilter(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 form
  • froot: filename root
  • years: iterator for multiple files
  • γ: GCM grid (meshArray type)
source
ECCOtour.matrixsaveinterpolationMethod
function 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 form
  • savefield: field to be saved to multiple files
  • frootin: filename root of input field
  • frootout: filename root of output field
  • years: iterator for multiple files
  • γ: GCM grid (meshArray type)
source
ECCOtour.matrixsprayMethod
matrixspray(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 form
  • rmfield: field to be removed from files
  • frootin: filename root of input field
  • frootout: filename root of output field
  • years: iterator for multiple files
  • γ: GCM grid (meshArray type)
source
ECCOtour.mdsio2dictMethod

function mdsio2dict(pathin,filein,γ)

Read native (i.e., mdsio) MITgcm files and return a Dictionary.

source
ECCOtour.mdsio2sigmaMethod
function 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 input
  • pathout::String: path of output
  • fileroots::String: beginning of file names in pathin
  • γ: grid description needed for preallocation
  • splorder::Integer: optional keyword argument, default = 3, order of spline
  • linearinterp::Logical: optional keyword argument, do linear interpolation?, default = false
  • eos::String: optional key argument for equation of state, default = "JMD95"
source
ECCOtour.mdsio2sigma1Method
function 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 input
  • pathout::String: path of output
  • fileroots::String: beginning of file names in pathin
  • γ: grid description needed for preallocation
  • splorder::Integer: optional keyword argument, default = 3, order of spline
  • linearinterp::Logical: optional keyword argument, do linear interpolation?, default = false
  • eos::String: optional key argument for equation of state, default = "JMD95"
source
ECCOtour.mdsio2sigma2Method
mdsio2sigma2(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
source
ECCOtour.nancountMethod
function 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
source
ECCOtour.patchmeanMethod
 function patchmean(x,area,ϕ,λ,ispatch,iswet)
 get weight for rectangle region.

Arguments

  • x: variable of interest
  • area: weighting
  • ϕ: lat
  • λ: lon
  • ispatch: true in patch of interest
  • iswet = function that is true if in ocean

Output

  • xbar: weighted filtered average of x
source
ECCOtour.position_labelMethod
function position_label(lon,lat)

Arguments

  • lon: longitude
  • lat: latitude

Output

  • lbl: LaTex string used for figure labels
source
ECCOtour.prereginterpMethod
function prereginterp(latgrid,longrid,γ)
prepare for regular interpolation
regular = onto rectangular 2d map

Arguments

  • latgrid: 1d array of latitude
  • longrid: 1d array of longitude
  • γ: GCM grid

Output

  • f,i,j,w: interpolation factors
source
ECCOtour.read_netcdfMethod
read_netcdf(fileName,fldName,mygrid)

Read model output from NetCDF files that are global and convert to MeshArray instance. ```

source
ECCOtour.reginterpMethod
function reginterp(fldin,nx,ny,f,i,j,w)
regular interpolation with precomputed factors
regular = onto rectangular 2d map

Arguments

  • fldin: gcmarray field of input
  • nx,ny: x,y dimension of output regular field
  • λ=(f,i,j,w): precomputed interpolation factors

Output

  • fldout: interpolated regular 2d field
source
ECCOtour.regional_maskMethod
function regional_mask(latpt,lonpt,latrect,lonrect,dlat,dlon)

Arguments

  • latpt: latitude grid
  • lonpt: longitude grid
  • latrect: tuple of latitude range of central rectangle
  • lonrect: tuple of longitude range of central rectangle
  • dlon: 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

source
ECCOtour.remove_climatologyMethod
function remove_climatology(x,xbar)
`x` = long monthly timeseries, starting in Jan
`xbar` = 12-month climatology, starting in Jan
remove climatology
source
ECCOtour.remove_seasonalMethod
function remove_seasonal(x,Ecycle,Fcycle,γ) natively
`x` = long monthly timeseries, starting in Jan

remove seasonal cycle of this timeseries
source
ECCOtour.remove_seasonalMethod
function 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 N
  • Ecycle::Matrix: matrix that operates on seasonal cycle parameters and return a seasonal cycle timeseries
  • Fcycle::Matrix: precomputed pseudo-inverse of Ecycle
  • γ::gcmgrid: MITgcm grid

Output

  • x′::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}: potential density referenced to p₀ minus 1000 kg/m³
source
ECCOtour.rotate_uvMethod

function rotate_uv(uvel,vvel,G) From Gael Forget, JuliaClimateNotebooks/Transport 3. Convert to Eastward/Northward transport 4. Display Subdomain Arrays (optional)

source
ECCOtour.searchdirMethod
function searchdir(path,key)
a useful one-line function

Arguments

  • path: directory to search for file
  • key: expression to search for in path
source
ECCOtour.seasonal_matricesFunction
function seasonal_matrices(fcycle,t,overtones=1)

Arguments

  • fcycle: frequency of seasonal cycle
  • t: time
  • overtones=1: optional argument for number of overtones

Output

  • E: matrix that solves E*parameters= seasonal cycle
  • F=E†: generalized inverse of E
source
ECCOtour.sigmaMethod
function sigma(θ,S,pz,p₀)
sigma values from SeaWaterDensity for gcmarrays

Arguments

  • θ::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}: potential temperature
  • S::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}: practical salinity
  • pz::Array{Float64,1} : vertical profile of standard pressures
  • p₀::Int64 : reference pressure for sigma value

Output

  • σ::MeshArrays.gcmarray{Float32,2,Array{Float32,2}}: potential density referenced to p₀ minus 1000 kg/m³
source
ECCOtour.sigma2gridMethod

function 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
source
ECCOtour.smoothMethod
function smooth(msk::MeshArrays.gcmarray,lengthscale)

Smooth a gcmarray with a lengthscale of `X` points

Based off Gael Forget, MeshArrays.jl
source
ECCOtour.time_labelMethod
function 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
source
ECCOtour.trend_matricesMethod
function trend_matrices(t)

Arguments

  • t: time

Output

  • E: matrix that solves E*parameters= timeseries
  • F=E†: generalized inverse of E
source
ECCOtour.trend_theta!Method
function trend_theta!(β,diagpathexp,path_out,tecco,γ,F)
get linear trend of potential temperature
from monthly-average ECCOv4r4 gcmarray fields

Arguments

  • β: gcmarray field of trends
  • diagpathexp: where to find files
  • tecco: time stamps of monthly output
  • γ: GCM grid information
  • F: linear estimator of trend

Output

  • β: updated trends
source
ECCOtour.trend_thetaMethod
function trend_theta(diagpathexp,path_out,tecco,γ,F)
get linear trend of potential temperature
from monthly-average ECCOv4r4 gcmarray fields

Arguments

  • diagpathexp: where to find files
  • tecco: time stamps of monthly output
  • γ: GCM grid information
  • F: linear estimator of trend

Output

  • β: updated trends
source
ECCOtour.vars2sigmaMethod
function 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 gcmarrays
  • p::Array{Float64,1} : vertical profile of standard pressures
  • siggrid: σ₁ surface values
  • γ: grid description needed for preallocation
  • splorder: 1-5, order of spline
  • linearinterp: optional logical

Output

  • varsσ::Dict{String,MeshArrays.gcmarray{T,N,Array{T,2}}}: dict of gcmarrays of variables on sigma1 surfaces
source
ECCOtour.velocity2center!Method
velocity2center!(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
source
ECCOtour.velocity2centerMethod
velocity2center(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
source
ECCOtour.wrapdistMethod
function wrapdist
Longitudinal distance in degrees
Passing the date line may be the shortest distance.
source
ECCOtour.write_varsMethod

function write_vars This version writes mdsio output on the native grid Missing: write the accompanying meta file

source
ECCOtour.writeregularpolesMethod

function writeregularpoles(vars::Dict{String,Array{Float32,2}},γ,pathout,filesuffix,filelog,λC,lonatts,ϕC,latatts,z,depthatts)

source
ECCOtour.θbudg2sigmaMethod

function θ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/V
  • pathin::String: path of input
  • pathout::String: path of output
  • fileroots::String: beginning of file names in pathin
  • γ: grid description needed for preallocation
  • splorder::Integer: optional keyword argument, default = 3, order of spline
  • linearinterp::Logical: optional keyword argument, do linear interpolation?, default = false
  • eos::String: optional key argument for equation of state, default = "JMD95"
source
ECCOtour.θbudg2sigma2Method

θ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
source
IsopycnalSurfaces.vars2sigma1Method
function 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 gcmarrays
  • p::Array{Float64,1} : vertical profile of standard pressures
  • sig1grid: σ₁ surface values
  • γ: grid description needed for preallocation
  • splorder: 1-5, order of spline
  • linearinterp: optional logical

Output

  • varsσ::Dict{String,MeshArrays.gcmarray{T,N,Array{T,2}}}: dict of gcmarrays of variables on sigma1 surfaces
source
Statistics.meanMethod
function 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 type
  • weight::MeshArrays.gcmarray{Float32,1,Array{Float32,2}}: weighting variable of gcmarray type
  • isgood::Function: returns true is a value to be used in the mean

Output

  • xbar::Float32: mean value (weighted and filtered)
source
Statistics.meanMethod
function 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 type
  • weight::MeshArrays.gcmarray{Float32,1,Array{Float32,2}}: weighting variable of gcmarray type
  • dryval::Float32: land value (doesn't work for NaN32)

Output

  • xbar::Float32: mean value (unweighted)
source
Statistics.meanMethod
function 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 type
  • dryval::T: land value (doesn't work for NaN32)

Output

  • xbar::T: mean value (unweighted)
source
Statistics.stdMethod
function 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 type
  • xbar::T: mean value
  • dryval::T: land value

Output

  • σx::T: standard deviation
source