MeshArrays.jl documentation

MeshArrays.jl documentation

MeshArrays.jl primarily defines composite types that embed inter-connected array collections within a struct and provides an exchange function that effectively transfers data between connected arrays. It was originally introduced, as gcmfaces.jl, in this JuliaCon-2018 presentation (see below for notebooks). Note: even though MeshArrays.jl is registered, documented, archived, and routinely tested, it is still regarded as a preliminary implementation.

Contents

Installation

using Pkg
Pkg.add("MeshArrays")
Pkg.test("MeshArrays")

Note: Julia's package manager is documented here within docs.julialang.org.

Use examples

The JuliaCon-2018 presentation relied on two Jupyter notebooks privided in the JuliaCon2018Notebooks repo and pre-defined grids that that can be downloaded as follows:

git clone https://github.com/gaelforget/GRID_CS32
git clone https://github.com/gaelforget/GRID_LLC90

At the command line, one can reproduce the same computation as follows:

using MeshArrays
isdir("GRID_LLC90") ? GridVariables=GCMGridLoad(GCMGridSpec("LLC90")) : 
					  GridVariables=GCMGridOnes("cs",6,100);
(Rini,Rend,DXCsm,DYCsm)= MeshArrays.demo2(GridVariables);

For plotting directions, see the help section (?MeshArrays.demo2). Additional demos are also provided (e.g., see ?MeshArrays.demo1, ?MeshArrays.demo3).

Main Features

MeshArrays.jl composite types contain array collections where arrays are typically inter-connected at their edges. It provides exchange() functions that transfer data between neighbor arrays so that the user can extend the domain of computation as needed e.g. for interpolation or to compute spatial derivatives.

The composite types specify how each array collection forms a mesh and provide information to allow exchange to dispatch to the appropriate method. Various configurations that are commonly used in Earth System Models are implemented using the concrete type called MeshArray. In fact MeshArray is an alias for other types that can thus used interchangeably (initially: gcmfaces or gcmarray). These types contain a gcmgrid instance that provides grid specifications (incl. location on disk). Embedded arrays, or meshes, each represent a subdomain within, e.g., an Earth System Model grid. For example, gcmarray represents a variable x as an array of 2D ragged arrays x.f.

This package derives from a previous Matlab / Octave package, gcmfaces introduced in Forget et al., 2015 (doi:10.5194/gmd-8-3071-2015). GCM is an acronym for General Circulation Model, or Global Climate Model (see this wikipedia entry), and faces can mean meshes, arrays, facets, or subdomains (x.f in a MeshArray instance x).

API index

API details

AbstractMeshArray{T, N}

Subtype of AbstractArray{T, N}

source
gcmgrid

gcmgrid data structure. Available constructors:

gcmgrid(path::String, class::String, nFaces::Int,
        fSize::Array{NTuple{2, Int},1}, ioSize::Array{Int64,2},
        ioPrec::Type, read::Function, write::Function)
source
GCMGridLoad(mygrid::gcmgrid)

Return a Dict of grid variables read from files located in mygrid.path (see ?GCMGridSpec).

Based on the MITgcm naming convention, grid variables are:

  • XC, XG, YC, YG, AngleCS, AngleSN, hFacC, hFacS, hFacW, Depth.
  • RAC, RAW, RAS, RAZ, DXC, DXG, DYC, DYG.
  • DRC, DRF, RC, RF (one-dimensional)
source
GCMGridOnes(grTp,nF,nP)

Define all-ones grid variables instead of using GCMGridSpec & GCMGridLoad.

source
GCMGridSpec(gridName)

Return a gmcgrid specification that provides grid files path, class, nFaces, ioSize, facesSize, ioPrec, & a read function (not yet) using hard-coded values for "LLC90", "CS32", "LL360" (for now).

source
LatitudeCircles(LatValues,GridVariables::Dict)

Compute integration paths that follow latitude circles

source
ThroughFlow(VectorField,IntegralPath,GridVariables::Dict)

Compute transport through an integration path

source
convergence(uFLD::MeshArray,vFLD::MeshArray)

Compute convergence of a vector field

source
exchange(fld::MeshArray)

Exchange / transfer data between neighboring arrays. Other methods are

exchange(fld::MeshArray,N::Integer)
exchange(u::MeshArray,v::MeshArray)
exchange(u::MeshArray,v::MeshArray,N::Integer)
source
MeshArrays.fijindMethod.
fijind(A::gcmfaces,ij::Int)

