MITgcm File Formats
The two main output formats of MITgcm gridded fields are called MDS
and MNC
. In addition, the standard output file
is a text file that records events during the model run, including its successful completion.
Standard Output
The standard output file
can be scanned to collect information about the model run.
MITgcm.monitor
— Functionfunction monitor(config::MITgcm_config)
Call scan_rundir
and show to REPL.
MITgcm.scan_rundir
— Functionscan_rundir(pth::String)
Scan a MITgcm run directory and then, if found, the standard output text file ("output.txt" or "STDOUT.0000") via scan_stdout
.
scan_rundir(config::MITgcm_config)
Scan a MITgcm run directory (joinpath(MC,"run")) and then, if found, the standard output text file ("output.txt" or "STDOUT.0000") via scan_stdout
.
MITgcm.scan_stdout
— Functionscan_stdout(filout::String)
Scan a MITgcm standard output text file ("output.txt" or "STDOUT.0000") and return a NamedTuple of information collected.
- packages : report of packages being compiled and used
- params_time : initial time, model duation, output frequency, etc
- params_grid : type of grid (Curvilinear, Cartesian, ...) and array sizes
- paramsfiles : type of output (usemdsio, use_mnc) and array size (ioSize)
- completed : true / false depending on the end of standard output file (filout)
MITgcm.read_available_diagnostics
— Functionread_available_diagnostics(fldname::String; filename="available_diagnostics.log")
Get the information for a particular variable fldname
from the available_diagnostics.log
text file generated by MITgcm
.
Input Files
Run-time parameters to MITgcm
are provided via text files. This package can read and write them in two formats (standard TOML
format, or the native MITgcm_namelist
format).
Other inputs (binary or NetCDF files) for gridded data, forcing fields, etc can be provided via an input_folder
or downloaded as shown in the setup_ECCO4!
.
MITgcm.read_toml
— Functionread_toml(toml_file::String)
Read toml parameter file into an OrderedDict with Symbol keys, consistent with tracked_parameters.toml
.
using MITgcm, TOML
pth=joinpath(dirname(pathof(MITgcm)),"..","examples","configurations")
toml_file=joinpath(pth,"tutorial_held_suarez_cs.toml")
params=read_toml(toml_file)
Writing parameters to file is straightforward. For example:
MC=MITgcm_config(configuration="tutorial_held_suarez_cs")
setup(MC)
open(tempname()*".toml", "w") do io
TOML.print(io, MC.inputs)
end
read_toml(config_name::Symbol)
Read toml parameter file specified by configuration name.
using MITgcm
params=read_toml(:OCCA2)
MITgcm.MITgcm_namelist
— TypeMITgcm_namelist(groups,params)
Data structure representing a MITgcm namelist file, such as data.pkg
, which contains
groups :: Array{Symbol,1} = Array{Symbol,1}(undef, 0)
params :: Array{OrderedDict{Symbol,Any},1} = Array{OrderedDict{Symbol,Any},1}(undef, 0)
with model parameters (params
) being organized in groups
as done in the files.
using MITgcm
fil=joinpath(MITgcm_path[1],"verification","advect_xy","run","data")
nml=read_namelist(fil)
MITgcm_namelist(nml.groups,nml.params)
MITgcm_namelist(groups=nml.groups,params=nml.params)
MITgcm_namelist(groups=nml.groups)
MITgcm.read_namelist
— Functionread_namelist(fil)
Read a MITgcm
namelist file in native format, parse it, and return as a NamedTuple
.
using MITgcm
testreport(MITgcm_config(configuration="advect_xy"))
fil=joinpath(MITgcm_path[1],"verification","advect_xy","run","data")
namelist=read_namelist(fil)
MITgcm.write_namelist
— Functionwrite_namelist(fil)
Save a MITgcm
namelist file (native format). In the example below, one is read from file, modified, and then saved to a new file using write_namelist.
using MITgcm
fil=joinpath(MITgcm_path[1],"verification","advect_xy","run","data")
nml=read(fil,MITgcm_namelist())
write(fil*"_new",nml)
or write_namelist(fil*"_new",namelist)
.
MITgcm.read_all_namelists
— Functionread_all_namelists(input_path)
Read all MITgcm
namelist files in input_path
, parse them, and return as a NamedTuple of NamedTuples.
using MITgcm; path0=default_path()
input_path=joinpath(path0,"verification","advect_xy","input")
params=read_all_namelists(input_path)
MITgcm.write_all_namelists
— Functionwrite_all_namelists(params,output_path=tempname())
Write all MITgcm
namelist to files in output_path
, from corresponding toml file
using MITgcm
params=read_toml(:OCCA2)
write_all_namelists(params)
MDS Files
MITgcm.read_mdsio
— Functionread_mdsio(fil::String)
Read a single MITgcm
MDSIO-type file (".data" binary + ".meta" text pair), and return as an Array
p="./hs94.cs-32x32x5/run/"
x=read_mdsio(p*"surfDiag.0000000020.002.001.data")
y=read_mdsio(p*"pickup.ckptA.002.001.data")
z=read_mdsio(p*"T.0000000000.002.001.data")
read_mdsio(fil::String,rec::Integer)
Read a single variable / record from a single MITgcm
MDSIO-type file, and return as an Array.
read_mdsio(fil,1)
read_mdsio(fil::String,nam::Symbol)
Read a single variable / record from a single MITgcm
MDSIO-type file, and return as an Array.
read_mdsio(fil,:THETA)
read_mdsio(pth::String,fil::String)
Read a set of MITgcm
MDSIO-type files (".data" binary + ".meta" text pair), combine, and return as an Array. Unlike read_mdsio(fil::String)
where fil
is one complete file name, this method will search within pth
for files that start with fil
.
p="./hs94.cs-32x32x5/run/"
x=read_mdsio(p,"surfDiag.0000000020")
y=read_mdsio(p,"pickup.ckptA")
z=read_mdsio(p,"T.0000000000")
MITgcm.read_meta
— Functionread_meta(metafile)
Read a MITgcm
metadata file, parse it, and return as a NamedTuple
p="./hs94.cs-32x32x5/run/"
meta=read_meta(p*"surfDiag.0000000020.002.001.meta")
pairs(meta)
meta.dimList
read_meta(pth::String,fil::String)
Read a MITgcm
metadata files, parse them, and return as an array of NamedTuple
p="./hs94.cs-32x32x5/run/"
meta=read_meta(p,"surfDiag.0000000020")
pairs(meta[end])
[meta[i].dimList for i in 1:length(meta)]
MNC Files
MITgcm.read_mnc
— Functionread_mnc(pth::String,fil::String,var::String)
Read variable var
from a set of MITgcm
MNC-type files (netcdf files), combine, and return as an Array. This method will search within pth
for files that start with fil
.
Grid Files
MITgcm.GridLoad_mdsio
— FunctionGridLoad_mdsio(myexp::MITgcm_config)
Load grid variables (XC, YC, Depth, etc) from model run directory (rundir
).
GridLoad_mdsio(rundir::String)
Load grid variables (XC, YC, Depth, etc) from model run directory (rundir
).
MITgcm.GridLoad_mnc
— FunctionGridLoad_mnc(γ::gcmgrid)
Load grid variabes (XC, YC, Depth) model run directory (joinpath(rundir,"mnc_test_0001")
).
GridLoad_mnc(myexp::MITgcm_config)
Load grid variables (XC, YC, Depth) from model run directory (joinpath(rundir,"mnc_test_0001")
).
GridLoad_mnc(rundir::String)
Load grid variables (XC, YC, Depth) from model run directory (joinpath(rundir,"mnc_test_0001")
).
MITgcm.ReadNativeGridFiles.GridLoad_native
— FunctionGridLoad_native(path,files,γ)
Load grid variables from native grid files.
path="GRID_LLC90"
files=["tile001.mitgrid","tile002.mitgrid","tile003.mitgrid","tile004.mitgrid","tile005.mitgrid"]
ioSize=[90 1170]
facesSize=[(90, 270), (90, 270), (90, 90), (270, 90), (270, 90)]
γ=gcmgrid(path,"LatLonCap",5,facesSize, ioSize, Float64, read, write)
Γ=GridLoad_native(path,files,γ)
or using another grid
path="llc_1080"
files=[ "llc_001_1080_3240.bin","llc_002_1080_3240.bin",
"llc_003_1080_1080.bin","llc_004_3240_1080.bin","llc_005_3240_1080.bin"]
ioSize=[1080 14040]
facesSize=[(1080, 3240), (1080, 3240), (1080, 1080), (3240, 1080), (3240, 1080)]
and for plotting
using GLMakie
#col=log10.(Γ.RAC); rng=(4,8)
#using MAT
#Depth=get_bathy(path,γ)
#col=write(Depth); rng=(-5000,5000)
XC=write(Γ.XC); YC=write(Γ.YC)
ii=findall((XC.>-80).&(XC.<-10).&(YC.>-10).&(YC.<60))
#ii=findall((XC.>-80).&(XC.<-60).&(YC.>15).&(YC.<35))
#ii=1:prod(γ.ioSize)
scatter(XC[ii],YC[ii],color=col[ii],colorrange=rng,markersize=0.1, markerspace = :data)
Other Output Files
MITgcm.read_flt
— Functionread_flt(dirIn::String,prec::DataType)
Read displacements from MITgcm/pkg/flt output file into a DataFrame.
MITgcm.read_bin
— Functionread_bin(fil::String,kt::Union{Int,Missing},kk::Union{Int,Missing},prec::DataType,mygrid::gcmgrid)
Read model output from binary file and convert to MeshArray. Other methods:
read_bin(fil::String,prec::DataType,mygrid::gcmgrid)
read_bin(fil::String,mygrid::gcmgrid)
MITgcm.read_nctiles
— Functionread_nctiles(fileName,fldName,mygrid; I, eccoVersion4Release4=false, verbose=false)
Read model output from NCTiles file and convert to MeshArray instance. Setting the keyword argument eccoVersion4Release4=true
allows read_nctiles
to read in ECCOv4r4 data which has a different file naming convention to previous versions.
mygrid=GridSpec("LatLonCap")
fileName="nctiles_grid/GRID"
Depth=read_nctiles(fileName,"Depth",mygrid)
hFacC=read_nctiles(fileName,"hFacC",mygrid)
hFacC=read_nctiles(fileName,"hFacC",mygrid,I=(:,:,1))