sfepy.discrete.problem module

class sfepy.discrete.problem.Problem(name, conf=None, functions=None, domain=None, fields=None, equations=None, auto_conf=True, active_only=True)[source]

Problem definition, the top-level class holding all data necessary to solve a problem.

It can be constructed from a ProblemConf instance using Problem.from_conf() or directly from a problem description file using Problem.from_conf_file()

For interactive use, the constructor requires only the equations, nls and ls keyword arguments, see below.

Parameters:
namestr

The problem name.

confProblemConf instance, optional

The ProblemConf describing the problem.

functionsFunctions instance, optional

The user functions for boundary conditions, materials, etc.

domainDomain instance, optional

The solution Domain.

fieldsdict, optional

The dictionary of Field instances.

equationsEquations instance, optional

The Equations to solve. This argument is required when auto_conf is True.

auto_confbool

If True, fields and domain are determined by equations.

active_onlybool

If True, the (tangent) matrices and residual vectors (right-hand sides) contain only active DOFs, see below.

Notes

The Problem is by default created with active_only set to True. Then the (tangent) matrices and residual vectors (right-hand sides) have reduced sizes and contain only the active DOFs, i.e., DOFs not constrained by EBCs or EPBCs.

Setting active_only to False results in full-size vectors and matrices. Then the matrix size non-zeros structure does not depend on the actual E(P)BCs applied. It must be False when using parallel PETSc solvers.

The active 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. However, the negative entries mean that the assembled matrices/residuals have zero values at positions corresponding to constrained DOFs.

The resulting linear system then provides a solution increment, that has to be added to the initial guess used to compute the residual, just like in the Newton iterations. The increment of the constrained DOFs is automatically zero.

When solving with a direct solver, the diagonal entries of a matrix at positions corresponding to constrained DOFs has to be set to ones, so that the matrix is not singular, see sfepy.discrete.evaluate.apply_ebc_to_matrix(), which is called automatically in sfepy.discrete.evaluate.Evaluator.eval_tangent_matrix(). It is not called automatically in Problem.evaluate(). Note that setting the diagonal entries to one might not be necessary with iterative solvers, as the zero matrix rows match the zero residual rows, i.e. if the reduced matrix would be regular, then the right-hand side (the residual) is orthogonal to the kernel of the matrix.

advance(ts=None)[source]
block_solve(state0=None, status=None, save_results=True, step_hook=None, post_process_hook=None, report_nls_status=False, log_nls_status=False, verbose=True)[source]

Call Problem.solve() sequentially for the individual matrix blocks of a block-triangular matrix. It is called by Problem.solve() if the ‘block_solve’ option is set to True.

clear_equations()[source]
copy(name=None)[source]

Make a copy of Problem.

create_evaluable(expression, try_equations=True, auto_init=False, preserve_caches=False, copy_materials=True, integrals=None, ebcs=None, epbcs=None, lcbcs=None, ts=None, functions=None, mode='eval', var_dict=None, strip_variables=True, extra_args=None, active_only=True, eterm_options=None, verbose=True, **kwargs)[source]

Create evaluable object (equations and corresponding variables) from the expression string. Convenience function calling create_evaluable() with defaults provided by the Problem instance self.

The evaluable can be repeatedly evaluated by calling eval_equations(), e.g. for different values of variables.

Parameters:
expressionstr

The expression to evaluate.

try_equationsbool

Try to get variables from self.equations. If this fails, variables can either be provided in var_dict, as keyword arguments, or are created automatically according to the expression.

auto_initbool

Set values of all variables to all zeros.

preserve_cachesbool

If True, do not invalidate evaluate caches of variables.

copy_materialsbool

Work with a copy of self.equations.materials instead of reusing them. Safe but can be slow.

integralsIntegrals instance, optional

The integrals to be used. Automatically created as needed if not given.

ebcsConditions instance, optional

The essential (Dirichlet) boundary conditions for ‘weak’ mode. If not given, self.ebcs are used.

epbcsConditions instance, optional

The periodic boundary conditions for ‘weak’ mode. If not given, self.epbcs are used.

lcbcsConditions instance, optional

The linear combination boundary conditions for ‘weak’ mode. If not given, self.lcbcs are used.

tsTimeStepper instance, optional

The time stepper. If not given, self.ts is used.

functionsFunctions instance, optional

The user functions for boundary conditions, materials etc. If not given, self.functions are used.

