High-level docstrings

This section lists the main interfaces for interacting with codes that use Jutul as their foundation.

Setting up models

Jutul.SimulationModelType
SimulationModel(domain, system; <kwarg>)

Instantiate a model for a given system discretized on the domain.

source
SimulationModel(g::JutulMesh, system; discretization = nothing, kwarg...)

Type that defines a simulation model - everything needed to solve the discrete equations.

The minimal setup requires a JutulMesh that defines topology together with a JutulSystem that imposes physical laws.

source
Jutul.MultiModelType
MultiModel(models)
MultiModel(models, :SomeLabel)

A model variant that is made up of many named submodels, each a fully realized SimulationModel.

models should be a NamedTuple or Dict{Symbol, JutulModel}.

source

Model API Getters

Jutul.get_secondary_variablesFunction
get_secondary_variables(model::SimulationModel)

Get the secondary variable definitions (as OrderedDict) for a given model.

Secondary variables are variables that can be computed from the primary variables together with the parameters.

source
Jutul.get_primary_variablesFunction
get_primary_variables(model::SimulationModel)

Get the primary variable definitions (as OrderedDict) for a given model.

Primary variables are sometimes referred to as solution variables or primary unknowns. The set of primary variables completely determines the state of the system together with the parameters.

source
Jutul.get_parametersFunction
get_parameters(model::SimulationModel)

Get the parameter definitions (as OrderedDict) for a given model.

Parameters are defined as static values in a forward simulation that combine with the primary variables to compute secondary variables and model equations.

source

Setters

Jutul.set_secondary_variables!Function
set_secondary_variables!(model, varname = vardef)
set_secondary_variables!(model, varname1 = vardef1, varname2 = vardef2)

Set a secondary variable with name varname to the definition vardef (adding if it does not exist)

source
Jutul.set_primary_variables!Function
set_primary_variables!(model, varname = vardef)
set_primary_variables!(model, varname1 = vardef1, varname2 = vardef2)

Set a primary variable with name varname to the definition vardef (adding if it does not exist)

source
Jutul.set_parameters!Function
set_parameters!(model, parname = pardef)

Set a parameter with name varname to the definition vardef (adding if it does not exist)

source

Various

Systems and domains

Jutul.JutulContextType

Abstract type for the context Jutul should execute in (matrix formats, memory allocation, etc.)

source
Jutul.DiscretizedDomainType
DiscretizedDomain(domain, disc = nothing)

A type for a discretized domain of some other domain or mesh. May contain one or more discretizations as-needed to write equations.

source

Grids/meshes and geometry

Jutul.CartesianMeshType
CartesianMesh(dims, [Δ, [origin]])

Create a Cartesian mesh with dimensions specified by the Tuple dims.

Arguments

  • dims::Tuple: Number of grid cells in each direction. For example, (nx, ny) will give a 2D grids with nx cells in the x-direction.
  • Δ::Tuple=Tuple(ones(length(dims))): Equal length to dims. First option: A

