load_data

moment_kinetics.load_data.get_run_info_no_setupMethod
get_run_info_no_setup(run_dir...; itime_min=1, itime_max=0, itime_skip=1, dfns=false,
                      initial_electron=false)
get_run_info_no_setup((run_dir, restart_index)...; itime_min=1, itime_max=0,
                      itime_skip=1, dfns=false, initial_electron=false)

Get file handles and other info for a single run

run_dir is either the directory to read output from (whose name should be the run_name), or a momentkinetics binary output file. If a file is passed, it is only used to infer the directory and `runname, so it is possible for example to pass a.moments.h5output file and alsodfns=trueand the.dfns.h5` file will be the one actually opened (as long as it exists).

restart_index can be given by passing a Tuple, e.g. ("runs/example", 42) as the positional argument. It specifies which restart to read if there are multiple restarts. If no restart_index is given or if nothing is passed, read all restarts and concatenate them. An integer value reads the restart with that index - -1 indicates the latest restart (which does not have an index).

Several runs can be loaded at the same time by passing multiple positional arguments. Each argument can be a String run_dir giving a directory to read output from or a Tuple (run_dir, restart_index) giving both a directory and a restart index (it is allowed to mix Strings and Tuples in a call).

By default load data from moments files, pass dfns=true to load from distribution functions files, or initial_electron=true and dfns=true to load from initial electron state files.

The itime_min, itime_max and itime_skip options can be used to select only a slice of time points when loading data. In makie_post_process these options are read from the input (if they are set) before get_run_info_no_setup() is called, so that the run_info returned can be passed to makie_post_processing.setup_makie_post_processing_input!(), to be used for defaults for the remaining options. If either itime_min or itime_max are ≤0, their values are used as offsets from the final time index of the run.

source
moment_kinetics.load_data.get_variableFunction
get_variable(run_info::Tuple, variable_name; kwargs...)
get_variable(run_info, variable_name; kwargs...)

Get an array (or Tuple of arrays, if run_info is a Tuple) of the data for variable_name from run_info.

Some derived variables need to be calculated from the saved output, not just loaded from file (with postproc_load_variable). This function takes care of that calculation, and handles the case where run_info is a Tuple (which postproc_load_data does not handle).

kwargs... are passed through to postproc_load_variable().

source
moment_kinetics.load_data.get_z_derivativeMethod
get_z_derivative(run_info, variable_name; kwargs...)

Get (i.e. load or calculate) variable_name from run_info and calculate its z-derivative. Returns the z-derivative

kwargs... are passed through to get_variable().

source
moment_kinetics.load_data.load_coordinate_dataMethod
load_coordinate_data(fid, name; printout=false, irank=nothing, nrank=nothing,
                     run_directory=nothing, ignore_MPI=true)

Load data for the coordinate name from a file-handle fid.

Returns (coord, spectral, chunk_size). coord is a coordinate object. spectral is the object used to implement the discretization in this coordinate. chunk_size is the size of chunks in this coordinate that was used when writing to the output file.

If printout is set to true a message will be printed when this function is called.

If irank and nrank are passed, then the coord and spectral objects returned will be set up for the parallelisation specified by irank and nrank, rather than the one implied by the output file.

Unless ignore_MPI=false is passed, the returned coordinates will be created without shared memory scratch arrays (ignore_MPI=true will be passed through to define_coordinate).

source
moment_kinetics.load_data.load_distributed_electron_pdf_sliceMethod

Read a slice of an electron distribution function

run_names is a tuple. If it has more than one entry, this means that there are multiple restarts (which are sequential in time), so concatenate the data from each entry together.

The slice to take is specified by the keyword arguments.

source
moment_kinetics.load_data.load_distributed_ion_pdf_sliceMethod

Read a slice of an ion distribution function

run_names is a tuple. If it has more than one entry, this means that there are multiple restarts (which are sequential in time), so concatenate the data from each entry together.

The slice to take is specified by the keyword arguments.

source
moment_kinetics.load_data.load_distributed_neutral_pdf_sliceMethod

Read a slice of a neutral distribution function

run_names is a tuple. If it has more than one entry, this means that there are multiple restarts (which are sequential in time), so concatenate the data from each entry together.

The slice to take is specified by the keyword arguments.

source
moment_kinetics.load_data.postproc_load_variableMethod
postproc_load_variable(run_info, variable_name; it=nothing, is=nothing,
                       ir=nothing, iz=nothing, ivperp=nothing, ivpa=nothing,
                       ivzeta=nothing, ivr=nothing, ivz=nothing)

Load a variable

run_info is the information about a run returned by makie_post_processing.get_run_info().

variable_name is the name of the variable to load.

The keyword arguments it, is, ir, iz, ivperp, ivpa, ivzeta, ivr, and ivz can be set to an integer or a range (e.g. 3:8 or 3:2:8) to select subsets of the data. Only the data for the subset requested will be loaded from the output file (mostly - when loading fields or moments from runs which used parallel_io = false, the full array will be loaded and then sliced).

source
moment_kinetics.load_data.read_Dict_from_sectionMethod
read_Dict_from_section(file_or_group, section_name; ignore_subsections=false)

Read information from section_name in file_or_group, returning a Dict.

By default, any subsections are included as nested Dicts. If ignore_subsections=true they are ignored.

source
moment_kinetics.load_data.read_distributed_zr_data!Method

Read data which is a function of (z,r,t) or (z,r,species,t)

run_names is a tuple. If it has more than one entry, this means that there are multiple restarts (which are sequential in time), so concatenate the data from each entry together.

source