sfepy.discrete.variables module

Classes of variables for equations/terms.

class sfepy.discrete.variables.FieldVariable(name, kind, field, order=None, primary_var_name=None, special=None, flags=None, history=None, **kwargs)[source]

A finite element field variable.

field .. field description of variable (borrowed)

apply_ebc(vec, offset=0, force_values=None)[source]

Apply essential (Dirichlet) and periodic boundary conditions to vector vec, starting at offset.

apply_ic(vec, offset=0, force_values=None)[source]

Apply initial conditions conditions to vector vec, starting at offset.

clear_evaluate_cache()[source]

Clear current evaluate cache.

create_output(vec=None, key=None, extend=True, fill_value=None, linearization=None)[source]

Convert the DOF vector to a dictionary of output data usable by Mesh.write().

Parameters:
vec : array, optional

An alternative DOF vector to be used instead of the variable DOF vector.

key : str, optional

The key to be used in the output dictionary instead of the variable name.

extend : bool

Extend the DOF values to cover the whole domain.

fill_value : float or complex

The value used to fill the missing DOF values if extend is True.

linearization : Struct or None

The linearization configuration for higher order approximations.

equation_mapping(bcs, var_di, ts, functions, problem=None, warn=False)[source]

Create the mapping of active DOFs from/to all DOFs.

Sets n_adof.

Returns:
active_bcs : set

The set of boundary conditions active in the current time.

evaluate(mode=’val’, region=None, integral=None, integration=None, step=0, time_derivative=None, is_trace=False, dt=None, bf=None)[source]

Evaluate various quantities related to the variable according to mode in quadrature points defined by integral.

The evaluated data are cached in the variable instance in evaluate_cache attribute.

Parameters:
mode : one of ‘val’, ‘grad’, ‘div’, ‘cauchy_strain’

The evaluation mode.

region : Region instance, optional

The region where the evaluation occurs. If None, the underlying field region is used.

integral : Integral instance, optional

The integral defining quadrature points in which the evaluation occurs. If None, the first order volume integral is created. Must not be None for surface integrations.

integration : ‘volume’, ‘surface’, ‘surface_extra’, or ‘point’

The term integration type. If None, it is derived from integral.

step : int, default 0

The time step (0 means current, -1 previous, …).

time_derivative : None or ‘dt’

If not None, return time derivative of the data, approximated by the backward finite difference.

is_trace : bool, default False

Indicate evaluation of trace of the variable on a boundary region.

dt : float, optional

The time step to be used if derivative is ‘dt’. If None, the dt attribute of the variable is used.

bf : Base function, optional

The base function to be used in ‘val’ mode.

Returns:
out : array

The 4-dimensional array of shape (n_el, n_qp, n_row, n_col) with the requested data, where n_row, n_col depend on mode.

evaluate_at(coors, mode=’val’, strategy=’general’, close_limit=0.1, get_cells_fun=None, cache=None, ret_cells=False, ret_status=False, ret_ref_coors=False, verbose=False)[source]

Evaluate the variable in the given physical coordinates. Convenience wrapper around Field.evaluate_at(), see its docstring for more details.

get_data_shape(integral, integration=’volume’, region_name=None)[source]

Get element data dimensions for given approximation.

Parameters:
integral : Integral instance

The integral describing used numerical quadrature.

integration : ‘volume’, ‘surface’, ‘surface_extra’, ‘point’ or ‘custom’

The term integration type.

region_name : str

The name of the region of the integral.

Returns:
data_shape : 5 ints

The (n_el, n_qp, dim, n_en, n_comp) for volume shape kind, (n_fa, n_qp, dim, n_fn, n_comp) for surface shape kind and (n_nod, 0, 0, 1, n_comp) for point shape kind.

Notes

  • n_el, n_fa = number of elements/facets
  • n_qp = number of quadrature points per element/facet
  • dim = spatial dimension
  • n_en, n_fn = number of element/facet nodes
  • n_comp = number of variable components in a point/node
  • n_nod = number of element nodes
get_dof_conn(dc_type, is_trace=False)[source]

Get active dof connectivity of a variable.

Notes

The primary and dual variables must have the same Region.

get_dof_info(active=False)[source]
get_element_diameters(cells, mode, square=False)[source]

Get diameters of selected elements.

get_field()[source]
get_full(r_vec, r_offset=0, force_value=None, vec=None, offset=0)[source]

Get the full DOF vector satisfying E(P)BCs from a reduced DOF vector.

Notes

The reduced vector starts in r_vec at r_offset. Passing a force_value overrides the EBC values. Optionally, vec argument can be provided to store the full vector (in place) starting at offset.

get_interp_coors(strategy=’interpolation’, interp_term=None)[source]

Get the physical coordinates to interpolate into, based on the strategy used.

get_mapping(region, integral, integration, get_saved=False, return_key=False)[source]

