High-level docstrings
This section lists the main interfaces for interacting with codes that use Jutul as their foundation.
Setting up models
Jutul.SimulationModel
— TypeSimulationModel(domain, system; <kwarg>)
Instantiate a model for a given system
discretized on the domain
.
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.
Jutul.MultiModel
— TypeMultiModel(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}
.
Model API
Getters
Jutul.get_secondary_variables
— Functionget_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.
Jutul.get_primary_variables
— Functionget_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
.
Jutul.get_parameters
— Functionget_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.
Jutul.get_variables
— Functionget_variables(model::SimulationModel)
Get all variable definitions (as OrderedDict
) for a given model
.
This is the union of get_secondary_variables
and get_primary_variables
.
Setters
Jutul.set_secondary_variables!
— Functionset_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)
Jutul.set_primary_variables!
— Functionset_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)
Jutul.set_parameters!
— Functionset_parameters!(model, parname = pardef)
Set a parameter with name varname
to the definition vardef
(adding if it does not exist)
Various
Jutul.number_of_degrees_of_freedom
— FunctionTotal number of degrees of freedom for a model, over all primary variables and all entities.
Jutul.number_of_values
— FunctionTotal number of values for a model, for a given type of variables over all entities
Systems and domains
Jutul.JutulSystem
— TypeAbstract type for the physical system to be solved.
Jutul.JutulContext
— TypeAbstract type for the context Jutul should execute in (matrix formats, memory allocation, etc.)
Jutul.JutulDomain
— TypeAbstract type for domains where equations can be defined
Jutul.DataDomain
— TypeDataDomain(domain::JutulDomain; property1 = p1, property2 = p2, ...)
A wrapper around a domain that allows for storing of entity-associated data.
Example:
# Grid with 6 cells and 7 interior faces
g = CartesianMesh((2, 3))
d = DataDomain(g)
d[:cell_vec] = rand(6) #ok, same as:
d[:cell_vec, Cells()] = rand(6) #ok
d[:cell_vec, Faces()] = rand(6) #not ok!
d[:face_vec, Faces()] = rand(7) #ok!
# Can also add general arrays if last dimension == entity dimension
d[:cell_vec, Cells()] = rand(10, 3, 6) #ok
# Can add general data too, but needs to be specified
d[:not_on_face_or_cell, nothing] = rand(3) # also ok
Jutul.DiscretizedDomain
— TypeDiscretizedDomain(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.
Jutul.JutulDiscretization
— TypeAbstract type for a Jutul discretization
Set up of system
Jutul.setup_state
— Functionsetup_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
Jutul.setup_parameters
— Functionsetup_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.
Jutul.setup_state_and_parameters
— Functionstate, 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)
Jutul.setup_forces
— Functionsetup_forces(model::JutulModel; force_name = force_value)
Set up forces for a given model. Keyword arguments varies depending on what the model supports.
Simulator interfaces
Used to run simulations once a model has been set up.
Jutul.simulate
— Functionsimulate(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 usingsetup_state
for themodel
in use.model::JutulModel
: model that describes the discretized system to solve, for exampleSimulationModel
orMultiModel
.timesteps::AbstractVector
: Vector of desired report steps. The simulator will perform time integration untilsum(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
: Eithernothing
(for no forces), a single set of forces fromsetup_forces(model)
or aVector
of such forces with equal length totimesteps
.restart=nothing
: If an integer is provided, the simulation will attempt to restart from that step. Requires thatoutput_path
is provided here or in theconfig
.config=simulator_config(model)
: ConfigurationDict
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
.
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.
Jutul.simulate!
— Functionsimulate!(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.
Jutul.solve_timestep!
— Functionsolve_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 solvedforces
: Driving forces for the time-stepmax_its
: Maximum number of steps/Newton iterations.config
: Configuration for solver (typically output fromsimulator_config
).
Note: This function is exported for fine-grained simulation workflows. The general simulate
interface is both easier to use and performs additional validation.
Configure a simulator:
Jutul.Simulator
— TypeSimulator(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.
Jutul.simulator_config
— Functionsimulator_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.
Jutul.JutulConfig
— TypeJutulConfig(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!
Jutul.add_option!
— Functionadd_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.
Sensitivities, adjoints and optimization
Jutul.solve_adjoint_sensitivities
— Functionsolve_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.
Jutul.solve_adjoint_sensitivities!
— Functionsolve_adjoint_sensitivities!(∇G, storage, states, state0, timesteps, G; forces)
Non-allocating version of solve_adjoint_sensitivities
.
Jutul.solve_numerical_sensitivities
— Functionsolve_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
.
Jutul.setup_parameter_optimization
— Functionsetup_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
. The objective should be on the form of sum over all steps, where each element of the sum is evaluated by model, state, dt, step_no, forces
.
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.
Generic optimization interface
Jutul.DictOptimization.DictParameters
— TypeDictParameters(parameters)
DictParameters(parameters::AbstractDict, setup_function = missing;
strict = true,
verbose = true,
active_type = Float64
)
Set up a DictParameters
object for optimization. Optionally, the setup function that takes an instance with the same keys as parameters
together with a step_info
dictionary can be provided. The setup function should return a JutulCase
set up from the parameters in the Dict.
Optional keyword arguments:
strict
: If true, the optimization will throw an error if any of the parameters are not set with at least one of the upper or lower bounds.verbose
: If true, the optimization will print information about the optimization process.active_type
: The type of the parameters that are considered active in the optimization. Defaults toFloat64
. This is used to determine which parameters are active and should be optimized. This means that all entries (and entries in nested dictionaries) of theparameters
dictionary must be of this type or an array with this type as element type.
Jutul.DictOptimization.free_optimization_parameter!
— Functionfree_optimization_parameter!(dopt, "parameter_name", rel_min = 0.01, rel_max = 100.0)
free_optimization_parameter!(dopt, ["dict_name", "parameter_name"], abs_min = -8.0, abs_max = 7.0)
Free an existing parameter for optimization in the DictParameters
object. This will allow the parameter to be optimized through a call to optimize
.
Nesting structures
If your DictParameters
has a nesting structure, you can use a vector of strings or symbols to specify the parameter name, e.g. ["dict_name", "parameter_name"]
to access the parameter located at ["dict_name"]["parameter_name"]
.
Setting limits
The limits can be set using the following keyword arguments:
abs_min
: Absolute minimum value for the parameter. If not set, no absolute minimum will be applied.abs_max
: Absolute maximum value for the parameter. If not set, no absolute maximum will be applied.rel_min
: Relative minimum value for the parameter. If not set, no relative minimum will be applied.rel_max
: Relative maximum value for the parameter. If not set, no relative maximum will be applied.
For either of these entries it is possible to pass either a scalar, or an array. If an array is passed, it must have the same size as the parameter being set.
Note that if dopt.strict
is set to true
, at least one of the upper or lower bounds must be set for free parameters. If dopt.strict
is set to false
, the bounds are optional and the DictParameters
object can be used to compute sensitivities, but the built-in optimization routine assumes that finite limits are set for all parameters.
Other keyword arguments
initial
: Initial value for the parameter. If not set, the current value indopt.parameters
will be used.scaler=missing
: Optional scaler for the parameter. If not set, no scaling will be applied. Available scalers are:log
,:exp
. The scaler will be appliedlumping=missing
: Optional lumping array for the parameter. If not set, no lumping will be applied. The lumping array should have the same size as the parameter and contain positive integers. The lumping array defines groups of indices that should be lumped together, i.e. the same value will be used for all indices in the same group. The lumping array should contain all integers from 1 to the maximum value in the array, and all indices in the same group should have the same value in the initial parameter, otherwise an error will be thrown.
Jutul.DictOptimization.freeze_optimization_parameter!
— Functionfreeze_optimization_parameter!(dopt, "parameter_name")
freeze_optimization_parameter!(dopt, ["dict_name", "parameter_name"])
freeze_optimization_parameter!(dopt::DictParameters, parameter_name, val = missing)
Freeze an optimization parameter in the DictParameters
object. This will remove the parameter from the optimization targets and set its value to val
if provided. Any limits/lumping/scaling settings for this parameter will be removed.
Jutul.DictOptimization.set_optimization_parameter!
— Functionset_optimization_parameter!(dopt::DictParameters, parameter_name, value)
Set a specific optimization parameter in the DictParameters
object. This function will update the value of the parameter in the dopt.parameters
dictionary.
Jutul.DictOptimization.optimize
— Functionoptimized_dict = optimize(dopt, objective)
optimize(dopt::DictParameters, objective, setup_fn = dopt.setup_function;
grad_tol = 1e-6,
obj_change_tol = 1e-6,
max_it = 25,
opt_fun = missing,
maximize = false,
simulator = missing,
config = missing,
solution_history = false,
backend_arg = (
use_sparsity = false,
di_sparse = true,
single_step_sparsity = false,
do_prep = true,
),
kwarg...
)
Optimize parameters defined in a DictParameters
object using the provided objective function. At least one variable has to be declared to be free using free_optimization_parameter!
prior to calling the optimizer.
Arguments
dopt::DictParameters
: Container with parameters to optimizeobjective
: The objective function to minimize (or maximize)setup_fn
: Function to set up the optimization problem. Defaults todopt.setup_function
Keyword Arguments
grad_tol
: Gradient tolerance for stopping criterionobj_change_tol
: Objective function change tolerance for stopping criterionmax_it
: Maximum number of iterationsopt_fun
: Optional custom optimization function. If missing, L-BFGS will be usedmaximize
: Set totrue
to maximize the objective instead of minimizingsimulator
: Optional simulator object used in forward simulationsconfig
: Optional configuration for the setupsolution_history
: Iftrue
, stores all intermediate solutionsbackend_arg
: Options for the autodiff backend:use_sparsity
: Enable sparsity detection for the objective functiondi_sparse
: Use sparse differentiationsingle_step_sparsity
: Enable single step sparsity detection (if sparsity does not change during timesteps)do_prep
: Perform preparation step
Returns
The optimized parameters as a dictionary.
Notes
- The function stores the optimization history and optimized parameters in the input
dopt
object. - If
solution_history
istrue
, intermediate solutions are stored indopt.history.solutions
. - The default optimization algorithm is L-BFGS with box constraints.
Jutul.DictOptimization.parameters_gradient
— Functionparameters_gradient(dopt::DictParameters, objective, setup_fn = dopt.setup_function)
Compute the gradient of the objective function with respect to the parameters defined in the DictParameters
object. This function will return the gradient as a dictionary with the same structure as the input parameters, where each entry is a vector of gradients for each parameter. Only gradients with respect to free parameters will be computed.
Linear solvers
Jutul.GenericKrylov
— TypeGenericKrylov(solver = :gmres; preconditioner = nothing; <kwarg>)
Solver that wraps Krylov.jl
with support for preconditioning.
Jutul.LUSolver
— TypeLUSolver(; 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.
Preconditioners
Jutul.AMGPreconditioner
— TypeAMG on CPU (Julia native)
Jutul.ILUZeroPreconditioner
— TypeILU(0) preconditioner on CPU
Jutul.JacobiPreconditioner
— TypeDamped Jacobi preconditioner on CPU
Jutul.SPAI0Preconditioner
— TypeSparse Approximate Inverse preconditioner of lowest order – SPAI(0)
Jutul.LUPreconditioner
— TypeFull LU factorization as preconditioner (intended for smaller subsystems)
Jutul.GroupWisePreconditioner
— TypeMulti-model preconditioners
Execution contexts
Jutul.DefaultContext
— TypeDefault context
Jutul.ParallelCSRContext
— TypeA context that uses a CSR sparse matrix format together with threads. Experimental.