modeone of ‘eval’, ‘el_avg’, ‘qp’, ‘weak’

The evaluation mode - ‘weak’ means the finite element assembling, ‘qp’ requests the values in quadrature points, ‘el_avg’ element averages and ‘eval’ means integration over each term region.

var_dictdict, optional

The variables (dictionary of (variable name) : (Variable instance)) to be used in the expression. Use this if the name of a variable conflicts with one of the parameters of this method.

strip_variablesbool

If False, the variables in var_dict or kwargs not present in the expression are added to the actual variables as a context.

extra_argsdict, optional

Extra arguments to be passed to terms in the expression.

active_onlybool

If True, in ‘weak’ mode, the (tangent) matrices and residual vectors (right-hand sides) contain only active DOFs.

eterm_optionsdict, optional

The einsum-based terms evaluation options.

verbosebool

If False, reduce verbosity.

**kwargskeyword arguments

Additional variables can be passed as keyword arguments, see var_dict.

Returns:
equationsEquations instance

The equations that can be evaluated.

variablesVariables instance

The corresponding variables. Set their values and use eval_equations().

Examples

problem is Problem instance.

>>> out = problem.create_evaluable('ev_integrate.i1.Omega(u)')
>>> equations, variables = out

vec is a vector of coefficients compatible with the field of ‘u’ - let’s use all ones.

>>> vec = nm.ones((variables['u'].n_dof,), dtype=nm.float64)
>>> variables['u'].set_data(vec)
>>> vec_qp = eval_equations(equations, variables, mode='qp')

Try another vector:

>>> vec = 3 * nm.ones((variables['u'].n_dof,), dtype=nm.float64)
>>> variables['u'].set_data(vec)
>>> vec_qp = eval_equations(equations, variables, mode='qp')
create_materials(mat_names=None)[source]

Create materials with names in mat_names. Their definitions have to be present in self.conf.materials.

Notes

This method does not change self.equations, so it should not have any side effects.

create_state()[source]
create_subproblem(var_names, known_var_names)[source]

Create a sub-problem with equations containing only terms with the given virtual variables.

Parameters:
var_nameslist

The list of names of virtual variables.

known_var_nameslist

The list of names of (already) known state variables.

Returns:
subpbProblem instance

The sub-problem.

create_variables(var_names=None)[source]

Create variables with names in var_names. Their definitions have to be present in self.conf.variables.

Notes

This method does not change self.equations, so it should not have any side effects.

eval_equations(names=None, preserve_caches=False, mode='eval', dw_mode='vector', term_mode=None, active_only=True, any_dof_conn=False, verbose=True)[source]

Evaluate (some of) the problem’s equations, convenience wrapper of eval_equations().

Parameters:
namesstr or sequence of str, optional

Evaluate only equations of the given name(s).

preserve_cachesbool

If True, do not invalidate evaluate caches of variables.

modeone of ‘eval’, ‘el_avg’, ‘qp’, ‘weak’

The evaluation mode - ‘weak’ means the finite element assembling, ‘qp’ requests the values in quadrature points, ‘el_avg’ element averages and ‘eval’ means integration over each term region.

dw_mode‘vector’ or ‘matrix’

The assembling mode for ‘weak’ evaluation mode.

term_modestr

The term call mode - some terms support different call modes and depending on the call mode different values are returned.

verbosebool

If False, reduce verbosity.

Returns:
outdict or result

The evaluation result. In ‘weak’ mode it is the vector or sparse matrix, depending on dw_mode. Otherwise, it is a dict of results with equation names as keys or a single result for a single equation.

evaluate(expression, try_equations=True, auto_init=False, preserve_caches=False, copy_materials=True, integrals=None, ebcs=None, epbcs=None, lcbcs=None, ts=None, functions=None, mode='eval', dw_mode='vector', term_mode=None, var_dict=None, strip_variables=True, ret_variables=False, any_dof_conn=False, active_only=True, eterm_options=None, verbose=True, extra_args=None, **kwargs)[source]

Evaluate an expression, convenience wrapper of Problem.create_evaluable() and eval_equations().

Parameters:
dw_mode‘vector’ or ‘matrix’

The assembling mode for ‘weak’ evaluation mode.

term_modestr

The term call mode - some terms support different call modes and depending on the call mode different values are returned.

ret_variablesbool

If True, return the variables that were created to evaluate the expression.

otherarguments

See docstrings of Problem.create_evaluable().

Returns:
outarray

