Documentation from Jutul.jl
JutulDarcy.jl
builds upon Jutul.jl
. You can use JutulDarcy.jl
without knowing the inner workings of Jutul.jl
, but if you want to dive under the hood these docstrings may be helpful.
Jutul.BlockMajorLayout Type
Same as EntityMajorLayout
, but the system is a sparse matrix where each entry is a small dense matrix.
For a test system with primary variables P, S and equations E1, E2 and two cells this will give a diagonal of length 2: [(∂E1/∂p)₁ (∂E1/∂S)₁ ; (∂E2/∂p)₁ (∂E2/∂S)₁] [(∂E1/∂p)₂ (∂E1/∂S)₂ ; (∂E2/∂p)₂ (∂E2/∂S)₂]
Jutul.BoundaryFaces Type
Entity for faces on the boundary (faces that are only connected to a single Cells
)
Jutul.Cells Type
Entity for Cells (closed volumes with averaged properties for a finite-volume solver)
Jutul.CoarseMesh Method
CoarseMesh(G::JutulMesh, p)
Construct a coarse mesh from a given JutulMesh
that can be converted to an UnstructuredMesh
instance. The second argument p
should be a partition Vector with one entry per cell in the original grid that assigns that cell to a coarse block. Should be one-indexed and the numbering should be sequential and contain at least one fine cell for each coarse index. This is tested by the function.
Jutul.CompactAutoDiffCache Type
Cache that holds an AD vector/matrix together with their positions.
Jutul.DataDomain Method
DataDomain(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 Type
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.
Jutul.EntityMajorLayout Type
Equations are grouped by entity, listing all equations and derivatives for entity 1 before proceeding to entity 2 etc.
For a test system with primary variables P, S and equations E1, E2 and two cells this will give the following ordering on the diagonal: (∂E1/∂p)₁, (∂E2/∂S)₁, (∂E1/∂p)₂, (∂E2/∂S)₂
Jutul.EquationMajorLayout Type
Equations are stored sequentially in rows, derivatives of same type in columns:
For a test system with primary variables P, S and equations E1, E2 and two cells this will give the following ordering on the diagonal: (∂E1/∂p)₁, (∂E1/∂p)₂, (∂E2/∂S)₁, (∂E2/∂S)₂
Jutul.FaceMap Type
Struct that contains mappings for a set of faces that are made up of nodes and are part of cells.
Jutul.FractionVariables Type
Abstract type for fraction variables (vector variables that sum up to unity over each entity).
By default, these are limited to the [0, 1] range through maximum_value
and minimum_value
default implementations.
Jutul.HelperSimulator Method
HelperSimulator(model::M, T = Float64; state0 = setup_state(model), executor::E = Jutul.default_executor()) where {M, E}
Construct a helper simulator that can be used to compute the residuals and/or accumulation terms for a given type T. Useful for coupling Jutul to other solvers and types of automatic differentiation.
Jutul.IndirectionMap Type
IndirectionMap(vals::Vector{V}, pos::Vector{Int}) where V
Create a indirection map that encodes a variable length dense vector.
pos
is assumed to be a Vector{Int} of length n+1 where n is the number of dense vectors that is encoded. The vals
array holds the entries for vector i in the range pos[i]:(pos[i+1]-1)
for fast lookup. Indexing into the indirection map with index k
will give a view into the values for vector k
.
Jutul.JutulAutoDiffCache Type
An AutoDiffCache is a type that holds both a set of AD values and a map into some global Jacobian.
Jutul.JutulCase Type
JutulCase(model, dt = [1.0], forces = setup_forces(model); state0 = nothing, parameters = nothing, kwarg...)
Set up a structure that holds the complete specification of a simulation case.
Jutul.JutulConfig Type
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!
Jutul.JutulContext Type
Abstract type for the context Jutul should execute in (matrix formats, memory allocation, etc.)
Jutul.JutulDiscretization Type
d = disc(i, Cells())
Ask discretization for entry i
when discretizing some equation on the chosen entity (e.g. Cells
)
Jutul.JutulMatrixLayout Type
Abstract type for matrix layouts. A layout determines how primary variables and equations are ordered in a sparse matrix representation. Note that this is different from the matrix format itself as it concerns the ordering itself: For example, if all equations for a single cell come in sequence, or if a single equation is given for all entities before the next equation is written.
Different layouts does not change the solution of the system, but different linear solvers support different layouts.
Jutul.JutulVariables Type
Abstract type for all variables in Jutul.
A variable is associated with a JutulEntity
through the associated_entity
function. A variable is local to that entity, and cannot depend on other entities. Variables are used by models to define:
primary variables: Sometimes referred to as degrees of freedom, primary unknowns or solution variables
parameters: Static quantities that impact the solution
secondary variables: Can be computed from a combination of other primary and secondary variables and parameters.
Jutul.LUPreconditioner Type
Full LU factorization as preconditioner (intended for smaller subsystems)
Jutul.LUSolver Type
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.
Jutul.MRSTWrapMesh Type
MRSTWrapMesh(G, N = nothing)
Mesh that adapts an exported MRST mesh to the Jutul interface. G
is assumed to be read directly from file using MAT.matread
. The raw exported grid can be found under the data
field.
Jutul.MultiModel Type
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}
.
Jutul.ParallelCSRContext Type
A context that uses a CSR sparse matrix format together with threads. Experimental.
Jutul.SPAI0Preconditioner Type
Sparse Approximate Inverse preconditioner of lowest order – SPAI(0)
Jutul.ScalarVariable Type
Abstract type for scalar variables (one entry per entity, e.g. pressure or temperature in each cell of a model)
Jutul.SimulationModel Method
SimulationModel(domain, system; <kwarg>)
Instantiate a model for a given system
discretized on the domain
.
Jutul.SimulationModel Method
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.Simulator Method
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.
Jutul.SparsityTracingWrapper Method
SparsityTracingWrapper(x::AbstractArray{T, N}) where {T, N}
Create a sparsity tracing wrapper for a numeric array. This wrapped array produces outputs that have the same value as the wrapped type, but contains a SparsityTracing seeded value with seed equal to the column index (if matrix) or linear index (if vector).
Jutul.TwoPointFiniteVolumeGeometry Type
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.
Jutul.UnstructuredMesh Method
UnstructuredMesh(g::CartesianMesh)
Convert CartesianMesh
instance to unstructured grid. Note that the mesh must be 2D and 3D for a 1-to-1 conversion. 1D meshes are implicitly converted to 2D.
Jutul.VectorVariables Type
Abstract type for vector variables (more than one entry per entity, for example saturations or displacements)
Jutul.absolute_increment_limit Method
Absolute allowable change for variable during a nonlinear update.
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.
Jutul.align_to_jacobian! Method
Update an equation so that it knows where to store its derivatives in the Jacobian representation.
Jutul.allocate_array_ad Method
allocate_array_ad(v::AbstractMatrix, ...)
Convert matrix to AD matrix.
Jutul.allocate_array_ad Method
allocate_array_ad(v::AbstractVector, ...)
Convert vector to AD vector.
Jutul.allocate_array_ad Method
allocate_array_ad(n[, m]; <keyword arguments>)
Allocate vector or matrix as AD with optionally provided context and a specified non-zero on the diagonal.
Arguments
n::Integer
: number of entries in vector, or number of rows ifm
is given.m::Integer
: number of rows (optional)
Keyword arguments
npartials = 1
: Number of partials derivatives to allocate for each elementdiag_pos = nothing
: Indices of where to put entities on the diagonal (if any)
Other keyword arguments are passed onto get_ad_entity_scalar
.
Examples:
Allocate a vector with a single partial:
julia> allocate_array_ad(2)
2-element Vector{ForwardDiff.Dual{nothing, Float64, 1}}:
Dual{nothing}(0.0,0.0)
Dual{nothing}(0.0,0.0)
Allocate a vector with two partials, and set the first to one:
julia> allocate_array_ad(2, diag_pos = 1, npartials = 2)
2-element Vector{ForwardDiff.Dual{nothing, Float64, 2}}:
Dual{nothing}(0.0,1.0,0.0)
Dual{nothing}(0.0,1.0,0.0)
Set up a matrix with two partials, where the first column has partials [1, 0] and the second [0, 1]:
julia> allocate_array_ad(2, 2, diag_pos = [1, 2], npartials = 2)
2×2 Matrix{ForwardDiff.Dual{nothing, Float64, 2}}:
Dual{nothing}(0.0,1.0,0.0) Dual{nothing}(0.0,1.0,0.0)
Dual{nothing}(0.0,0.0,1.0) Dual{nothing}(0.0,0.0,1.0)
Jutul.apply_forces! Method
Apply a set of forces to all equations. Equations that don't support a given force will just ignore them, thanks to the power of multiple dispatch.
Jutul.apply_forces_to_equation! Method
apply_forces_to_equation!(diag_part, storage, model, eq, eq_s, force, time)
Update an equation with the effect of a force. The default behavior for any force we do not know about is to assume that the force does not impact this particular equation.
Jutul.as_value Method
Create a mapped array that produces only the values when indexed.
Only useful for AD arrays, otherwise it does nothing.
Jutul.associated_entity Method
The entity a variable is associated with, and can hold partial derivatives with respect to.
Jutul.cell_dims Method
cell_dims(g, pos)::Tuple
Get physical cell dimensions of cell with index pos
for grid g
.
Jutul.cell_index Method
cell_index(g, pos)
Get linear (scalar) index of mesh cell from provided IJK tuple pos
.
Jutul.check_amgcl_availability Method
check_amgcl_availability(; throw = true)
Check if AMGCLWrap extension is available. If throw=true
this wil be an error, otherwise a Boolean indicating if the extension is available will be returned.
Jutul.check_plotting_availability Method
check_plotting_availability(; throw = true, interactive = false)
Check if plotting through at least one Makie
backend is available in the Julia session (after package has been loaded by for example using GLMakie
). The argument throw
can be used to control if this function acts as a programmatic check (throw=false
) there the return value indicates availability, or if an error message is to be printed telling the user how to get plotting working (throw=true
)
An additional check for specifically interactive
plots can also be added.
Jutul.compress_timesteps Function
compress_timesteps(timesteps, forces = nothing; max_step = Inf)
Compress a set of timesteps and forces to the largest possible steps that still covers the same interval and changes forces at exactly the same points in time, while being limited to a maximum size of max_step
.
Jutul.compress_timesteps Method
compress_timesteps(case::JutulCase; max_step = Inf)
Compress time steps for a Jutul case. See compress_timesteps
for the general case.
Jutul.compute_boundary_trans Method
compute_boundary_trans(d::DataDomain, perm)
Compute the boundary half face transmissibilities for perm. The input perm
can either be the symbol of some data defined on Cells()
, a vector of numbers for each cell or a matrix with number of columns equal to the number of cells.
Jutul.compute_face_trans Method
compute_face_trans(g::DataDomain, perm)
Compute face trans for the interior faces. The input perm
can either be the symbol of some data defined on Cells()
, a vector of numbers for each cell or a matrix with number of columns equal to the number of cells.
Jutul.compute_half_face_trans Method
compute_half_face_trans(g::DataDomain, perm)
Compute half-face trans for the interior faces. The input perm
can either be the symbol of some data defined on Cells()
, a vector of numbers for each cell or a matrix with number of columns equal to the number of cells.
Jutul.convergence_criterion Method
convergence_criterion(model, storage, eq, eq_s, r; dt = 1)
Get the convergence criterion values for a given equation. Can be checked against the corresponding tolerances.
Arguments
model
: model that generated the current equation.storage
: global simulator storage.eq::JutulEquation
: equation implementation currently being checkedeq_s
: storage foreq
where values are contained.r
: the local residual part corresponding to this model, as a matrix with column index equaling entity index
Jutul.convert_from_si Method
convert_from_si(value, unit_name::Union{Symbol, String})
Convert value
from SI representation to the unit in unit_symbol
.
Examples
julia> convert_from_si(3600.0, :hour) # Get 3600 s represented as hours
1.0
Jutul.convert_state_ad Function
Convert a state containing variables as arrays of doubles to a state where those arrays contain the same value as Dual types. The dual type is currently taken from ForwardDiff.
Jutul.convert_to_si Method
convert_to_si(value, unit_name::String)
Convert value
to SI representation from value in the unit given by unit_symbol
.
Examples
julia> convert_to_si(1.0, :hour) # Get 1 hour represented as seconds
3600.0
Jutul.data_domain_to_parameters_gradient Method
data_domain_to_parameters_gradient(model, parameter_gradient; dp_dd = missing, config = nothing)
Make a data_domain copy that contains the gradient of some objective with respect to the fields in the data_domain, assuming that the parameters were initialized directly from the data_domain via (setup_parameters
)[@ref].
Jutul.declare_pattern Method
Give out source, target arrays of equal length for a given equation attached to the given model.
Jutul.declare_sparsity Function
Give out I, J arrays of equal length for a given equation attached to the given model.
Jutul.degrees_of_freedom_per_entity Method
Number of independent primary variables / degrees of freedom per computational entity.
Jutul.descalarize_primary_variable! Method
descalarize_primary_variable!(dest_array, model, V, var::Jutul.ScalarVariable, index)
Descalarize a primary variable, overwriting dest_array at entity index
. The AD status of entries in dest_array
will be retained.
Jutul.descalarize_primary_variables! Function
descalarize_primary_variables!(state, model, V, pvars::NamedTuple = (; pairs(model.primary_variables)...), ind = eachindex(V))
Replace valeus in state
by the scalarized values found in V.
Jutul.expand_to_ministeps Method
substates, dt, report_index = expand_to_ministeps(states, reports)
Get states and timesteps at the finest stored resolution. Output lengths depend on if output_substates
option to simulator was enabled.
Jutul.expand_to_ministeps Method
substates, dt, report_index = expand_to_ministeps(result::SimResult)
Get states and timesteps at the finest stored resolution. Output lengths depend on if output_substates
option to simulator was enabled.
Jutul.extra_debug_output! Method
extra_debug_output!(report, storage, model, config, iteration, dt)
Add extra debug output to report during a nonlinear iteration.
Jutul.extract_submesh Method
extract_submesh(g::UnstructuredMesh, cells)
Extract a subgrid for a given mesh and a iterable of cells
to keep.
Jutul.find_enclosing_cell Method
find_enclosing_cell(G::UnstructuredMesh{D}, pt::SVector{D, T},
normals::AbstractVector{SVector{D, T}},
face_centroids::AbstractVector{SVector{D, T}},
boundary_normals::AbstractVector{SVector{D, T}},
boundary_centroids::AbstractVector{SVector{D, T}},
cells = 1:number_of_cells(G)
) where {D, T}
Find enclosing cell of a point. This can be a bit expensive for larger meshes. Recommended to use the more high level find_enclosing_cells
instead.
Jutul.find_enclosing_cells Method
find_enclosing_cells(G, traj; geometry = tpfv_geometry(G), n = 25)
Find the cell indices of cells in the mesh G
that are intersected by a given trajectory traj
. traj
can be either a matrix with equal number of columns as dimensions in G (i.e. three columns for 3D) or a Vector
of SVector
instances with the same length.
The optional argument geometry
is used to define the centroids and normals used in the calculations. You can precompute this if you need to perform many searches. The keyword argument n
can be used to set the number of discretizations in each segment.
use_boundary
is by default set to false
. If set to true, the boundary faces of cells are treated more rigorously when picking exactly what cells are cut by a trajectory, but this requires that the boundary normals are oriented outwards, which is currently not the case for all meshes from downstream packages.
Examples:
# 3D mesh
G = CartesianMesh((4, 4, 5), (100.0, 100.0, 100.0))
trajectory = [
50.0 25.0 1;
55 35.0 25;
65.0 40.0 50.0;
70.0 70.0 90.0
]
cells = Jutul.find_enclosing_cells(G, trajectory)
# Optional plotting, requires Makie:
fig, ax, plt = Jutul.plot_mesh_edges(G)
plot_mesh!(ax, G, cells = cells, alpha = 0.5, transparency = true)
lines!(ax, trajectory, linewidth = 10)
fig
# 2D mesh
G = CartesianMesh((50, 50), (1.0, 2.0))
trajectory = [
0.1 0.1;
0.2 0.4;
0.3 1.2
]
fig, ax, plt = Jutul.plot_mesh_edges(G)
cells = Jutul.find_enclosing_cells(G, trajectory)
# Plotting, needs Makie
fig, ax, plt = Jutul.plot_mesh_edges(G)
plot_mesh!(ax, G, cells = cells, alpha = 0.5, transparency = true)
lines!(ax, trajectory[:, 1], trajectory[:, 2], linewidth = 3)
fig
Jutul.get_1d_interpolator Method
get_1d_interpolator(xs, ys; <keyword arguments>)
Get a 1D interpolator F(x) ≈ y
for a table xs, ys
that by default does constant extrapolation
Arguments
xs
: sorted list of parameter points.ys
: list of function values with equal length toxs
method=LinearInterpolant
: constructor for the interpolation. Defaults toLinearInterpolant
which does simple linear interpolation.cap_endpoints = true
: Add values so that the endpoints are capped (constant extrapolation). Otherwise, the extrapolation will match the method.cap_start = cap_endpoints
: Fine-grained version of cap_endpoints for the start of the interval only (extrapolation forx < xs[1]
)cap_end = cap_endpoints
:Fine-grained version of cap_endpoints for the end of the interval only (extrapolation forx > xs[end]
)
Additional keyword arguments are passed onto the interpolator constructor.
Jutul.get_2d_interpolator Method
get_2d_interpolator(xs, ys, fs; method = BilinearInterpolant, cap_endpoints = true)
For xs
of length nx
and ys
of length ny
generate a 2D interpolation for values given as a nx
by ny
matrix. By default cap_endpoints=true
, and constant extrapolation is used. Fine-grined control over extrapolation can be achieved by setting the keywords arguments cap_x = (cap_low_x, cap_high_x)
and analogously for cap_y
.
Jutul.get_ad_entity_scalar Method
get_ad_entity_scalar(v::Real, npartials, diag_pos = nothing; <keyword_arguments>)
Get scalar with partial derivatives as AD instance.
Arguments
v::Real
: Value of AD variable.npartials
: Number of partial derivatives each AD instance holds.diag_pos
= nothing: Position(s) of where to set 1 as the partial derivative instead of zero.
Keyword arguments
tag = nothing
: Tag for AD instance. Two AD values of the different tag cannot interoperate to avoid perturbation confusion (see ForwardDiff documentation).
Jutul.get_dependencies Method
Get dependencies of variable when viewed as a secondary variable. Normally autogenerated with @jutul_secondary
Jutul.get_diagonal_entries Method
get_diagonal_entries(eq::JutulEquation, eq_s)
Get the diagonal entries of a cache, i.e. the entries where entity type and index equals that of the governing equation.
Note: Be very careful about modifications to this array, as this is a view into the internal AD buffers and it is very easy to create inconsistent Jacobians.
Jutul.get_entity_tag Method
get_entity_tag(basetag, entity)
Combine a base tag (which can be nothing) with a entity to get a tag that captures base tag + entity tag for use with AD initialization.
Jutul.get_entries Method
Get entries of autodiff cache. Entries are AD vectors that hold values and derivatives.
Jutul.get_entries Method
Get the entries of the main autodiff cache for an equation.
Note: This only gets the .equation field's entries.
Jutul.get_mesh_entity_tag Method
get_mesh_entity_tag(met::JutulMesh, entity::JutulEntity, tag_group::Symbol, tag_value = missing; throw = true)
Get the indices tagged for entity
in group tag_group
, optionally for the specific tag_value
. If ismissing(tag_value)
, the Dict containing the tag group will be returned.
Jutul.get_parameters Method
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.
Jutul.get_primary_variable_ordered_entities Method
get_primary_variable_ordered_entities(model::SimulationModel)
Get only the entities where primary variables are present, sorted by their order in the primary variables.
Jutul.get_primary_variables Method
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
.
Jutul.get_secondary_variables Method
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.
Jutul.get_sparse_arguments Method
get_sparse_arguments(storage, model)
Get the [SparsePattern
]@ref for the Jacobian matrix of a given simulator storage and corresponding model.
Jutul.get_tstr Function
get_tstr(dT, lim = 3)
Get formatted time string of dT
given in seconds, limited to lim
number of units.
Jutul.get_variable Method
get_variable(model::SimulationModel, name::Symbol)
Get implementation of variable or parameter with name name
for the model.
Jutul.get_variables Method
get_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
.
Jutul.initialize_extra_state_fields! Method
initialize_extra_state_fields!(state, model::JutulModel)
Add model-dependent changing variables that need to be in state, but are never AD variables themselves (for example status flags).
Jutul.initialize_storage! Method
initialize_storage!(storage, model::JutulModel; initialize_state0 = true)
Initialize the already allocated storage at the beginning of a simulation. Use this to e.g. set up extra stuff in state0 needed for initializing the simulation loop.
Jutul.initialize_variable_value Method
Initializer for the value of non-scalar primary variables
Jutul.interpolation_constant_lookup Function
interpolation_constant_lookup(X, constant_dx = missing)
Generate a lookup table for linear interpolation when dx is evenly spaced.
Note: Setting constant_dx=true
can lead to incorrect interpolations if the data is not evenly spaced.
Jutul.jutul_output_path Function
pth = jutul_output_path(name = missing; subfolder = "jutul", basedir = missing, create = true)
Get path for output. The final path will be found in /basedir/<subfolder/name. If subfolder=missing
, the path will be set to /basedir/name instead. name
will be autogenerated if not provided.
Pass the optional input create = false
to avoid making the directory. To globally set the default output dir, set ENV["JUTUL_OUTPUT_PATH"]``to your desired
basedir``.
Jutul.linear_timestep_selection Function
linear_timestep_selection(x, x0, x1, dt0, dt1)
Produce linear estimate of timestep dt
for some value x
from observed observations. If the observations have the same x
or dt
values, a simple scaling based on the x1
value is used.
Jutul.load_balanced_endpoint Method
load_balanced_endpoint(block_index, nvals, nblocks)
Endpoint for interval block_index
that subdivides nvals
into nblocks
in a load balanced manner. This is done by adding one element to the first set of blocks whenever possible.
Jutul.load_balanced_interval Method
load_balanced_interval(b, n, m)
Create UnitRange for block b ∈ [1, m] for interval of total length n
Jutul.local_ad Method
local_ad(state::T, index::I, ad_tag::∂T) where {T, I<:Integer, ∂T}
Create local_ad for state for index I of AD tag of type ad_tag
Jutul.local_residual_view Method
local_residual_view(r_buf, model, eq, equation_offset)
Get a matrix view of the residual so that, independent of ordering, the column index corresponds to the entity index for the given equation eq
starting at equation_offset
in the global residual buffer r_buf
.
Jutul.merge_step_report_errors Method
merge_step_report_errors(data; fn = max)
Merge step reports errors of the same type using a pair wise reduction (default: max)
Jutul.model_accumulation! Function
model_accumulation!(acc, sim::HelperSimulator, x, dt = 1.0;
forces = setup_forces(sim.model),
update_secondary = true,
kwarg...
)
Compute the accumulation term into Vector acc.
Jutul.model_residual! Function
model_residual!(r, sim, x, x0 = missing, dt = 1.0;
forces = setup_forces(sim.model),
include_accumulation = true,
kwarg...
)
Fill in the model residual into Vector r.
Jutul.model_residual Method
model_residual(sim::HelperSimulator, x, y = missing; kwarg...)
Out of place version of model_residual!
Jutul.number_of_cells Method
number_of_cells(D::Union{DataDomain, DiscretizedDomain})
Get the number of cells in a DataDomain
or DiscretizedDomain
.
Jutul.number_of_degrees_of_freedom Method
Total number of degrees of freedom for a model, over all primary variables and all entities.
Jutul.number_of_entities Method
Get the number of entities (e.g. the number of cells) that the equation is defined on.
Jutul.number_of_entities Method
Number of entities (e.g. Cells, Faces) a variable is defined on. By default, each primary variable exists on all cells of a discretized domain
Jutul.number_of_entities Method
Number of entities for vector stored in state (just the number of elements)
Jutul.number_of_entities Method
Number of entities for matrix stored in state (convention is number of columns)
Jutul.number_of_equations_per_entity Method
Get the number of equations per entity. For example, mass balance of two components will have two equations per grid cell (= entity)
Jutul.number_of_faces Method
number_of_faces(D::Union{DataDomain, DiscretizedDomain})
Get the number of faces in a DataDomain
or DiscretizedDomain
.
Jutul.number_of_half_faces Method
number_of_half_faces(D::Union{DataDomain, DiscretizedDomain})
Get the number of half-faces in a DataDomain
or DiscretizedDomain
.
Jutul.number_of_partials_per_entity Method
number_of_partials_per_entity(model::SimulationModel, entity::JutulEntity)
Get the number of local partial derivatives per entity in a model
for a given JutulEntity
. This is the sum of degrees_of_freedom_per_entity
for all primary variables defined on entity
.
Jutul.number_of_values Function
Total number of values for a model, for a given type of variables over all entities
Jutul.numerical_eltype Method
numerical_eltype(x::AbstractArray{T}) where T
Get the numerical eltype (i.e. the inner type of the element type that could potentially be AD)
Jutul.numerical_type Method
numerical_type(::T) where T
Get the numerical eltype (i.e. the inner type of the element type that could potentially be AD). This function should be overloaded if you have a custom type that wraps a numeric/potentially AD type.
Jutul.parameters_jacobian_wrt_data_domain Method
parameters_jacobian_wrt_data_domain(model; copy = true, config = nothing)
Compute the (sparse) Jacobian of parameters with respect to data_domain values (i.e. floating point values). Optionally, config
can be passed to allow vectorize_variables
to only include a subset of the parameters.
Jutul.partition Function
partition(N::AbstractMatrix, num_coarse, weights = ones(size(N, 2)); partitioner = MetisPartitioner(), groups = nothing, n = maximum(N), group_by_weights = false, buffer_group = true)
Partition based on neighborship (with optional groups kept contigious after partitioning)
Jutul.partition_hypergraph Function
partition_hypergraph(g, n::Int, partitioner = MetisPartitioner(); expand = true)
Partition a hypergraph from setup_partitioner_hypergraph
using a given partitioner. If the optional expand
parameter is set to true the result will be expanded to the full graph (i.e. where groups are not condensed).
Jutul.physical_representation Method
physical_representation(x)
Get the physical representation of an object. The physical representation is usually some kind of mesh or domain that represents a physical domain.
Jutul.physical_representation Method
physical_representation(x::DataDomain)
Get the underlying physical representation (domain or mesh) that is wrapped.
Jutul.physical_representation Method
physical_representation(x::DiscretizedDomain)
Get the underlying physical representation (domain or mesh) that was discretized.
Jutul.physical_representation Method
physical_representation(m::SimulationModel)
Get the underlying physical representation for the model (domain or mesh)
Jutul.pick_next_timestep Method
pick_next_timestep(sel::IterationTimestepSelector, sim, config, dt_prev, dT, forces, reports, current_reports, step_index, new_step)
Pick the next time-step for IterationTimestepSelector
. This function uses the number of iterations from previous timesteps to estimate the relationship between the last and the new time step.
Jutul.plot_mesh! Method
plot_mesh!(ax, mesh)
Mutating version of plot_mesh
that plots into an existing Makie Axis
instance.
Jutul.plot_mesh Method
plot_mesh(mesh)
plot_mesh(mesh;
cells = nothing,
faces = nothing,
boundaryfaces = nothing,
outer = false,
color = :lightblue,
)
Plot a mesh
with uniform colors. Optionally, indices cells
, faces
or boundaryfaces
can be passed to limit the plotting to a specific selection of entities.
Jutul.plot_mesh_edges! Method
plot_mesh_edges!(ax, mesh; kwarg...)
Plot the edges of all cells on the exterior of a mesh into existing Makie Axis
ax
.
Jutul.plot_mesh_edges Method
plot_mesh_edges(mesh; kwarg...)
Plot the edges of all cells on the exterior of a mesh.
Jutul.prepare_step_storage Method
prepare_step_storage(storage, model, ::Missing)
Initialize storage for prepare_step_handler.
Jutul.read_results Method
states, reports = read_results(pth; read_states = true, read_reports = true)
Read results from a given output_path
provded to simulate
or simulator_config
.
Jutul.relative_increment_limit Method
Relative allowable change for variable during a nonlinear update. A variable with value |x| and relative limit 0.2 cannot change more than |x|*0.2.
Jutul.replace_variables! Method
replace_variables!(model, throw = true, varname = vardef, varname2 = vardef2)
Replace one or more variables that already exists in the model (primary, secondary or parameters) with a new definition.
Arguments
model
: instance where variables is to be replacedvarname=vardef::JutulVariables
: replace variable withvarname
byvardef
throw=true
: throw an error if the named variable definition is not found in primary or secondary, otherwise silently return themodel
.
Jutul.scalarize_primary_variable Method
scalarize_primary_variable(model, source_vec, var::Jutul.ScalarVariable, index)
Scalarize a primary variable. For scalars, this means getting the value itself.
Jutul.scalarize_primary_variables Function
scalarize_primary_variables(model, state, pvars = model.primary_variables)
Create a vector where each entry corresponds to a tuple of values that minimally defines the given variables. All variables must belong to the same type of entity. This is checked by this function.
Jutul.scalarize_primary_variables! Method
scalarize_primary_variables!(V::Vector{T}, model, state, pvars::NamedTuple) where T
Scalarize into array. See scalarize_primary_variables
for more details.
Jutul.scalarized_primary_variable_type Method
scalarized_primary_variable_type(model, var::Jutul.ScalarVariable)
Get the type of a scalarized numerical variable (=Float64 for variables that are already represented as scalars)
Jutul.set_default_tolerances Method
set_default_tolerances(model)
Set default tolerances for the nonlinear convergence check of the governing equations.
Jutul.set_parameters! Method
set_parameters!(model, parname = pardef)
Set a parameter with name varname
to the definition vardef
(adding if it does not exist)
Jutul.set_primary_variables! Method
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)
Jutul.set_secondary_variables! Method
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)
Jutul.setup_adjoint_storage Method
setup_adjoint_storage(model; state0 = setup_state(model), parameters = setup_parameters(model))
Set up storage for use with solve_adjoint_sensitivities!
.
Jutul.setup_forces Method
setup_forces(model::JutulModel; force_name = force_value)
Set up forces for a given model. Keyword arguments varies depending on what the model supports.
Jutul.setup_parameter_optimization Method
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: F_o(x) -> evaluate objective dF_o(dFdx, x) -> evaluate gradient of objective, mutating dFdx (may trigger evaluation of F_o) F_and_dF(F, dFdx, x) -> evaluate F and/or dF. Value of nothing will mean that the corresponding entry is skipped.
Jutul.setup_parameters Method
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.
Jutul.setup_partitioner_hypergraph Method
setup_partitioner_hypergraph(N::Matrix{Int};
num_nodes::Int = maximum(N),
num_edges::Int = size(N, 2),
node_weights::Vector{Int} = ones(Int, num_nodes),
edge_weights::Vector{Int} = ones(Int, num_edges),
groups = [Int[]]
)
Set up a hypergraph structure for a given neighborship matrix. N
should be a matrix with two rows, with one pair of cells in each column. Optionally node and edge weights can be provided. If a list of groups are provided, these nodes will be accumulated together in the hypergraph.
Jutul.setup_state! Function
setup_state!(state, model::JutulModel, init_values::AbstractDict = Dict())
Initialize primary variables and other state fields, given initial values as a Dict
Jutul.setup_state Method
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
Jutul.setup_state_and_parameters Method
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)
Jutul.setup_storage! Method
setup_storage!(storage, model::JutulModel; setup_linearized_system = true,
setup_equations = true,
state0 = setup_state(model),
parameters = setup_parameters(model),
tag = nothing,
state0_ad = false,
state_ad = true,
kwarg...)
Allocate storage for a given model. The storage consists of all dynamic quantities used in the simulation. The default implementation allocates properties, equations and linearized system.
Jutul.setup_storage Method
setup_storage(model::JutulModel; kwarg...)
Allocate storage for the model. You should overload setup_storage! if you have a custom definition.
Jutul.si_unit Method
si_unit(u::Union{String, Symbol})
Get the multiplicative SI unit conversion factor for a single unit. The return value is given so that x*si_unit(:name)
will convert x
to the SI representation of the unit with the given name.
Examples
julia> si_unit(:day) # Get days represented as seconds
86400.0
Jutul.si_units Method
u1_val = si_units(u1)
meter = si_units(:meter)
meter, hour = si_units(:meter, :hour)
Get multiplicative SI unit conversion factors for multiple units simultaneously. The return value will be a Tuple
of values, one for each input argument. Each input arguments can be either a String
a Symbol
.
Examples
julia> meter, hour = si_units(:meter, :hour)
(1.0, 3600.0)
Jutul.simulate! Method
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.
Jutul.simulate Method
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 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
.
Jutul.simulate Method
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.simulator_config Method
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.
Jutul.solve_adjoint_sensitivities! Method
solve_adjoint_sensitivities!(∇G, storage, states, state0, timesteps, G; forces = setup_forces(model))
Non-allocating version of solve_adjoint_sensitivities
.
Jutul.solve_adjoint_sensitivities Method
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.
Jutul.solve_numerical_sensitivities Method
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
.
Jutul.solve_timestep! Method
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 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.
Jutul.successful_reports Function
successful_reports(old_reports, current_reports, step_index, n = 1)
Get the n
last successful solve reports from all previous reports (old_reports) and the current ministep set.
Jutul.synchronize Method
Synchronize backend after allocations.
Some backends may require notification that storage has been allocated.
Jutul.tpfv_geometry Function
tpfv_geometry(g)
Generate two-point finite-volume geometry for a given grid, if supported.
See also TwoPointFiniteVolumeGeometry
.
Jutul.transfer Method
Transfer v to the representation expected by a given context.
For the defalt context, the transfer function does nothing. For other context such as the CUDA version, it may convert integers and floats to other types (e.g. Float32) and Arrays to CuArrays.
You will likely have to implement some transfer operators for your own types if you want to simulate with a non-default context.
Jutul.unsafe_reinterpret Method
unsafe_reinterpret(Vt, v, n)
Unsafely reinterpret v as a n length vector of value type Vt
Jutul.update_equations! Function
update_equations!(storage, model, dt = nothing)
Update the governing equations using the current set of primary variables, parameters and secondary variables. Does not fill linearized system.
Jutul.update_equations_and_apply_forces! Method
update_equations_and_apply_forces!(storage, model, dt, forces; time = NaN)
Update the model equations and apply boundary conditions and forces. Does not fill linearized system.
Jutul.update_linearized_system! Function
update_linearized_system!(storage, model::JutulModel; <keyword arguments>)
Update the linearized system with the current set of equations.
Jutul.update_linearized_system_equation! Method
Update a linearized system based on the values and derivatives in the equation.
Jutul.update_parameter_before_step! Method
update_parameter_before_step!(prm_val, prm, storage, model, dt, forces)
Update parameters before time-step. Used for hysteretic parameters.
Jutul.update_secondary_variable! Function
Update a secondary variable. Normally autogenerated with @jutul_secondary
Jutul.update_state_dependents! Method
update_state_dependents!(storage, model, dt, forces; time = NaN, update_secondary = true)
Perform updates of everything that depends on the state: A full linearization for the current primary variables.
This includes properties, governing equations and the linearized system itself.
Jutul.update_values! Method
update_values!(x, dx)
Replace values (for non-Real types, direct assignment)
Jutul.update_values! Method
update_values!(x, dx)
Replace values of x
in-place by y
, leaving x
with the values of y and the partials of x
.
Jutul.values_per_entity Method
Number of values held by a primary variable. Normally this is equal to the number of degrees of freedom, but some special primary variables are most conveniently defined by having N values and N-1 independent variables.
Jutul.variable_scale Method
Define a "typical" numerical value for a variable to scale the linear system entries.
Jutul.write_reports_to_mat_format Function
write_reports_to_mat_format(reports::Vector, path = jutul_output_path(); name = "report", config = missing, verbose = false)
Write the reports to MAT files named "report_1", "report_2", ... to the given path.
Jutul.@jutul_secondary Macro
Designate the function as updating a secondary variable.
A generic evaluator is then defined, together with a function for getting the dependencies of that function upon the state. This is most easily documented with an example. If we define the following function annotated with the macro when updating the array containing the values of MyVarType
realized for some model:
@jutul_secondary function some_fn!(target, var::MyVarType, model, a, b, c, ix)
for i in ix
target[i] = a[i] + b[i] / c[i]
end
end
The purpose of the macro is to translate this into two functions. The first defines for the dependencies of the function with respect to the fields of the state (primary variables, secondary variables and parameters):
function get_dependencies(var::MyVarType, model)
return (:a, :b, :c)
end
The second function defines a generic version that takes in state, and automatically expands the set of dependencies into getfield
calls.
function update_secondary_variable!(array_target, var::MyVarType, model, state, ix)
some_fn!(array_target, var, model, state.a, state.b, state.c, ix)
end
Note that the input names of arguments 4 to end-1 matter, as these will be fetched from state, exactly as written.