Get the reference element mapping of the underlying field.

See also

sfepy.discrete.fem.fields.Field.get_mapping

get_reduced(vec, offset=0, follow_epbc=False)[source]

Get the reduced DOF vector, with EBC and PBC DOFs removed.

Notes

The full vector starts in vec at offset. If ‘follow_epbc’ is True, values of EPBC master DOFs are not simply thrown away, but added to the corresponding slave DOFs, just like when assembling. For vectors with state (unknown) variables it should be set to False, for assembled vectors it should be set to True.

get_state_in_region(region, reshape=True, step=0)[source]

Get DOFs of the variable in the given region.

Parameters:
region : Region

The selected region.

reshape : bool

If True, reshape the DOF vector to a 2D array with the individual components as columns. Otherwise a 1D DOF array of the form [all DOFs in region node 0, all DOFs in region node 1, …] is returned.

step : int, default 0

The time step (0 means current, -1 previous, …).

Returns:
out : array

The selected DOFs.

has_same_mesh(other)[source]
Returns:
flag : int

The flag can be either ‘different’ (different meshes), ‘deformed’ (slightly deformed same mesh), or ‘same’ (same).

invalidate_evaluate_cache(step=0)[source]

Invalidate variable data in evaluate cache for time step given by step (0 is current, -1 previous, …).

This should be done, for example, prior to every nonlinear solver iteration.

save_as_mesh(filename)[source]

Save the field mesh and the variable values into a file for visualization. Only the vertex values are stored.

set_data_from_qp(data_qp, integral, step=0)[source]

Set DOFs of variable using values in quadrature points corresponding to the given integral.

set_from_mesh_vertices(data)[source]

Set the variable using values at the mesh vertices.

set_from_other(other, strategy=’projection’, close_limit=0.1)[source]

Set the variable using another variable. Undefined values (e.g. outside the other mesh) are set to numpy.nan, or extrapolated.

Parameters:
strategy : ‘projection’ or ‘interpolation’

The strategy to set the values: the L^2 orthogonal projection (not implemented!), or a direct interpolation to the nodes (nodal elements only!)

Notes

If the other variable uses the same field mesh, the coefficients are set directly.

setup_initial_conditions(ics, di, functions, warn=False)[source]

Setup of initial conditions.

time_update(ts, functions)[source]

Store time step, set variable data for variables with the setter function.

class sfepy.discrete.variables.Variable(name, kind, order=None, primary_var_name=None, special=None, flags=None, **kwargs)[source]
advance(ts)[source]

Advance in time the DOF state history. A copy of the DOF vector is made to prevent history modification.

static from_conf(key, conf, fields)[source]
get_dual()[source]

Get the dual variable.

Returns:
var : Variable instance

The primary variable for non-state variables, or the dual variable for state variables.

get_initial_condition()[source]
get_primary()[source]

Get the corresponding primary variable.

Returns:
var : Variable instance

The primary variable, or self for state variables or if primary_var_name is None, or None if no other variables are defined.

get_primary_name()[source]
init_data(step=0)[source]

Initialize the dof vector data of time step step to zeros.

init_history()[source]

Initialize data of variables with history.

is_complex()[source]
is_finite(step=0, derivative=None, dt=None)[source]
is_kind(kind)[source]
is_parameter()[source]
is_real()[source]
is_state()[source]
is_state_or_parameter()[source]
is_virtual()[source]
static reset()[source]
set_constant(val)[source]

Set the variable to a constant value.

set_data(data=None, indx=None, step=0, preserve_caches=False)[source]

Set data (vector of DOF values) of the variable.

Parameters:
data : array

The vector of DOF values.

indx : int, optional

If given, data[indx] is used.

step : int, optional

The time history step, 0 (default) = current.

preserve_caches : bool

If True, do not invalidate evaluate caches of the variable.

time_update(ts, functions)[source]

Implemented in subclasses.

class sfepy.discrete.variables.Variables(variables=None)[source]

Container holding instances of Variable.

advance(ts)[source]
apply_ebc(vec, force_values=None)[source]

Apply essential (Dirichlet) and periodic boundary conditions defined for the state variables to vector vec.

apply_ic(vec, force_values=None)[source]

Apply initial conditions defined for the state variables to vector vec.

check_vector_size(vec, stripped=False)[source]

Check whether the shape of the DOF vector corresponds to the total number of DOFs of the state variables.

Parameters:
vec : array

The vector of DOF values.

stripped : bool

If True, the size of the DOF vector should be reduced, i.e. without DOFs fixed by boundary conditions.

create_state_vector()[source]
create_stripped_state_vector()[source]
equation_mapping(ebcs, epbcs, ts, functions, problem=None, active_only=True)[source]

Create the mapping of active DOFs from/to all DOFs for all state variables.

Parameters:
ebcs : Conditions instance