The result of the evaluation.

variablesVariables instance

The variables that were created to evaluate the expression. Only provided if ret_variables is True.

static from_conf(conf, init_fields=True, init_equations=True, init_solvers=True)[source]
static from_conf_file(conf_filename, required=None, other=None, init_fields=True, init_equations=True, init_solvers=True)[source]
get_default_ts(t0=None, t1=None, dt=None, n_step=None, step=None)[source]
get_dim(get_sym=False)[source]

Returns mesh dimension, symmetric tensor dimension (if get_sym is True).

get_ebc_indices()[source]

Get indices of E(P)BC-constrained DOFs in the full global state vector.

get_evaluator(reuse=False)[source]

Either create a new Evaluator instance (reuse == False), or return an existing instance, created in a preceding call to Problem.init_solvers().

get_initial_state(vec=None)[source]

Create a zero state and apply initial conditions.

get_integrals(names=None)[source]

Get integrals, initialized from problem configuration if available.

Parameters:
nameslist, optional

If given, only the named integrals are returned.

Returns:
integralsIntegrals instance

The requested integrals.

get_ls()[source]
get_materials()[source]
get_mesh_coors(actual=False)[source]
get_meshes_from_region(reg_names)[source]
get_nls()[source]
get_nls_functions()[source]

Returns functions to be used by a nonlinear solver to evaluate the nonlinear function value (the residual) and its gradient (the tangent matrix) corresponding to the problem equations.

Returns:
funfunction

The function fun(x) for computing the residual.

fun_gradfunction

The function fun_grad(x) for computing the tangent matrix.

iter_hookfunction

The optional (user-defined) function to be called before each nonlinear solver iteration iteration.

get_output_name(suffix=None, extra=None, mode=None)[source]

Return default output file name, based on the output directory, output format, step suffix and mode. If present, the extra string is put just before the output format suffix.

get_restart_filename(ts=None)[source]

If restarts are allowed in problem definition options, return the restart file name, based on the output directory and time step.

get_solver()[source]
get_solver_conf(name)[source]
get_timestepper()[source]
get_tss()[source]
get_tss_functions(update_bcs=True, update_materials=True, save_results=True, step_hook=None, post_process_hook=None)[source]

Get the problem-dependent functions required by the time-stepping solver during the solution process.

Parameters:
update_bcsbool, optional

If True, update the boundary conditions in each prestep_fun call.

update_materialsbool, optional

If True, update the values of material parameters in each prestep_fun call.

save_resultsbool, optional

If True, save the results in each poststep_fun call.

step_hookcallable, optional

The optional user-defined function that is called in each poststep_fun call before saving the results.

post_process_hookcallable, optional

The optional user-defined function that is passed in each poststep_fun to Problem.save_state().

Returns:
init_funcallable

The initialization function called before the actual time-stepping.

prestep_funcallable

The function called in each time (sub-)step prior to the nonlinear solver call.

poststep_funcallable

The function called at the end of each time step.

get_variables(auto_create=False)[source]
init_solvers(status=None, ls_conf=None, nls_conf=None, tsc_conf=None, ts_conf=None, force=False)[source]

Create and initialize solver instances.

Parameters:
statusdict-like, IndexedStruct, optional

The user-supplied object to hold the time-stepping/nonlinear solver convergence statistics.

ls_confStruct, optional

The linear solver options.

nls_confStruct, optional

The nonlinear solver options.

tsc_confStruct, optional

The time step contoller options.

ts_confStruct, optional

The time stepping solver options.

forcebool

If True, re-create the solver instances even if they already exist in self.nls attribute.

init_time(ts)[source]
is_linear()[source]
load_restart(filename, ts=None)[source]

Load the current state and time step from a restart file.

Alternatively, a regular output file in the HDF5 format can be used in place of the restart file. In that case the restart is only approximate, because higher order field DOFs (if any) were stripped out. Files with the adaptive linearization are not supported. Use with caution!

Parameters:
filenamestr

The restart file name.

tsTimeStepper instance, optional

The time stepper. If not given, a default one is created. Otherwise, it is modified in place.

Returns:
variablesVariables instance

The loaded variables.

refine_uniformly(level)[source]

Refine the mesh uniformly level-times.

Notes

This operation resets almost everything (fields, equations, …) - it is roughly equivalent to creating a new Problem instance with the refined mesh.

remove_bcs()[source]

Convenience function to remove boundary conditions.