Tuple of scalars where each entry is the length of each cell in that direction. For example, specifying (Δx, Δy) for a uniform grid with each grid cell having area ofΔx*Δy. Second option: ATuple` of vectors where each entry contains the cell sizes in the direction.

  • origin=zeros(length(dims)): The origin of the first corner in the grid.

Examples

Generate a uniform 3D mesh that discretizes a domain of 2 by 3 by 5 units with 3 by 5 by 2 cells:

julia> CartesianMesh((3, 5, 2), (2.0, 3.0, 5.0))
CartesianMesh (3D) with 3x5x2=30 cells

Generate a non-uniform 2D mesh:

julia> CartesianMesh((2, 3), ([1.0, 2.0], [0.1, 3.0, 2.5]))
CartesianMesh (2D) with 2x3x1=6 cells
source
Jutul.TwoPointFiniteVolumeGeometryType
TwoPointFiniteVolumeGeometry(neighbors, areas, volumes, normals, cell_centers, face_centers)

Store two-point geometry information for a given list of neighbors specified as a 2 by n matrix where n is the number of faces such that face i connectes cells N[1, i] and N[2, i].

The two-point finite-volume geometry contains the minimal set of geometry information required to compute standard finite-volume discretizations.

source
Jutul.number_of_cellsFunction
number_of_cells(D::Union{DataDomain, DiscretizedDomain})

Get the number of cells in a DataDomain or DiscretizedDomain.

source
number_of_cells(g)::Integer

Get the number of cells in a mesh.

source

Set up of system

Jutul.setup_stateFunction
setup_state(model::JutulModel, name1 = value1, name2 = value2)

Set up a state for a given model with values for the primary variables defined in the model. Normally all primary variables must be initialized in this way.

Arguments

  • name=value: The name of the primary variable together with the value(s) used to initialize the primary variable.

A scalar (or short vector of the right size for VectorVariables) will be repeated over the entire domain, while a vector (or matrix for VectorVariables) with length (number of columns for VectorVariables) equal to the entity count (for example, number of cells for a cell variable) will be used directly.

Note: You likely want to overload [setup_state!]@ref for a custom model instead of setup_state

source
Jutul.setup_parametersFunction
setup_parameters(model::JutulModel; name = value)

Set up a parameter storage for a given model with values for the parameter defined in the model.

Arguments

  • name=value: The name of the parameter together with the value(s) of the parameter.

A scalar (or short vector of the right size for VectorVariables) will be repeated over the entire domain, while a vector (or matrix for VectorVariables) with length (number of columns for VectorVariables) equal to the entity count (for example, number of cells for a cell variable) will be used directly.

source
Jutul.setup_state_and_parametersFunction
state, prm = setup_state_and_parameters(model, init)

Simultaneously set up state and parameters from a single init file (typically a Dict containing values that might either be initial values or parameters)

source
Jutul.setup_forcesFunction
setup_forces(model::JutulModel; force_name = force_value)

Set up forces for a given model. Keyword arguments varies depending on what the model supports.

source

Simulator interfaces

Used to run simulations once a model has been set up.

Jutul.simulateFunction
simulate(state0, model, timesteps, parameters = setup_parameters(model))
simulate(state0, model, timesteps, info_level = 3)
simulate(state0, model, timesteps; <keyword arguments>)

Simulate a set of timesteps with model for the given initial state0 and optionally specific parameters. Additional keyword arguments are passed onto simulator_config and simulate!. This interface is primarily for convenience, as all storage for the simulator is allocated upon use and discared upon return. If you want to perform multiple simulations with the same model it is advised to instead instantiate Simulator and combine it with simulate!.

Arguments

  • state0::Dict: initial state, typically created using setup_state for the model in use.
  • model::JutulModel: model that describes the discretized system to solve, for example SimulationModel or MultiModel.
  • timesteps::AbstractVector: Vector of desired report steps. The simulator will perform time integration until sum(timesteps) is reached, providing outputs at the end of each report step.
  • parameters=setup_parameters(model): Optional overrides the default parameters for the model.
  • forces=nothing: Either nothing (for no forces), a single set of forces from setup_forces(model) or a Vector of such forces with equal length to timesteps.
  • restart=nothing: If an integer is provided, the simulation will attempt to restart from that step. Requires that output_path is provided here or in the config.
  • config=simulator_config(model): Configuration Dict that holds many fine grained settings for output, linear solver, time-steps, outputs etc.

Additional arguments are passed onto simulator_config.

See also simulate!, Simulator, SimulationModel, simulator_config.

source
simulate(state0, sim::JutulSimulator, timesteps::AbstractVector; parameters = nothing, kwarg...)

Simulate a set of timesteps with simulator for the given initial state0 and optionally specific parameters.

source
Jutul.simulate!Function
simulate!(sim::JutulSimulator, timesteps::AbstractVector;
    forces = nothing,
    config = nothing,
    initialize = true,
    restart = nothing,
    state0 = nothing,
    parameters = nothing,
    kwarg...
)

Non-allocating (or perhaps less allocating) version of simulate!.

Arguments

  • initialize=true: Perform internal updates as if this is the first time

See also simulate for additional supported input arguments.

source
Jutul.solve_timestep!Function
solve_timestep!(sim, dT, forces, max_its, config; <keyword arguments>)

Internal function for solving a single time-step with fixed driving forces.

Arguments

  • sim: Simulator instance.
  • dT: time-step to be solved
  • forces: Driving forces for the time-step
  • max_its: Maximum number of steps/Newton iterations.
  • config: Configuration for solver (typically output from simulator_config).

Note: This function is exported for fine-grained simulation workflows. The general simulate interface is both easier to use and performs additional validation.

source

Configure a simulator:

Jutul.SimulatorType
Simulator(model; <kwarg>)

Set up a simulator object for a model that can be used by simulate!. To avoid manually instantiating the simulator, the non-mutating simulate interface can be used instead.

source
Jutul.simulator_configFunction
simulator_config(sim; info_level = 3, linear_solver = GenericKrylov())

Set up a simulator configuration object that can be passed onto simulate!.

There are many options available to configure a given simulator. The best way to get an overview of these possible configuration options is to instatiate the config without any arguments and inspecting the resulting table by calling simulator_config(sim) in the REPL.

source
Jutul.JutulConfigType

JutulConfig(name = nothing)

A configuration object that acts like a Dict{Symbol,Any} but contains additional data to limit the valid keys and values to those added by add_option!

source
Jutul.add_option!Function
add_option!(opts::JutulConfig, :my_cool_option, 3, "My option has this brief description")

Add an option to existing JutulConfig structure. Additional currently undocumented keyword arguments can be used to restrict valid types and values.

source

Sensitivities, adjoints and optimization

Jutul.solve_adjoint_sensitivitiesFunction
solve_adjoint_sensitivities(model, states, reports_or_timesteps, G; extra_timing = false, state0 = setup_state(model), forces = setup_forces(model), raw_output = false, kwarg...)

Compute sensitivities of model parameter with name target for objective function G.

The objective function is at the moment assumed to be a sum over all states on the form: obj = Σₙ G(model, state, dt_n, n, forces_for_step_n)

Solves the adjoint equations: For model equations F the gradient with respect to parameters p is ∇ₚG = Σₙ (∂Fₙ / ∂p)ᵀ λₙ where n ∈ [1, N]. Given Lagrange multipliers λₙ from the adjoint equations (∂Fₙ / ∂xₙ)ᵀ λₙ = - (∂J / ∂xₙ)ᵀ - (∂Fₙ₊₁ / ∂xₙ)ᵀ λₙ₊₁ where the last term is omitted for step n = N and G is the objective function.

source
Jutul.solve_adjoint_sensitivities!Function
solve_adjoint_sensitivities!(∇G, storage, states, state0, timesteps, G; forces = setup_forces(model))

Non-allocating version of solve_adjoint_sensitivities.

source
Jutul.solve_numerical_sensitivitiesFunction
solve_numerical_sensitivities(model, states, reports, G, target;
                                            forces = setup_forces(model),
                                            state0 = setup_state(model),
                                            parameters = setup_parameters(model),
                                            epsilon = 1e-8)

Compute sensitivities of model parameter with name target for objective function G.

This method uses numerical perturbation and is primarily intended for testing of solve_adjoint_sensitivities.

source
Jutul.setup_parameter_optimizationFunction
setup_parameter_optimization(model, state0, param, dt, forces, G, opt_cfg = optimization_config(model, param);
                                                        grad_type = :adjoint,
                                                        config = nothing,
                                                        print = 1,
                                                        copy_case = true,
                                                        param_obj = false,
                                                        kwarg...)

Set up function handles for optimizing the case defined by the inputs to simulate together with a per-timestep objective function G.

Generally calling either of the functions will mutate the data Dict. The options are: Fo(x) -> evaluate objective dFo(dFdx, x) -> evaluate gradient of objective, mutating dFdx (may trigger evaluation of Fo) Fand_dF(F, dFdx, x) -> evaluate F and/or dF. Value of nothing will mean that the corresponding entry is skipped.

source

Linear solvers

Jutul.GenericKrylovType
GenericKrylov(solver = :gmres; preconditioner = nothing; <kwarg>)

Solver that wraps Krylov.jl with support for preconditioning.

source
Jutul.LUSolverType
LUSolver(; reuse_memory = true, check = true, max_size = 50000)

Direct solver that calls lu directly. Direct solvers are highly accurate, but are costly in terms of memory usage and execution speed for larger systems.

source

Preconditioners

Execution contexts