Compute face and local indices (f,j,k) from global index (ij).

(needed in other types?)

source
MeshArrays.findtilesFunction.
findtiles(ni,nj,grid="llc90")

Return a MeshArray map of tile indices for tile size ni,nj

source
gradient(inFLD::MeshArray,GridVariables::Dict)

Compute spatial derivatives. Other methods:

gradient(inFLD::MeshArray,GridVariables::Dict,doDIV::Bool)
gradient(inFLD::MeshArray,iDXC::MeshArray,iDYC::MeshArray)
source
MeshArrays.maskMethod.
mask(fld::MeshArray, val::Number)

Replace non finite values with val. Other methods:

mask(fld::MeshArray)
mask(fld::MeshArray, val::Number, noval::Number)
source
nFacesEtc(a::gcmarray)

Return nFaces, n3 (1 in 2D case; >1 otherwise)

source
MeshArrays.smoothMethod.
smooth(FLD::MeshArray,DXCsm::MeshArray,DYCsm::MeshArray,GridVariables::Dict)

Smooth out scales below DXCsm / DYCsm via diffusion

source
gcmarray{T, N}

gcmarray data structure. Available constructors:

gcmarray{T,N}(grid::gcmgrid,f::Array{Array{T,N},1},
         fSize::Array{NTuple{N, Int}},fIndex::Array{Int,1})
gcmarray(grid::gcmgrid,
         fSize::Array{NTuple{N, Int}},fIndex::Array{Int,1})

gcmarray(grid::gcmgrid,::Type{T})
gcmarray(grid::gcmgrid,::Type{T},n3::Int)
source
gcmfaces{T, N}

gcmfaces data structure. Available constructors:

gcmfaces{T,N}(grid::gcmgrid,f::Array{Array{T,N},1},
         fSize::Array{NTuple{N, Int}}, aSize::NTuple{N,Int})
gcmfaces(grid::gcmgrid,v1::Array{Array{T,N},1}) where {T,N}

gcmfaces(grid::gcmgrid,::Type{T},
         fSize::Array{NTuple{N, Int}}, aSize::NTuple{N,Int}) where {T,N}
gcmfaces(grid::gcmgrid)
gcmfaces()
source
gcmsubset{T, N}

gcmsubset data structure for subsets of gcmfaces. Available constructors:

gcmsubset{T,N}(grid::gcmgrid,f::Array{Array{T,N},1},
               fSize::Array{NTuple{N, Int}},aSize::NTuple{N, Int},
               i::Array{Array{T,N},1},iSize::Array{NTuple{N, Int}})
gcmsubset(grid::gcmgrid,::Type{T},fSize::Array{NTuple{N, Int}},
          aSize::NTuple{N,Int},dims::NTuple{N,Int}) where {T,N}
source
Base.readMethod.
read(fil::String,x::MeshArray)

Read binary file to MeshArray. Other methods:

read(xx::Array,x::MeshArray) #from Array
source
Base.writeMethod.
write(fil::String,x::MeshArray)

Write MeshArray to binary file. Other methods:

write(xx::Array,x::MeshArray) #to Array
source
MeshArrays.demo1Method.
demo1(gridChoice::String)

Demonstrate basic functionalities (load grid, arithmetic, exchange, gradient, etc.). Call sequence:

!isdir("GRID_LLC90") ? error("missing files") : nothing

(D,Dexch,Darr,DD)=MeshArrays.demo1("LLC90");
source
MeshArrays.demo2Method.
demo2()

Demonstrate higher level functions using smooth. Call sequence:

(Rini,Rend,DXCsm,DYCsm)=MeshArrays.demo2();

include(joinpath(dirname(pathof(MeshArrays)),"gcmfaces_plot.jl"))
qwckplot(Rini,"raw noise")
qwckplot(Rend,"smoothed noise")
source
MeshArrays.demo3Method.
demo3()

Demonstrate ocean transport computations. Call sequence:

(UV,LC,Tr)=MeshArrays.demo3();
using Plots; plot(Tr/1e6,title="meridional transport")

include(joinpath(dirname(pathof(MeshArrays)),"gcmfaces_plot.jl"))
qwckplot(UV["U"][:,:,1,1],"U component (note varying face orientations)")
qwckplot(UV["V"][:,:,1,1],"V component (note varying face orientations)")
source
getindexetc(A::gcmarray, I::Vararg{_}) where {T,N}

Same as getindex but also returns the face size and index

(needed in other types?)

source