reset()[source]
save_ebc(filename, ebcs=None, epbcs=None, force=True, default=0.0)[source]

Save essential boundary conditions as state variables.

Parameters:
filenamestr

The output file name.

ebcsConditions instance, optional

The essential (Dirichlet) boundary conditions. If not given, self.conf.ebcs are used.

epbcsConditions instance, optional

The periodic boundary conditions. If not given, self.conf.epbcs are used.

forcebool

If True, sequential nonzero values are forced to individual ebcs so that the conditions are visible even when zero.

defaultfloat

The default constant value of state vector.

save_regions(filename_trunk, region_names=None)[source]

Save regions as meshes.

Parameters:
filename_trunkstr

The output filename without suffix.

region_nameslist, optional

If given, only the listed regions are saved.

save_regions_as_groups(filename_trunk, region_names=None)[source]

Save regions in a single mesh but mark them by using different element/node group numbers.

See Domain.save_regions_as_groups() for more details.

Parameters:
filename_trunkstr

The output filename without suffix.

region_nameslist, optional

If given, only the listed regions are saved.

save_restart(filename, ts=None)[source]

Save the current state and time step to a restart file.

Parameters:
filenamestr

The restart file name.

tsTimeStepper instance, optional

The time stepper. If not given, a default one is created.

Notes

Does not support terms with internal state.

save_state(filename, state=None, out=None, fill_value=None, post_process_hook=None, linearization=None, split_results_by=None, **kwargs)[source]
Parameters:
split_results_byNone, ‘region’, ‘variable’

If ‘region’ or ‘variable’, data of each region/variable are stored in a separate file. If None, it is set to the application option value.

linearizationStruct or None

The linearization configuration for higher order approximations. If its kind is ‘adaptive’, split_results_by is assumed ‘variable’.

select_bcs(ebc_names=None, epbc_names=None, lcbc_names=None, create_matrix=False, any_dof_conn=None)[source]
select_materials(material_names, only_conf=False)[source]
select_variables(variable_names, only_conf=False)[source]
set_bcs(ebcs=None, epbcs=None, lcbcs=None)[source]

Update boundary conditions.

set_conf_solvers(conf_solvers=None, options=None)[source]

Choose which solvers should be used. If solvers are not set in options, use the ones named ls, nls, ts and optionally tsc. If such solver names do not exist, use the first of each required solver kind listed in conf_solvers.

set_default_state(vec=None)[source]

Return variables with an initialized state.

A convenience function that obtains the problem equations’ variables, initializes the state ones with zeros (default) or using vec and then returns the variables.

set_equations(conf_equations=None, user=None, keep_solvers=False, make_virtual=False)[source]

Set equations of the problem using the equations problem description entry.

Fields and Regions have to be already set.

set_equations_instance(equations, keep_solvers=False)[source]

Set equations of the problem to equations.

set_fields(conf_fields=None)[source]
set_ics(ics=None)[source]

Set the initial conditions to use.

set_linear(is_linear)[source]
set_materials(conf_materials=None)[source]

Set definition of materials.

set_mesh_coors(coors, update_fields=False, actual=False, clear_all=True, extra_dofs=False)[source]

Set mesh coordinates.

Parameters:
coorsarray

The new coordinates.

update_fieldsbool

If True, update also coordinates of fields.

actualbool

If True, update the actual configuration coordinates, otherwise the undeformed configuration ones.

set_output_dir(output_dir=None)[source]

Set the directory for output files.

The directory is created if it does not exist.

set_regions(conf_regions=None, conf_materials=None, functions=None, allow_empty=False)[source]
set_solver(solver, status=None)[source]

Set a time-stepping or nonlinear solver to be used in Problem.solve() call.

Parameters:
solverNonlinearSolver or TimeSteppingSolver instance

The nonlinear or time-stepping solver.

statusdict-like, optional

The user-supplied object to hold the solver convergence statistics.

Notes

A copy of the solver is used, and the nonlinear solver functions are set to those returned by Problem.get_nls_functions(), if not set already. If a nonlinear solver is set, a default StationarySolver instance is created automatically as the time-stepping solver. Also sets self.ts attribute.

If self.conf.options.auto_transform_equations is True (the default is False), the problem equations are automatically transformed to a form suitable for the given solver. Implemented for ElastodynamicsBaseTS-based solvers. If it is False, solver.var_names have to be defined.

set_variables(conf_variables=None)[source]

Set definition of variables.