The essential (Dirichlet) boundary conditions.

epbcs : Conditions instance

The periodic boundary conditions.

ts : TimeStepper instance

The time stepper.

functions : Functions instance

The user functions for boundary conditions.

problem : Problem instance, optional

The problem that can be passed to user functions as a context.

active_only : bool

If True, the active DOF info self.adi uses the reduced (active DOFs only) numbering. Otherwise it is the same as self.di.

Returns:
active_bcs : set

The set of boundary conditions active in the current time.

static from_conf(conf, fields)[source]

This method resets the variable counters for automatic order!

get_dual_names()[source]

Get names of pairs of dual variables.

Returns:
duals : dict

The dual names as virtual name : state name pairs.

get_indx(var_name, stripped=False, allow_dual=False)[source]
get_lcbc_operator()[source]
get_matrix_shape()[source]
get_state_part_view(state, var_name, stripped=False)[source]
get_state_parts(vec=None)[source]

Return parts of a state vector corresponding to individual state variables.

Parameters:
vec : array, optional

The state vector. If not given, then the data stored in the variables are returned instead.

Returns:
out : dict

The dictionary of the state parts.

has_ebc(vec, force_values=None)[source]
has_virtuals()[source]
init_history()[source]
iter_state(ordered=True)[source]

Link state variables with corresponding virtual variables, and assign link to self to each variable instance.

Usually, when solving a PDE in the weak form, each state variable has a corresponding virtual variable.

make_full_vec(svec, force_value=None, vec=None)[source]

Make a full DOF vector satisfying E(P)BCs from a reduced DOF vector.

Parameters:
svec : array

The reduced DOF vector.

force_value : float, optional

Passing a force_value overrides the EBC values.

vec : array, optional

If given, the buffer for storing the result (zeroed).

Returns:
vec : array

The full DOF vector.

set_adof_conns(adof_conns)[source]

Set all active DOF connectivities to self as well as relevant sub-dicts to the individual variables.

set_data(data, step=0, ignore_unknown=False, preserve_caches=False)[source]

Set data (vectors of DOF values) of variables.

Parameters:
data : array

The state vector or dictionary of {variable_name : data vector}.

step : int, optional

The time history step, 0 (default) = current.

ignore_unknown : bool, optional

Ignore unknown variable names if data is a dict.

preserve_caches : bool

If True, do not invalidate evaluate caches of variables.

set_data_from_state(var_names, state, var_names_state)[source]

Set variables with names in var_names from state variables with names in var_names_state using DOF values in the state vector state.

set_state_part(state, part, var_name, stripped=False)[source]
setup_dof_info(make_virtual=False)[source]

Setup global DOF information.

setup_dtype()[source]

Setup data types of state variables - all have to be of the same data type, one of nm.float64 or nm.complex128.

setup_initial_conditions(ics, functions)[source]
setup_lcbc_operators(lcbcs, ts=None, functions=None)[source]

Prepare linear combination BC operator matrix and right-hand side vector.

setup_ordering()[source]

Setup ordering of variables.

state_to_output(vec, fill_value=None, var_info=None, extend=True, linearization=None)[source]

Convert a state vector to a dictionary of output data usable by Mesh.write().

strip_state_vector(vec, follow_epbc=False, svec=None)[source]

Get the reduced DOF vector, with EBC and PBC DOFs removed.

Notes

If ‘follow_epbc’ is True, values of EPBC master dofs are not simply thrown away, but added to the corresponding slave dofs, just like when assembling. For vectors with state (unknown) variables it should be set to False, for assembled vectors it should be set to True.

time_update(ts, functions, verbose=True)[source]
sfepy.discrete.variables.create_adof_conn(eq, conn, dpn, offset)[source]

Given a node connectivity, number of DOFs per node and equation mapping, create the active dof connectivity.

Locally (in a connectivity row), the DOFs are stored DOF-by-DOF (u_0 in all local nodes, u_1 in all local nodes, …).

Globally (in a state vector), the DOFs are stored node-by-node (u_0, u_1, …, u_X in node 0, u_0, u_1, …, u_X in node 1, …).

sfepy.discrete.variables.create_adof_conns(conn_info, var_indx=None, active_only=True, verbose=True)[source]

Create active DOF connectivities for all variables referenced in conn_info.

If a variable has not the equation mapping, a trivial mapping is assumed and connectivity with all DOFs active is created.

DOF connectivity key is a tuple (primary variable name, region name, type, is_trace flag).

Notes

If active_only is False, the DOF connectivities contain all DOFs, with the E(P)BC-constrained ones stored as -1 - <DOF number>, so that the full connectivities can be reconstructed for the matrix graph creation.

sfepy.discrete.variables.expand_basis(basis, dpn)[source]

Expand basis for variables with several components (DOFs per node), in a way compatible with create_adof_conn(), according to dpn (DOF-per-node count).