setup_default_output(conf=None, options=None)[source]

Provide default values to Problem.setup_output() from conf.options and options.

setup_hooks(options=None)[source]

Setup various hooks (user-defined functions), as given in options.

Supported hooks:

  • matrix_hook

    • check/modify tangent matrix in each nonlinear solver iteration

  • nls_iter_hook

    • called prior to every iteration of nonlinear solver, if the solver supports that

    • takes the Problem instance (self) as the first argument

setup_output(output_filename_trunk=None, output_dir=None, output_format=None, file_format=None, float_format=None, split_results_by=None, linearization=None)[source]

Sets output options to given values, or uses the defaults for each argument that is None.

solve(state0=None, status=None, force_values=None, var_data=None, update_bcs=True, update_materials=True, save_results=True, step_hook=None, post_process_hook=None, post_process_hook_final=None, report_nls_status=False, log_nls_status=False, verbose=True)[source]

Solve the problem equations by calling the top-level solver.

Before calling this function the top-level solver has to be set, see Problem.set_solver(). Also, the boundary conditions and the initial conditions (for time-dependent problems) has to be set, see Problem.set_bcs(), Problem.set_ics().

Parameters:
state0array, optional

If given, the initial state - then the initial conditions stored in the Problem instance are ignored. By default, the initial state is created and the initial conditions are applied automatically.

statusdict-like, optional

The user-supplied object to hold the solver convergence statistics.

force_valuesdict of floats or float, optional

If given, the supplied values override the values of the essential boundary conditions.

var_datadict, optional

A dictionary of {variable_name : data vector} used to initialize parameter variables.

update_bcsbool, optional

If True, update the boundary conditions in each prestep_fun call. See Problem.get_tss_functions().

update_materialsbool, optional

If True, update the values of material parameters in each prestep_fun call. See Problem.get_tss_functions().

save_resultsbool, optional

If True, save the results in each poststep_fun call. See Problem.get_tss_functions().

step_hookcallable, optional

The optional user-defined function that is called in each poststep_fun call before saving the results. See Problem.get_tss_functions().

post_process_hookcallable, optional

The optional user-defined function that is passed in each poststep_fun to Problem.save_state(). See Problem.get_tss_functions().

post_process_hook_finalcallable, optional

The optional user-defined function that is called after the top-level solver returns.

report_nls_status: bool, optional

If True, print summary non-linear solver info.

log_nls_status: bool, optional

If True, log non-linear solver info into a CSV file.

Returns:
variablesVariables

The variables with the final time step state.

time_update(ts=None, ebcs=None, epbcs=None, lcbcs=None, functions=None, create_matrix=False, is_matrix=True, any_dof_conn=None)[source]
try_presolve(mtx)[source]
update_equations(ts=None, ebcs=None, epbcs=None, lcbcs=None, functions=None, create_matrix=False, is_matrix=True, any_dof_conn=None)[source]

Update equations for current time step.

The tangent matrix graph is automatically recomputed if the set of active essential or periodic boundary conditions changed w.r.t. the previous time step.

Parameters:
tsTimeStepper instance, optional

The time stepper. If not given, self.ts is used.

ebcsConditions instance, optional

The essential (Dirichlet) boundary conditions. If not given, self.ebcs are used.

epbcsConditions instance, optional

The periodic boundary conditions. If not given, self.epbcs are used.

lcbcsConditions instance, optional

The linear combination boundary conditions. If not given, self.lcbcs are used.

functionsFunctions instance, optional

The user functions for boundary conditions, materials, etc. If not given, self.functions are used.

create_matrixbool

If True, force the matrix graph computation.

is_matrixbool

If False, the matrix is not created. Has precedence over create_matrix.

any_dof_connbool or None, default False

If True, all DOF connectivities are used to pre-allocate the matrix graph. If False, only cell region connectivities are used. If None, the value is, if available, taken from conf.options.

update_materials(ts=None, mode='normal', verbose=True)[source]

Update materials used in equations.

Parameters:
tsTimeStepper instance

The time stepper.

mode‘normal’, ‘update’ or ‘force’

The update mode, see Material.time_update().

verbosebool

If False, reduce verbosity.

update_time_stepper(ts)[source]
sfepy.discrete.problem.make_is_save(options)[source]

Given problem options, return a callable that determines whether to save results of a time step.

sfepy.discrete.problem.prepare_matrix(problem, state)[source]

Pre-assemble tangent system matrix.