cashocs#

cashocs is a shape optimization and optimal control software for python.

cashocs is based on the finite element package FEniCS and uses its high-level unified form language UFL to treat general PDE constrained optimization problems, in particular, shape optimization and optimal control problems.

The documentation for cashocs can be found here.

Functions

compute_mesh_quality(mesh[, quality_type, ...])

This computes the mesh quality of a given mesh.

convert(input_file[, output_file, mode, quiet])

Converts the input mesh file to a xdmf mesh file for cashocs to work with.

create_dirichlet_bcs(function_space, value, ...)

Create several Dirichlet boundary conditions at once.

import_mesh(mesh_file[, comm])

Imports a mesh file for use with cashocs / FEniCS.

interpolate_levelset_function_to_cells(...)

Interpolates jumping values to the mesh cells based on the levelset function.

interval_mesh([n, start, end, partitions, comm])

Creates an 1D interval mesh starting at x=0 to x=length.

linear_solve(linear_form, u, bcs[, ...])

Solves a linear problem.

load_config(path)

Loads a config object from a config file.

newton_solve(nonlinear_form, u, bcs[, ...])

Solves a nonlinear problem with Newton's method.

picard_iteration(form_list, u_list, bcs_list)

Solves a system of coupled PDEs via a Picard iteration.

regular_box_mesh([n, start_x, start_y, ...])

Creates a mesh corresponding to a rectangle or cube.

regular_mesh([n, length_x, length_y, ...])

Creates a mesh corresponding to a rectangle or cube.

set_log_level(level)

Determines the log level of cashocs.

Classes

ConstrainedOptimalControlProblem(...[, ...])

An optimal control problem with additional (in-)equality constraints.

ConstrainedShapeOptimizationProblem(...[, ...])

A shape optimization problem with additional (in-)equality constraints.

EqualityConstraint(variable_function, target)

Models an (additional) equality constraint.

Functional()

Base class for all cost functionals.

InequalityConstraint(variable_function[, ...])

Models an (additional) inequality constraint.

IntegralFunctional(form)

A functional which is given by the integral of form.

Interpolator(origin_space, target_space)

Efficient interpolation between two function spaces.

LogLevel()

Stores the various log levels of cashocs.

MinMaxFunctional(integrand[, lower_bound, ...])

Cost functional involving a maximum of 0 and a integral term squared.

OptimalControlProblem(state_forms, bcs_list, ...)

Implements an optimal control problem.

ScalarTrackingFunctional(integrand, ...[, ...])

Tracking cost functional for scalar quantities arising due to integration.

ShapeOptimizationProblem(state_forms, ...[, ...])

A shape optimization problem.

TopologyOptimizationProblem(state_forms, ...)

A topology optimization problem.

class cashocs.ConstrainedOptimalControlProblem(state_forms, bcs_list, cost_functional_form, states, controls, adjoints, constraint_list, config=None, riesz_scalar_products=None, control_constraints=None, initial_guess=None, ksp_options=None, adjoint_ksp_options=None, gradient_ksp_options=None, control_bcs_list=None, preconditioner_forms=None, pre_callback=None, post_callback=None, linear_solver=None, adjoint_linear_solver=None, newton_linearizations=None)[source]#

Bases: ConstrainedOptimizationProblem

An optimal control problem with additional (in-)equality constraints.

Initializes self.

Parameters:
  • state_forms (Union[ufl.Form, List[ufl.Form]]) – The weak form of the state equation (user implemented). Can be either a single UFL form, or a (ordered) list of UFL forms.

  • bcs_list (Union[fenics.DirichletBC, List[fenics.DirichletBC], List[List[fenics.DirichletBC]]]) – The list of fenics.DirichletBC objects describing Dirichlet (essential) boundary conditions. If this is None, then no Dirichlet boundary conditions are imposed.

  • cost_functional_form (Union[List[_typing.CostFunctional], _typing.CostFunctional]) – UFL form of the cost functional. Can also be a list of summands of the cost functional

  • states (Union[fenics.Function, List[fenics.Function]]) – The state variable(s), can either be a fenics.Function, or a list of these.

  • controls (Union[fenics.Function, List[fenics.Function]]) – The control variable(s), can either be a fenics.Function, or a list of these.

  • adjoints (Union[fenics.Function, List[fenics.Function]]) – The adjoint variable(s), can either be a fenics.Function, or a (ordered) list of these.

  • constraint_list (Union[_typing.Constraint, List[_typing.Constraint]]) – (A list of) additional equality and inequality constraints for the problem.

  • config (Optional[io.Config]) – config: The config file for the problem, generated via cashocs.load_config(). Alternatively, this can also be None, in which case the default configurations are used, except for the optimization algorithm. This has then to be specified in the solve method. The default is None.

  • riesz_scalar_products (Optional[Union[ufl.Form, List[ufl.Form]]]) – The scalar products of the control space. Can either be None, a single UFL form, or a (ordered) list of UFL forms. If None, the \(L^2(\Omega)\) product is used (default is None).

  • control_constraints (Optional[List[List[Union[float, fenics.Function]]]]) – Box constraints posed on the control, None means that there are none (default is None). The (inner) lists should contain two elements of the form [u_a, u_b], where u_a is the lower, and u_b the upper bound.

  • initial_guess (Optional[List[fenics.Function]]) – List of functions that act as initial guess for the state variables, should be valid input for fenics.assign(). Defaults to None, which means a zero initial guess.

  • ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to solve the state systems. If this is None, then the direct solver mumps is used (default is None).

  • adjoint_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to solve the adjoint systems. If this is None, then the same options as for the state systems are used (default is None).

  • gradient_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to compute the (shape) gradient. If this is None, either a direct or an iterative method is used (depending on the configuration, section OptimizationRoutine, key gradient_method).

  • control_bcs_list (Optional[Union[fenics.DirichletBC, List[fenics.DirichletBC], List[List[fenics.DirichletBC]]]]) – A list of boundary conditions for the control variables. This is passed analogously to bcs_list. Default is None.

  • preconditioner_forms (Optional[List[ufl.Form]]) – The list of forms for the preconditioner. The default is None, so that the preconditioner matrix is the same as the system matrix.

  • pre_callback (Optional[Callable]) – A function (without arguments) that will be called before each solve of the state system

  • post_callback (Optional[Callable]) – A function (without arguments) that will be called after the computation of the gradient.

  • linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the linear systems arising from the discretized PDE.

  • adjoint_linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the (linear) adjoint system.

  • newton_linearizations (Optional[Union[ufl.Form, List[ufl.Form]]]) – A (list of) UFL forms describing which (alternative) linearizations should be used for the (nonlinear) state equations when solving them (with Newton’s method). The default is None, so that the Jacobian of the supplied state forms is used.

inject_post_callback(function)#

Changes the a-posteriori callback of the OptimizationProblem.

Parameters:

function (Optional[Union[Callable[[], None], Callable[[_typing.OptimizationProblem], None]]]) – A custom function without arguments, which will be called after the computation of the gradient(s)

Return type:

None

inject_pre_callback(function)#

Changes the a-priori callback of the OptimizationProblem.

Parameters:

function (Optional[Union[Callable[[], None], Callable[[_typing.OptimizationProblem], None]]]) – A custom function without arguments, which will be called before each solve of the state system

Return type:

None

inject_pre_post_callback(pre_function, post_function)#

Changes the a-priori (pre) and a-posteriori (post) callbacks of the problem.

Parameters:
  • pre_function (Optional[Union[Callable[[], None], Callable[[_typing.OptimizationProblem], None]]]) – A function without arguments, which is to be called before each solve of the state system

  • post_function (Optional[Union[Callable[[], None], Callable[[_typing.OptimizationProblem], None]]]) – A function without arguments, which is to be called after each computation of the (shape) gradient

Return type:

None

solve(method='Augmented Lagrangian', tol=0.01, max_iter=25, inner_rtol=None, inner_atol=None, constraint_tol=None, mu_0=None, lambda_0=None)#

Solves the constrained optimization problem.

Parameters:
  • method (Literal['Augmented Lagrangian', 'AL', 'Quadratic Penalty', 'QP']) – The solution algorithm, either an augmented Lagrangian method (“Augmented Lagrangian”, “AL”) or quadratic penalty method (“Quadratic Penalty”, “QP”)

  • tol (float) – An overall tolerance to be used in the algorithm. This will set the relative tolerance for the inner optimization problems to tol. Default is 1e-2.

  • max_iter (int) – Maximum number of iterations for the outer solver. Default is 25.

  • inner_rtol (float | None) – Relative tolerance for the inner problem. Default is None,

  • used. (so that inner_rtol = tol is)

  • inner_atol (float | None) – Absolute tolerance for the inner problem. Default is None, so that inner_atol = tol/10 is used.

  • constraint_tol (float | None) – The tolerance for the constraint violation, which is desired. If this is None (default), then this is specified as tol/10.

  • mu_0 (float | None) – Initial value of the penalty parameter. Default is None, which means that mu_0 = 1 is used.

  • lambda_0 (List[float] | None) – Initial guess for the Lagrange multipliers. Default is None, which corresponds to a zero guess.

Return type:

None

total_constraint_violation()#

Computes the total constraint violation.

Returns:

The 2-norm of the total constraint violation.

Return type:

float

solver: solvers.AugmentedLagrangianMethod | solvers.QuadraticPenaltyMethod#
class cashocs.ConstrainedShapeOptimizationProblem(state_forms, bcs_list, cost_functional_form, states, adjoints, boundaries, constraint_list, config=None, shape_scalar_product=None, initial_guess=None, ksp_options=None, adjoint_ksp_options=None, gradient_ksp_options=None, preconditioner_forms=None, pre_callback=None, post_callback=None, linear_solver=None, adjoint_linear_solver=None, newton_linearizations=None)[source]#

Bases: ConstrainedOptimizationProblem

A shape optimization problem with additional (in-)equality constraints.

Initializes self.

Parameters:
  • state_forms (Union[ufl.Form, List[ufl.Form]]) – The weak form of the state equation (user implemented). Can be either a single UFL form, or a (ordered) list of UFL forms.

  • bcs_list (Union[fenics.DirichletBC, List[fenics.DirichletBC], List[List[fenics.DirichletBC]], None]) – The list of fenics.DirichletBC objects describing Dirichlet (essential) boundary conditions. If this is None, then no Dirichlet boundary conditions are imposed.

  • cost_functional_form (Union[List[_typing.CostFunctional], _typing.CostFunctional]) – UFL form of the cost functional. Can also be a list of summands of the cost functional

  • states (Union[fenics.Function, List[fenics.Function]]) – The state variable(s), can either be a fenics.Function, or a list of these.

  • adjoints (Union[fenics.Function, List[fenics.Function]]) – The adjoint variable(s), can either be a fenics.Function, or a (ordered) list of these.

  • boundaries (fenics.MeshFunction) – A fenics.MeshFunction that indicates the boundary markers.

  • constraint_list (Union[_typing.Constraint, List[_typing.Constraint]]) – (A list of) additional equality and inequality constraints for the problem.

  • config (Optional[io.Config]) – config: The config file for the problem, generated via cashocs.load_config(). Alternatively, this can also be None, in which case the default configurations are used, except for the optimization algorithm. This has then to be specified in the solve method. The default is None.

  • shape_scalar_product (Optional[ufl.Form]) – The bilinear form for computing the shape gradient (or gradient deformation). This has to use fenics.TrialFunction and fenics.TestFunction objects to define the weak form, which have to be in a fenics.VectorFunctionSpace of continuous, linear Lagrange finite elements. Moreover, this form is required to be symmetric.

  • initial_guess (Optional[List[fenics.Function]]) – List of functions that act as initial guess for the state variables, should be valid input for fenics.assign(). Defaults to None, which means a zero initial guess.

  • ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to solve the state systems. If this is None, then the direct solver mumps is used (default is None).

  • adjoint_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to solve the adjoint systems. If this is None, then the same options as for the state systems are used (default is None).

  • gradient_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to compute the (shape) gradient. If this is None, either a direct or an iterative method is used (depending on the configuration, section OptimizationRoutine, key gradient_method).

  • preconditioner_forms (Optional[List[ufl.Form]]) – The list of forms for the preconditioner. The default is None, so that the preconditioner matrix is the same as the system matrix.

  • pre_callback (Optional[Callable]) – A function (without arguments) that will be called before each solve of the state system

  • post_callback (Optional[Callable]) – A function (without arguments) that will be called after the computation of the gradient.

  • linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the linear systems arising from the discretized PDE.

  • adjoint_linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the (linear) adjoint system.

  • newton_linearizations (Optional[Union[ufl.Form, List[ufl.Form]]]) – A (list of) UFL forms describing which (alternative) linearizations should be used for the (nonlinear) state equations when solving them (with Newton’s method). The default is None, so that the Jacobian of the supplied state forms is used.

inject_post_callback(function)#

Changes the a-posteriori callback of the OptimizationProblem.

Parameters:

function (Optional[Union[Callable[[], None], Callable[[_typing.OptimizationProblem], None]]]) – A custom function without arguments, which will be called after the computation of the gradient(s)

Return type:

None

inject_pre_callback(function)#

Changes the a-priori callback of the OptimizationProblem.

Parameters:

function (Optional[Union[Callable[[], None], Callable[[_typing.OptimizationProblem], None]]]) – A custom function without arguments, which will be called before each solve of the state system

Return type:

None

inject_pre_post_callback(pre_function, post_function)#

Changes the a-priori (pre) and a-posteriori (post) callbacks of the problem.

Parameters:
  • pre_function (Optional[Union[Callable[[], None], Callable[[_typing.OptimizationProblem], None]]]) – A function without arguments, which is to be called before each solve of the state system

  • post_function (Optional[Union[Callable[[], None], Callable[[_typing.OptimizationProblem], None]]]) – A function without arguments, which is to be called after each computation of the (shape) gradient

Return type:

None

solve(method='Augmented Lagrangian', tol=0.01, max_iter=25, inner_rtol=None, inner_atol=None, constraint_tol=None, mu_0=None, lambda_0=None)#

Solves the constrained optimization problem.

Parameters:
  • method (Literal['Augmented Lagrangian', 'AL', 'Quadratic Penalty', 'QP']) – The solution algorithm, either an augmented Lagrangian method (“Augmented Lagrangian”, “AL”) or quadratic penalty method (“Quadratic Penalty”, “QP”)

  • tol (float) – An overall tolerance to be used in the algorithm. This will set the relative tolerance for the inner optimization problems to tol. Default is 1e-2.

  • max_iter (int) – Maximum number of iterations for the outer solver. Default is 25.

  • inner_rtol (float | None) – Relative tolerance for the inner problem. Default is None,

  • used. (so that inner_rtol = tol is)

  • inner_atol (float | None) – Absolute tolerance for the inner problem. Default is None, so that inner_atol = tol/10 is used.

  • constraint_tol (float | None) – The tolerance for the constraint violation, which is desired. If this is None (default), then this is specified as tol/10.

  • mu_0 (float | None) – Initial value of the penalty parameter. Default is None, which means that mu_0 = 1 is used.

  • lambda_0 (List[float] | None) – Initial guess for the Lagrange multipliers. Default is None, which corresponds to a zero guess.

Return type:

None

total_constraint_violation()#

Computes the total constraint violation.

Returns:

The 2-norm of the total constraint violation.

Return type:

float

solver: solvers.AugmentedLagrangianMethod | solvers.QuadraticPenaltyMethod#
class cashocs.EqualityConstraint(variable_function, target, measure=None)[source]#

Bases: Constraint

Models an (additional) equality constraint.

Initializes self.

Parameters:
  • variable_function (Union[ufl.Form, ufl_expr.Expr]) – Either a UFL Form (when we have a scalar / integral constraint) or an ufl expression (when we have a pointwise constraint), which models the part that is to be constrained.

  • target (float) – The target (rhs) of the equality constraint.

  • measure (Optional[fenics.Measure]) – A measure indicating where a pointwise constraint should be satisfied.

constraint_violation()[source]#

Computes the constraint violation for the problem.

Returns:

The computed violation

Return type:

float

lower_bound: fenics.Function | float#
min_max_term: MinMaxFunctional#
multiplier: fenics.Function#
target: float#
upper_bound: fenics.Function | float#
class cashocs.Functional[source]#

Bases: ABC

Base class for all cost functionals.

Initialize the functional.

abstract coefficients()[source]#

Computes the ufl coefficients which are used in the functional.

Returns:

The set of used coefficients.

Return type:

Tuple[fenics.Function]

abstract derivative(argument, direction)[source]#

Computes the derivative of the functional w.r.t. argument towards direction.

Parameters:
  • argument (ufl_legacy.core.expr.Expr) – The argument, w.r.t. which the functional is differentiated

  • direction (ufl_legacy.core.expr.Expr) – The direction into which the derivative is computed

Returns:

A form of the resulting derivative

Return type:

ufl_legacy.Form

abstract evaluate()[source]#

Evaluates the functional.

Returns:

The current value of the functional.

Return type:

float

abstract scale(scaling_factor)[source]#

Scales the functional by a scalar.

Parameters:

scaling_factor (float | int) – The scaling factor used to scale the functional

Return type:

None

abstract update()[source]#

Updates the functional after solving the state equation.

Return type:

None

class cashocs.InequalityConstraint(variable_function, lower_bound=None, upper_bound=None, measure=None)[source]#

Bases: Constraint

Models an (additional) inequality constraint.

Initializes self.

Parameters:
  • variable_function (Union[ufl.Form, ufl_expr.Expr]) – Either a UFL Form (when we have a scalar / integral constraint) or an ufl expression (when we have a pointwise constraint), which models the part that is to be constrained

  • lower_bound (fenics.Function | float) – The lower bound for the inequality constraint

  • upper_bound (fenics.Function | float) – The upper bound for the inequality constraint

  • measure (Optional[fenics.Measure]) – A measure indicating where a pointwise constraint should be satisfied.

constraint_violation()[source]#

Computes the constraint violation for the problem.

Returns:

The computed violation

Return type:

float

lower_bound: fenics.Function | float#
min_max_term: MinMaxFunctional#
multiplier: fenics.Function#
target: float#
upper_bound: fenics.Function | float#
class cashocs.IntegralFunctional(form)[source]#

Bases: Functional

A functional which is given by the integral of form.

Initializes self.

Parameters:

form (ufl.Form) – The form of the integrand, which is to be calculated for evaluating the functional.

coefficients()[source]#

Computes the ufl coefficients which are used in the functional.

Returns:

The set of used coefficients.

Return type:

Tuple[fenics.Function]

derivative(argument, direction)[source]#

Computes the derivative of the functional w.r.t. argument towards direction.

Parameters:
  • argument (ufl_legacy.core.expr.Expr) – The argument, w.r.t. which the functional is differentiated

  • direction (ufl_legacy.core.expr.Expr) – The direction into which the derivative is computed

Returns:

A form of the resulting derivative

Return type:

ufl_legacy.Form

evaluate()[source]#

Evaluates the functional.

Returns:

The current value of the functional.

Return type:

float

scale(scaling_factor)[source]#

Scales the functional by a scalar.

Parameters:

scaling_factor (float | int) – The scaling factor used to scale the functional

Return type:

None

update()[source]#

Updates the functional after solving the state equation.

Return type:

None

class cashocs.Interpolator(origin_space, target_space)[source]#

Bases: object

Efficient interpolation between two function spaces.

This is very useful, if multiple interpolations have to be carried out between the same spaces, which is made significantly faster by computing the corresponding matrix. The function spaces can even be defined on different meshes.

Notes

This class only works properly for continuous Lagrange elements and constant, discontinuous Lagrange elements.

Examples

Here, we consider interpolating from CG1 elements to CG2 elements

import fenics
import cashocs

mesh, _, _, _, _, _ = cashocs.regular_mesh(25)
V1 = fenics.FunctionSpace(mesh, 'CG', 1)
V2 = fenics.FunctionSpace(mesh, 'CG', 2)

expr = fenics.Expression('sin(2*pi*x[0])', degree=1)
u = fenics.interpolate(expr, V1)

interp = cashocs._utils.Interpolator(V1, V2)
interp.interpolate(u)

Initializes self.

Parameters:
  • origin_space (fenics.FunctionSpace) – The function space whose objects shall be interpolated.

  • target_space (fenics.FunctionSpace) – The space into which they shall be interpolated.

interpolate(u)[source]#

Interpolates function to target space.

The function has to belong to the origin space, i.e., the first argument of __init__, and it is interpolated to the destination space, i.e., the second argument of __init__. There is no need to call set_allow_extrapolation on the function (this is done automatically due to the method).

Parameters:

u (fenics.Function) – The function that shall be interpolated.

Returns:

The result of the interpolation.

Return type:

fenics.Function

class cashocs.LogLevel[source]#

Bases: object

Stores the various log levels of cashocs.

CRITICAL = 50#
DEBUG = 10#
ERROR = 40#
INFO = 20#
WARNING = 30#
class cashocs.MinMaxFunctional(integrand, lower_bound=None, upper_bound=None, mu=1.0, lambd=0.0)[source]#

Bases: Functional

Cost functional involving a maximum of 0 and a integral term squared.

Initializes self.

Parameters:
  • integrand (ufl.Form)

  • lower_bound (Optional[Union[float, int]])

  • upper_bound (Optional[Union[float, int]])

  • mu (Union[float, int])

  • lambd (Union[float, int])

coefficients()[source]#

Computes the ufl coefficients which are used in the functional.

Returns:

The set of used coefficients.

Return type:

Tuple[fenics.Function]

derivative(argument, direction)[source]#

Computes the derivative of the functional w.r.t. argument towards direction.

Parameters:
  • argument (ufl_legacy.core.expr.Expr) – The argument, w.r.t. which the functional is differentiated

  • direction (ufl_legacy.core.expr.Expr) – The direction into which the derivative is computed

Returns:

A form of the resulting derivative

Return type:

ufl_legacy.Form

evaluate()[source]#

Evaluates the functional.

Returns:

The current value of the functional.

Return type:

float

scale(scaling_factor)[source]#

Scales the functional by a scalar.

Parameters:

scaling_factor (float | int) – The scaling factor used to scale the functional

Return type:

None

update()[source]#

Updates the functional after solving the state equation.

Return type:

None

class cashocs.OptimalControlProblem(state_forms, bcs_list, cost_functional_form, states, controls, adjoints, config=None, riesz_scalar_products=None, control_constraints=None, initial_guess=None, ksp_options=None, adjoint_ksp_options=None, gradient_ksp_options=None, desired_weights=None, control_bcs_list=None, preconditioner_forms=None, pre_callback=None, post_callback=None, linear_solver=None, adjoint_linear_solver=None, newton_linearizations=None)[source]#

Bases: OptimizationProblem

Implements an optimal control problem.

This class is used to define an optimal control problem, and also to solve it subsequently. For a detailed documentation, see the examples in the tutorial. For easier input, when considering single (state or control) variables, these do not have to be wrapped into a list. Note, that in the case of multiple variables these have to be grouped into ordered lists, where state_forms, bcs_list, states, adjoints have to have the same order (i.e. [y1, y2] and [p1, p2], where p1 is the adjoint of y1 and so on).

Initializes self.

Parameters:
  • state_forms (Union[List[ufl.Form], ufl.Form]) – The weak form of the state equation (user implemented). Can be either a single UFL form, or a (ordered) list of UFL forms.

  • bcs_list (Union[List[List[fenics.DirichletBC]], List[fenics.DirichletBC], fenics.DirichletBC]) – The list of fenics.DirichletBC objects describing Dirichlet (essential) boundary conditions. If this is None, then no Dirichlet boundary conditions are imposed.

  • cost_functional_form (Union[List[_typing.CostFunctional], _typing.CostFunctional]) – UFL form of the cost functional. Can also be a list of summands of the cost functional

  • states (Union[List[fenics.Function], fenics.Function]) – The state variable(s), can either be a fenics.Function, or a list of these.

  • controls (Union[List[fenics.Function], fenics.Function]) – The control variable(s), can either be a fenics.Function, or a list of these.

  • adjoints (Union[List[fenics.Function], fenics.Function]) – The adjoint variable(s), can either be a fenics.Function, or a (ordered) list of these.

  • config (Optional[io.Config]) – The config file for the problem, generated via cashocs.load_config(). Alternatively, this can also be None, in which case the default configurations are used, except for the optimization algorithm. This has then to be specified in the solve method. The default is None.

  • riesz_scalar_products (Optional[Union[List[ufl.Form], ufl.Form]]) – The scalar products of the control space. Can either be None, a single UFL form, or a (ordered) list of UFL forms. If None, the \(L^2(\Omega)\) product is used (default is None).

  • control_constraints (Optional[List[List[Union[float, fenics.Function]]]]) – Box constraints posed on the control, None means that there are none (default is None). The (inner) lists should contain two elements of the form [u_a, u_b], where u_a is the lower, and u_b the upper bound.

  • initial_guess (Optional[List[fenics.Function]]) – List of functions that act as initial guess for the state variables, should be valid input for fenics.assign(). Defaults to None, which means a zero initial guess.

  • ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to solve the state systems. If this is None, then the direct solver mumps is used (default is None).

  • adjoint_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to solve the adjoint systems. If this is None, then the same options as for the state systems are used (default is None).

  • gradient_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to compute the (shape) gradient. If this is None, either a direct or an iterative method is used (depending on the configuration, section OptimizationRoutine, key gradient_method).

  • desired_weights (Optional[List[float]]) – A list of values for scaling the cost functional terms. If this is supplied, the cost functional has to be given as list of summands. The individual terms are then scaled, so that term i has the magnitude of desired_weights[i] for the initial iteration. In case that desired_weights is None, no scaling is performed. Default is None.

  • control_bcs_list (Optional[Union[List[List[fenics.DirichletBC]], List[fenics.DirichletBC], fenics.DirichletBC]]) – A list of boundary conditions for the control variables. This is passed analogously to bcs_list. Default is None.

  • preconditioner_forms (Optional[Union[List[ufl.Form], ufl.Form]]) – The list of forms for the preconditioner. The default is None, so that the preconditioner matrix is the same as the system matrix.

  • pre_callback (Optional[Callable]) – A function (without arguments) that will be called before each solve of the state system

  • post_callback (Optional[Callable]) – A function (without arguments) that will be called after the computation of the gradient.

  • linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the linear systems arising from the discretized PDE.

  • adjoint_linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the (linear) adjoint system.

  • newton_linearizations (Optional[Union[ufl.Form, List[ufl.Form]]]) – A (list of) UFL forms describing which (alternative) linearizations should be used for the (nonlinear) state equations when solving them (with Newton’s method). The default is None, so that the Jacobian of the supplied state forms is used.

Examples

Examples how to use this class can be found in the tutorial.

compute_adjoint_variables()#

Solves the adjoint system.

This can be used for debugging purposes and solver validation. Updates / overwrites the user input for the adjoint variables. The solution of the corresponding state system needed to determine the adjoints is carried out automatically.

Return type:

None

compute_gradient()[source]#

Solves the Riesz problem to determine the gradient.

This can be used for debugging, or code validation. The necessary solutions of the state and adjoint systems are carried out automatically.

Returns:

A list consisting of the (components) of the gradient.

Return type:

List[fenics.Function]

compute_state_variables()#

Solves the state system.

This can be used for debugging purposes and to validate the solver. Updates and overwrites the user input for the state variables.

Return type:

None

gradient_test(u=None, h=None, rng=None)[source]#

Performs a Taylor test to verify correctness of the computed gradient.

Parameters:
  • u (List[fenics.Function] | None) – The point, at which the gradient shall be verified. If this is None, then the current controls of the optimization problem are used. Default is None.

  • h (List[fenics.Function] | None) – The direction(s) for the directional (Gâteaux) derivative. If this is None, one random direction is chosen. Default is None.

  • rng (RandomState | None) – A numpy random state for calculating a random direction.

Returns:

The convergence order from the Taylor test. If this is (approximately) 2 or larger, everything works as expected.

Return type:

float

initialize_solve_parameters(algorithm=None, rtol=None, atol=None, max_iter=None)#

Solves the optimization problem by the method specified in the config file.

Parameters:
  • algorithm (str | None) – Selects the optimization algorithm. Valid choices are 'gradient_descent' or 'gd' for a gradient descent method, 'conjugate_gradient', 'nonlinear_cg', 'ncg' or 'cg' for nonlinear conjugate gradient methods, and 'lbfgs' or 'bfgs' for limited memory BFGS methods. This overwrites the value specified in the config file. If this is None, then the value in the config file is used. Default is None. In addition, for optimal control problems, one can use 'newton' for a truncated Newton method.

  • rtol (float | None) – The relative tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • atol (float | None) – The absolute tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • max_iter (int | None) – The maximum number of iterations the optimization algorithm can carry out before it is terminated. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

Return type:

None

Notes

If either rtol or atol are specified as arguments to the .solve call, the termination criterion changes to:

  • a purely relative one (if only rtol is specified), i.e.,

\[|| \nabla J(u_k) || \leq \texttt{rtol} || \nabla J(u_0) ||.\]
  • a purely absolute one (if only atol is specified), i.e.,

\[|| \nabla J(u_K) || \leq \texttt{atol}.\]
  • a combined one if both rtol and atol are specified, i.e.,

\[|| \nabla J(u_k) || \leq \texttt{atol} + \texttt{rtol} || \nabla J(u_0) ||\]
inject_post_callback(function)#

Changes the a-posteriori callback of the OptimizationProblem to function.

Parameters:

function (Callable | None) – A custom function without arguments, which will be called after the computation of the gradient(s)

Return type:

None

inject_pre_callback(function)#

Changes the a-priori callback of the OptimizationProblem to function.

Parameters:

function (Callable | None) – A custom function without arguments, which will be called before each solve of the state system

Return type:

None

inject_pre_post_callback(pre_function, post_function)#

Changes the a-priori (pre) and a-posteriori (post) callbacks.

Parameters:
  • pre_function (Callable | None) – A function without arguments, which is to be called before each solve of the state system

  • post_function (Callable | None) – A function without arguments, which is to be called after each computation of the (shape) gradient

Return type:

None

shift_cost_functional(shift=0.0)#

Shifts the cost functional by a constant.

Parameters:

shift (float) – The constant, by which the cost functional is shifted.

Return type:

None

solve(algorithm=None, rtol=None, atol=None, max_iter=None)[source]#

Solves the optimization problem by the method specified in the config file.

Updates / overwrites states, controls, and adjoints according to the optimization method, i.e., the user-input fenics.Function() objects.

Parameters:
  • algorithm (str | None) – Selects the optimization algorithm. Valid choices are 'gradient_descent' or 'gd' for a gradient descent method, 'conjugate_gradient', 'nonlinear_cg', 'ncg' or 'cg' for nonlinear conjugate gradient methods, 'lbfgs' or 'bfgs' for limited memory BFGS methods, and 'newton' for a truncated Newton method. This overwrites the value specified in the config file. If this is None, then the value in the config file is used. Default is None.

  • rtol (float | None) – The relative tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • atol (float | None) – The absolute tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • max_iter (int | None) – The maximum number of iterations the optimization algorithm can carry out before it is terminated. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

Return type:

None

Notes

If either rtol or atol are specified as arguments to the .solve call, the termination criterion changes to: - a purely relative one (if only rtol is specified), i.e.,

\[|| \nabla J(u_k) || \leq \texttt{rtol} || \nabla J(u_0) ||.\]
  • a purely absolute one (if only atol is specified), i.e.,

\[|| \nabla J(u_K) || \leq \texttt{atol}.\]
  • a combined one if both rtol and atol are specified, i.e.,

\[|| \nabla J(u_k) || \leq \texttt{atol} + \texttt{rtol} || \nabla J(u_0) ||\]
supply_adjoint_forms(adjoint_forms, adjoint_bcs_list)#

Overwrites the computed weak forms of the adjoint system.

This allows the user to specify their own weak forms of the problems and to use cashocs merely as a solver for solving the optimization problems.

Parameters:
  • adjoint_forms (ufl_legacy.Form | List[ufl_legacy.Form]) – The UFL forms of the adjoint system(s).

  • adjoint_bcs_list (fenics.DirichletBC | List[fenics.DirichletBC] | List[List[fenics.DirichletBC]]) – The list of Dirichlet boundary conditions for the adjoint system(s).

Return type:

None

supply_custom_forms(derivatives, adjoint_forms, adjoint_bcs_list)[source]#

Overrides both adjoint system and derivatives with user input.

This allows the user to specify both the derivatives of the reduced cost functional and the corresponding adjoint system, and allows them to use cashocs as a solver.

Parameters:
  • derivatives (ufl_legacy.Form | List[ufl_legacy.Form]) – The derivatives of the reduced (!) cost functional w.r.t. the control variables.

  • adjoint_forms (ufl_legacy.Form | List[ufl_legacy.Form]) – The UFL forms of the adjoint system(s).

  • adjoint_bcs_list (fenics.DirichletBC | List[fenics.DirichletBC] | List[List[fenics.DirichletBC]]) – The list of Dirichlet boundary conditions for the adjoint system(s).

Return type:

None

supply_derivatives(derivatives)[source]#

Overwrites the derivatives of the reduced cost functional w.r.t. controls.

This allows users to implement their own derivatives and use cashocs as a solver library only.

Parameters:
  • derivatives (ufl_legacy.Form | List[ufl_legacy.Form]) – The derivatives of the reduced (!) cost functional w.r.t.

  • variables. (the control)

Return type:

None

adjoint_problem: _pde_problems.AdjointProblem#
config: io.Config#
cost_functional_list: List[_typing.CostFunctional]#
form_handler: _forms.ControlFormHandler#
gradient: List[fenics.Function]#
gradient_problem: _pde_problems.ControlGradientProblem#
hessian_problem: _pde_problems.HessianProblem#
initial_guess: List[fenics.Function] | None#
optimization_variable_abstractions: ova.OptimizationVariableAbstractions#
output_manager: io.OutputManager#
reduced_cost_functional: cost_functional.ReducedCostFunctional#
solver: optimization_algorithms.OptimizationAlgorithm | topology_optimization_algorithm.TopologyOptimizationAlgorithm#
state_problem: _pde_problems.StateProblem#
uses_custom_scalar_product: bool = False#
class cashocs.ScalarTrackingFunctional(integrand, tracking_goal, weight=1.0)[source]#

Bases: Functional

Tracking cost functional for scalar quantities arising due to integration.

Initializes self.

Parameters:
  • integrand (ufl.Form) – The integrand of the functional

  • tracking_goal (Union[float, int, ctypes.c_float, ctypes.c_double]) – A real number, which the integral of the integrand should track. Note, that when a ctypes object is passed, the float is assumed to be mutable and the tracking_goal is updated every iteration.

  • weight (Union[float, int]) – A real number which gives the scaling factor for this functional

coefficients()[source]#

Computes the ufl coefficients which are used in the functional.

Returns:

The set of used coefficients.

Return type:

Tuple[fenics.Function]

derivative(argument, direction)[source]#

Computes the derivative of the functional w.r.t. argument towards direction.

Parameters:
  • argument (ufl_legacy.core.expr.Expr) – The argument, w.r.t. which the functional is differentiated

  • direction (ufl_legacy.core.expr.Expr) – The direction into which the derivative is computed

Returns:

A form of the resulting derivative

Return type:

ufl_legacy.Form

evaluate()[source]#

Evaluates the functional.

Returns:

The current value of the functional.

Return type:

float

scale(scaling_factor)[source]#

Scales the functional by a scalar.

Parameters:

scaling_factor (float | int) – The scaling factor used to scale the functional

Return type:

None

update()[source]#

Updates the functional after solving the state equation.

Return type:

None

class cashocs.ShapeOptimizationProblem(state_forms, bcs_list, cost_functional_form, states, adjoints, boundaries, config=None, shape_scalar_product=None, initial_guess=None, ksp_options=None, adjoint_ksp_options=None, gradient_ksp_options=None, desired_weights=None, temp_dict=None, initial_function_values=None, preconditioner_forms=None, pre_callback=None, post_callback=None, linear_solver=None, adjoint_linear_solver=None, newton_linearizations=None)[source]#

Bases: OptimizationProblem

A shape optimization problem.

This class is used to define a shape optimization problem, and to solve it subsequently. For a detailed documentation, we refer to the tutorial. For easier input, when consider single (state or control) variables, these do not have to be wrapped into a list. Note, that in the case of multiple variables these have to be grouped into ordered lists, where state_forms, bcs_list, states, adjoints have to have the same order (i.e. [y1, y2] and [p1, p2], where p1 is the adjoint of y1 and so on).

Initializes self.

Parameters:
  • state_forms (Union[List[ufl.Form], ufl.Form]) – The weak form of the state equation (user implemented). Can be either a single UFL form, or a (ordered) list of UFL forms.

  • bcs_list (Union[List[List[fenics.DirichletBC]], List[fenics.DirichletBC], fenics.DirichletBC]) – The list of fenics.DirichletBC objects describing Dirichlet (essential) boundary conditions. If this is None, then no Dirichlet boundary conditions are imposed.

  • cost_functional_form (Union[List[_typing.CostFunctional], _typing.CostFunctional]) – UFL form of the cost functional. Can also be a list of summands of the cost functional

  • states (Union[List[fenics.Function], fenics.Function]) – The state variable(s), can either be a fenics.Function, or a list of these.

  • adjoints (Union[List[fenics.Function], fenics.Function]) – The adjoint variable(s), can either be a fenics.Function, or a (ordered) list of these.

  • boundaries (fenics.MeshFunction) – A fenics.MeshFunction that indicates the boundary markers.

  • config (Optional[io.Config]) – The config file for the problem, generated via cashocs.load_config(). Alternatively, this can also be None, in which case the default configurations are used, except for the optimization algorithm. This has then to be specified in the solve method. The default is None.

  • shape_scalar_product (Optional[ufl.Form]) – The bilinear form for computing the shape gradient (or gradient deformation). This has to use fenics.TrialFunction and fenics.TestFunction objects to define the weak form, which have to be in a fenics.VectorFunctionSpace of continuous, linear Lagrange finite elements. Moreover, this form is required to be symmetric.

  • initial_guess (Optional[List[fenics.Function]]) – List of functions that act as initial guess for the state variables, should be valid input for fenics.assign(). Defaults to None, which means a zero initial guess.

  • ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to solve the state systems. If this is None, then the direct solver mumps is used (default is None).

  • adjoint_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to solve the adjoint systems. If this is None, then the same options as for the state systems are used (default is None).

  • gradient_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to compute the (shape) gradient. If this is None, either a direct or an iterative method is used (depending on the configuration, section OptimizationRoutine, key gradient_method).

  • desired_weights (Optional[List[float]]) – A list of values for scaling the cost functional terms. If this is supplied, the cost functional has to be given as list of summands. The individual terms are then scaled, so that term i has the magnitude of desired_weights[i] for the initial iteration. In case that desired_weights is None, no scaling is performed. Default is None.

  • temp_dict (Optional[Dict]) – This is a private parameter of the class, required for remeshing. This parameter must not be set by the user and should be ignored. Using this parameter may result in unintended side effects.

  • initial_function_values (Optional[List[float]]) – This is a privatve parameter of the class, required for remeshing. This parameter must not be set by the user and should be ignored. Using this parameter may result in unintended side effects.

  • preconditioner_forms (Optional[Union[List[ufl.Form], ufl.Form]]) – The list of forms for the preconditioner. The default is None, so that the preconditioner matrix is the same as the system matrix.

  • pre_callback (Optional[Callable]) – A function (without arguments) that will be called before each solve of the state system

  • post_callback (Optional[Callable]) – A function (without arguments) that will be called after the computation of the gradient.

  • linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the linear systems arising from the discretized PDE.

  • adjoint_linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the (linear) adjoint system.

  • newton_linearizations (Optional[Union[ufl.Form, List[ufl.Form]]]) – A (list of) UFL forms describing which (alternative) linearizations should be used for the (nonlinear) state equations when solving them (with Newton’s method). The default is None, so that the Jacobian of the supplied state forms is used.

compute_adjoint_variables()#

Solves the adjoint system.

This can be used for debugging purposes and solver validation. Updates / overwrites the user input for the adjoint variables. The solution of the corresponding state system needed to determine the adjoints is carried out automatically.

Return type:

None

compute_shape_gradient()[source]#

Solves the Riesz problem to determine the shape gradient.

This can be used for debugging, or code validation. The necessary solutions of the state and adjoint systems are carried out automatically.

Returns:

A list containing the shape gradient.

Return type:

List[fenics.Function]

compute_state_variables()#

Solves the state system.

This can be used for debugging purposes and to validate the solver. Updates and overwrites the user input for the state variables.

Return type:

None

get_vector_field()[source]#

Returns the TestFunction for defining shape derivatives.

Returns:

The TestFunction object.

Return type:

dolfin.function.argument.Argument

gradient_test(h=None, rng=None)[source]#

Performs a Taylor test to verify that the computed shape gradient is correct.

Parameters:
  • h (fenics.Function | None) – The direction used to compute the directional derivative. If this is None, then a random direction is used (default is None).

  • rng (RandomState | None) – A numpy random state for calculating a random direction

Returns:

The convergence order from the Taylor test. If this is (approximately) 2 or larger, everything works as expected.

Return type:

float

initialize_solve_parameters(algorithm=None, rtol=None, atol=None, max_iter=None)#

Solves the optimization problem by the method specified in the config file.

Parameters:
  • algorithm (str | None) – Selects the optimization algorithm. Valid choices are 'gradient_descent' or 'gd' for a gradient descent method, 'conjugate_gradient', 'nonlinear_cg', 'ncg' or 'cg' for nonlinear conjugate gradient methods, and 'lbfgs' or 'bfgs' for limited memory BFGS methods. This overwrites the value specified in the config file. If this is None, then the value in the config file is used. Default is None. In addition, for optimal control problems, one can use 'newton' for a truncated Newton method.

  • rtol (float | None) – The relative tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • atol (float | None) – The absolute tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • max_iter (int | None) – The maximum number of iterations the optimization algorithm can carry out before it is terminated. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

Return type:

None

Notes

If either rtol or atol are specified as arguments to the .solve call, the termination criterion changes to:

  • a purely relative one (if only rtol is specified), i.e.,

\[|| \nabla J(u_k) || \leq \texttt{rtol} || \nabla J(u_0) ||.\]
  • a purely absolute one (if only atol is specified), i.e.,

\[|| \nabla J(u_K) || \leq \texttt{atol}.\]
  • a combined one if both rtol and atol are specified, i.e.,

\[|| \nabla J(u_k) || \leq \texttt{atol} + \texttt{rtol} || \nabla J(u_0) ||\]
inject_post_callback(function)#

Changes the a-posteriori callback of the OptimizationProblem to function.

Parameters:

function (Callable | None) – A custom function without arguments, which will be called after the computation of the gradient(s)

Return type:

None

inject_pre_callback(function)#

Changes the a-priori callback of the OptimizationProblem to function.

Parameters:

function (Callable | None) – A custom function without arguments, which will be called before each solve of the state system

Return type:

None

inject_pre_post_callback(pre_function, post_function)#

Changes the a-priori (pre) and a-posteriori (post) callbacks.

Parameters:
  • pre_function (Callable | None) – A function without arguments, which is to be called before each solve of the state system

  • post_function (Callable | None) – A function without arguments, which is to be called after each computation of the (shape) gradient

Return type:

None

shift_cost_functional(shift=0.0)#

Shifts the cost functional by a constant.

Parameters:

shift (float) – The constant, by which the cost functional is shifted.

Return type:

None

solve(algorithm=None, rtol=None, atol=None, max_iter=None)[source]#

Solves the optimization problem by the method specified in the config file.

Parameters:
  • algorithm (str | None) – Selects the optimization algorithm. Valid choices are 'gradient_descent' or 'gd' for a gradient descent method, 'conjugate_gradient', 'nonlinear_cg', 'ncg' or 'cg' for nonlinear conjugate gradient methods, and 'lbfgs' or 'bfgs' for limited memory BFGS methods. This overwrites the value specified in the config file. If this is None, then the value in the config file is used. Default is None.

  • rtol (float | None) – The relative tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • atol (float | None) – The absolute tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • max_iter (int | None) – The maximum number of iterations the optimization algorithm can carry out before it is terminated. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

Return type:

None

Notes

If either rtol or atol are specified as arguments to the .solve call, the termination criterion changes to: - a purely relative one (if only rtol is specified), i.e.,

\[|| \nabla J(u_k) || \leq \texttt{rtol} || \nabla J(u_0) ||.\]
  • a purely absolute one (if only atol is specified), i.e.,

\[|| \nabla J(u_K) || \leq \texttt{atol}.\]
  • a combined one if both rtol and atol are specified, i.e.,

\[|| \nabla J(u_k) || \leq \texttt{atol} + \texttt{rtol} || \nabla J(u_0) ||\]
supply_adjoint_forms(adjoint_forms, adjoint_bcs_list)#

Overwrites the computed weak forms of the adjoint system.

This allows the user to specify their own weak forms of the problems and to use cashocs merely as a solver for solving the optimization problems.

Parameters:
  • adjoint_forms (ufl_legacy.Form | List[ufl_legacy.Form]) – The UFL forms of the adjoint system(s).

  • adjoint_bcs_list (fenics.DirichletBC | List[fenics.DirichletBC] | List[List[fenics.DirichletBC]]) – The list of Dirichlet boundary conditions for the adjoint system(s).

Return type:

None

supply_custom_forms(shape_derivative, adjoint_forms, adjoint_bcs_list)[source]#

Overrides both adjoint system and shape derivative with user input.

This allows the user to specify both the shape_derivative of the reduced cost functional and the corresponding adjoint system, and allows them to use cashocs as a solver.

Parameters:
  • shape_derivative (ufl_legacy.Form) – The shape derivative of the reduced (!) cost functional.

  • adjoint_forms (ufl_legacy.Form | List[ufl_legacy.Form]) – The UFL forms of the adjoint system(s).

  • adjoint_bcs_list (fenics.DirichletBC | List[fenics.DirichletBC] | List[List[fenics.DirichletBC]]) – The list of Dirichlet boundary conditions for the adjoint system(s).

Return type:

None

supply_shape_derivative(shape_derivative)[source]#

Overrides the shape derivative of the reduced cost functional.

This allows users to implement their own shape derivative and use cashocs as a solver library only.

Parameters:

shape_derivative (ufl_legacy.Form) – The shape_derivative of the reduced cost functional.

Return type:

None

adjoint_problem: _pde_problems.AdjointProblem#
config: io.Config#
cost_functional_list: List[_typing.CostFunctional]#
form_handler: _forms.ShapeFormHandler#
gradient: List[fenics.Function]#
gradient_problem: _typing.GradientProblem#
hessian_problem: _pde_problems.HessianProblem#
initial_guess: List[fenics.Function] | None#
optimization_variable_abstractions: ova.OptimizationVariableAbstractions#
output_manager: io.OutputManager#
reduced_cost_functional: cost_functional.ReducedCostFunctional#
solver: optimization_algorithms.OptimizationAlgorithm | topology_optimization_algorithm.TopologyOptimizationAlgorithm#
state_problem: _pde_problems.StateProblem#
uses_custom_scalar_product: bool = False#
class cashocs.TopologyOptimizationProblem(state_forms, bcs_list, cost_functional_form, states, adjoints, levelset_function, topological_derivative_neg, topological_derivative_pos, update_levelset, config=None, riesz_scalar_products=None, initial_guess=None, ksp_options=None, adjoint_ksp_options=None, gradient_ksp_options=None, desired_weights=None, preconditioner_forms=None, pre_callback=None, post_callback=None, linear_solver=None, adjoint_linear_solver=None, newton_linearizations=None)[source]#

Bases: OptimizationProblem

A topology optimization problem.

This class is used to define a topology optimization problem, and to solve it subsequently. For a detailed documentation, we refer to the tutorial. For easier input, when considering single (state or control) variables, these do not have to be wrapped into a list. Note, that in the case of multiple variables these have to be grouped into ordered lists, where state_forms, bcs_list, states, adjoints have to have the same order (i.e. [y1, y2] and [p1, p2], where p1 is the adjoint of y1 and so on).

Initializes the topology optimization problem.

Parameters:
  • state_forms (list[ufl.Form] | ufl.Form) – The weak form of the state equation (user implemented). Can be either a single UFL form, or a (ordered) list of UFL forms.

  • bcs_list (list[list[fenics.DirichletBC]] | list[fenics.DirichletBC] | fenics.DirichletBC) – The list of fenics.DirichletBC objects describing Dirichlet (essential) boundary conditions. If this is None, then no Dirichlet boundary conditions are imposed.

  • cost_functional_form (list[_typing.CostFunctional] | _typing.CostFunctional) – UFL form of the cost functional. Can also be a list of summands of the cost functional

  • states (list[fenics.Function] | fenics.Function) – The state variable(s), can either be a fenics.Function, or a list of these.

  • adjoints (list[fenics.Function] | fenics.Function) – The adjoint variable(s), can either be a fenics.Function, or a (ordered) list of these.

  • levelset_function (fenics.Function) – A fenics.Function which represents the levelset function.

  • topological_derivative_neg (fenics.Function | ufl.Form) – The topological derivative inside the domain, where the levelset function is negative.

  • topological_derivative_pos (fenics.Function | ufl.Form) – The topological derivative inside the domain, where the levelset function is positive.

  • update_levelset (Callable) – A python function (without arguments) which is called to update the coefficients etc. when the levelset function is changed.

  • config (io.Config | None) – The config file for the problem, generated via cashocs.create_config(). Alternatively, this can also be None, in which case the default configurations are used, except for the optimization algorithm. This has then to be specified in the solve method. The default is None.

  • riesz_scalar_products (list[ufl.Form] | ufl.Form | None) – The scalar products of the control space. Can either be None or a single UFL form. If it is None, the \(L^2(\Omega)\) product is used (default is None).

  • initial_guess (list[fenics.Function] | None) – List of functions that act as initial guess for the state variables, should be valid input for fenics.assign(). Defaults to None, which means a zero initial guess.

  • ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of strings corresponding to command line options for PETSc, used to solve the state systems. If this is None, then the direct solver mumps is used (default is None).

  • adjoint_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of strings corresponding to command line options for PETSc, used to solve the adjoint systems. If this is None, then the same options as for the state systems are used (default is None).

  • gradient_ksp_options (Optional[Union[_typing.KspOption, List[_typing.KspOption]]]) – A list of dicts corresponding to command line options for PETSc, used to compute the (shape) gradient. If this is None, either a direct or an iterative method is used (depending on the configuration, section OptimizationRoutine, key gradient_method).

  • desired_weights (list[float] | None) – A list of values for scaling the cost functional terms. If this is supplied, the cost functional has to be given as list of summands. The individual terms are then scaled, so that term i has the magnitude of desired_weights[i] for the initial iteration. In case that desired_weights is None, no scaling is performed. Default is None.

  • preconditioner_forms (Optional[Union[List[ufl.Form], ufl.Form]]) – The list of forms for the preconditioner. The default is None, so that the preconditioner matrix is the same as the system matrix.

  • pre_callback (Optional[Callable]) – A function (without arguments) that will be called before each solve of the state system

  • post_callback (Optional[Callable]) – A function (without arguments) that will be called after the computation of the gradient.

  • linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the linear systems arising from the discretized PDE.

  • adjoint_linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the (linear) adjoint system.

  • newton_linearizations (Optional[Union[ufl.Form, List[ufl.Form]]]) – A (list of) UFL forms describing which (alternative) linearizations should be used for the (nonlinear) state equations when solving them (with Newton’s method). The default is None, so that the Jacobian of the supplied state forms is used.

compute_adjoint_variables()#

Solves the adjoint system.

This can be used for debugging purposes and solver validation. Updates / overwrites the user input for the adjoint variables. The solution of the corresponding state system needed to determine the adjoints is carried out automatically.

Return type:

None

compute_state_variables()#

Solves the state system.

This can be used for debugging purposes and to validate the solver. Updates and overwrites the user input for the state variables.

Return type:

None

gradient_test()[source]#

Gradient test for topology optimization - not implemented.

Return type:

float

initialize_solve_parameters(algorithm=None, rtol=None, atol=None, max_iter=None)#

Solves the optimization problem by the method specified in the config file.

Parameters:
  • algorithm (str | None) – Selects the optimization algorithm. Valid choices are 'gradient_descent' or 'gd' for a gradient descent method, 'conjugate_gradient', 'nonlinear_cg', 'ncg' or 'cg' for nonlinear conjugate gradient methods, and 'lbfgs' or 'bfgs' for limited memory BFGS methods. This overwrites the value specified in the config file. If this is None, then the value in the config file is used. Default is None. In addition, for optimal control problems, one can use 'newton' for a truncated Newton method.

  • rtol (float | None) – The relative tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • atol (float | None) – The absolute tolerance used for the termination criterion. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

  • max_iter (int | None) – The maximum number of iterations the optimization algorithm can carry out before it is terminated. Overwrites the value specified in the config file. If this is None, the value from the config file is taken. Default is None.

Return type:

None

Notes

If either rtol or atol are specified as arguments to the .solve call, the termination criterion changes to:

  • a purely relative one (if only rtol is specified), i.e.,

\[|| \nabla J(u_k) || \leq \texttt{rtol} || \nabla J(u_0) ||.\]
  • a purely absolute one (if only atol is specified), i.e.,

\[|| \nabla J(u_K) || \leq \texttt{atol}.\]
  • a combined one if both rtol and atol are specified, i.e.,

\[|| \nabla J(u_k) || \leq \texttt{atol} + \texttt{rtol} || \nabla J(u_0) ||\]
inject_post_callback(function)#

Changes the a-posteriori callback of the OptimizationProblem to function.

Parameters:

function (Callable | None) – A custom function without arguments, which will be called after the computation of the gradient(s)

Return type:

None

inject_pre_callback(function)#

Changes the a-priori callback of the OptimizationProblem to function.

Parameters:

function (Callable | None) – A custom function without arguments, which will be called before each solve of the state system

Return type:

None

inject_pre_post_callback(pre_function, post_function)#

Changes the a-priori (pre) and a-posteriori (post) callbacks.

Parameters:
  • pre_function (Callable | None) – A function without arguments, which is to be called before each solve of the state system

  • post_function (Callable | None) – A function without arguments, which is to be called after each computation of the (shape) gradient

Return type:

None

plot_shape()[source]#

Visualize the current shape in a plot.

Return type:

None

shift_cost_functional(shift=0.0)#

Shifts the cost functional by a constant.

Parameters:

shift (float) – The constant, by which the cost functional is shifted.

Return type:

None

solve(algorithm=None, rtol=None, atol=None, max_iter=None, angle_tol=None)[source]#

Solves the optimization problem.

Parameters:
  • algorithm (str | None) – Selects the optimization algorithm. Valid choices are 'gradient_descent' or 'gd' for a gradient descent method, 'conjugate_gradient', 'nonlinear_cg', 'ncg' or 'cg' for nonlinear conjugate gradient methods, 'lbfgs' or 'bfgs' for limited memory BFGS methods, 'sphere_combination' for Euler’s method on the spehere, and 'convex_combination' for a convex combination approach.

  • rtol (float | None) – The relative tolerance used for the termination criterion (i.e. the norm of the projected topological gradient). If this is None, then the value provided in the config file is used. Default is None.

  • atol (float | None) – The absolute tolerance for the termination criterion (i.e. the norm of the projected topological gradient). If this is None, then the value provided in the config file is used. Default is None.

  • max_iter (int | None) – The maximum number of iterations the optimization algorithm can carry out before it is terminated. If this is None, then the value provided in the config file is used. Default is None.

  • angle_tol (float | None) – The absolute tolerance for the angle between topological derivative and levelset function. If this is None, then the value provided in the config file is used. Default is None.

Return type:

None

supply_adjoint_forms(adjoint_forms, adjoint_bcs_list)#

Overwrites the computed weak forms of the adjoint system.

This allows the user to specify their own weak forms of the problems and to use cashocs merely as a solver for solving the optimization problems.

Parameters:
  • adjoint_forms (ufl_legacy.Form | List[ufl_legacy.Form]) – The UFL forms of the adjoint system(s).

  • adjoint_bcs_list (fenics.DirichletBC | List[fenics.DirichletBC] | List[List[fenics.DirichletBC]]) – The list of Dirichlet boundary conditions for the adjoint system(s).

Return type:

None

adjoint_problem: _pde_problems.AdjointProblem#
config: io.Config#
cost_functional_list: List[_typing.CostFunctional]#
form_handler: _forms.ControlFormHandler#
gradient: List[fenics.Function]#
gradient_problem: _pde_problems.ControlGradientProblem#
hessian_problem: _pde_problems.HessianProblem#
initial_guess: List[fenics.Function] | None#
optimization_variable_abstractions: ova.OptimizationVariableAbstractions#
output_manager: io.OutputManager#
reduced_cost_functional: cost_functional.ReducedCostFunctional#
solver: topology_optimization_algorithm.TopologyOptimizationAlgorithm#
state_problem: _pde_problems.StateProblem#
uses_custom_scalar_product: bool = False#
cashocs.compute_mesh_quality(mesh, quality_type='min', quality_measure='skewness')[source]#

This computes the mesh quality of a given mesh.

Parameters:
  • mesh (fenics.Mesh) – The mesh whose quality shall be computed.

  • quality_type (str) – The type of measurement for the mesh quality, either minimum quality or average quality over all mesh cells, default is ‘min’.

  • quality_measure (str) – The type of quality measure which is used to compute the quality measure, default is ‘skewness’

Returns:

The quality of the mesh, in the interval \([0,1]\), where 0 is the worst, and 1 the best possible quality.

Return type:

float

cashocs.convert(input_file, output_file=None, mode='physical', quiet=False)[source]#

Converts the input mesh file to a xdmf mesh file for cashocs to work with.

Parameters:
  • input_file (str) – A gmsh .msh file.

  • output_file (str | None) – The name of the output .xdmf file or None. If this is None, then a file name.msh will be converted to name.xdmf, i.e., the name of the input file stays the same

  • quiet (bool) – A boolean flag which silences the output.

  • mode (str) – The mode which is used to define the subdomains and boundaries. Should be one of ‘physical’ (the default), ‘geometrical’, or ‘none’.

Return type:

None

cashocs.create_dirichlet_bcs(function_space, value, boundaries, idcs, **kwargs)[source]#

Create several Dirichlet boundary conditions at once.

Wraps multiple Dirichlet boundary conditions into a list, in case they have the same value but are to be defined for multiple boundaries with different markers. Particularly useful for defining homogeneous boundary conditions.

Parameters:
  • function_space (fenics.FunctionSpace) – The function space onto which the BCs should be imposed on.

  • value (fenics.Constant | fenics.Expression | fenics.Function | float | Tuple[float]) – The value of the boundary condition. Has to be compatible with the function_space, so that it could also be used as fenics.DirichletBC(function_space, value, ...).

  • boundaries (fenics.MeshFunction) – The fenics.MeshFunction object representing the boundaries.

  • idcs (List[int | str] | int | str) – A list of indices / boundary markers that determine the boundaries onto which the Dirichlet boundary conditions should be applied to. Can also be a single entry for a single boundary. If your mesh file is named, then you can also use the names of the boundaries to define the boundary conditions.

  • **kwargs (Any) – Keyword arguments for fenics.DirichletBC

Returns:

A list of DirichletBC objects that represent the boundary conditions.

Return type:

List[fenics.DirichletBC]

Examples

Generate homogeneous Dirichlet boundary conditions for all 4 sides of the unit square

import fenics
import cashocs

mesh, _, _, _, _, _ = cashocs.regular_mesh(25)
V = fenics.FunctionSpace(mesh, 'CG', 1)
bcs = cashocs.create_dirichlet_bcs(V, fenics.Constant(0), boundaries,
    [1,2,3,4])
cashocs.import_mesh(mesh_file, comm=None)[source]#

Imports a mesh file for use with cashocs / FEniCS.

This function imports a mesh file. The mesh file can either be a Gmsh mesh file or an xdmf mesh file that was generated by GMSH and converted to .xdmf with the function cashocs.convert(). If there are Physical quantities specified in the GMSH file, these are imported to the subdomains and boundaries output of this function and can also be directly accessed via the measures, e.g., with dx(1), ds(1), etc.

Parameters:
  • mesh_file (str) – The location of the mesh file in .xdmf or .msh file format.

  • comm (Optional[MPI.Comm]) – MPI communicator that is to be used for creating the mesh.

Returns:

A tuple (mesh, subdomains, boundaries, dx, ds, dS), where mesh is the imported FEM mesh, subdomains is a mesh function for the subdomains, boundaries is a mesh function for the boundaries, dx is a volume measure, ds is a surface measure, and dS is a measure for the interior facets.

Return type:

_typing.MeshTuple

Notes

In case the boundaries in the Gmsh .msh file are not only marked with numbers (as physical groups), but also with names (i.e. strings), these strings can be used with the integration measures dx and ds returned by this method. E.g., if one specified the following in a 2D Gmsh .geo file

Physical Surface("domain", 1) = {i,j,k};

where i,j,k are representative for some integers, then this can be used in the measure dx (as we are 2D) as follows. The command

dx(1)

is completely equivalent to

dx("domain")

and both can be used interchangeably.

cashocs.interpolate_levelset_function_to_cells(levelset_function, val_neg, val_pos, cell_function)[source]#

Interpolates jumping values to the mesh cells based on the levelset function.

Parameters:
  • levelset_function (fenics.Function) – The levelset function representing the domain.

  • val_neg (float) – The value inside the domain (levelset_function < 0).

  • val_pos (float) – The value outside the domain (levelset_function > 0).

  • cell_function (fenics.Function) – The piecewise constant cell function, into which the result is written.

Return type:

None

cashocs.interval_mesh(n=10, start=0.0, end=1.0, partitions=None, comm=None)[source]#

Creates an 1D interval mesh starting at x=0 to x=length.

This function creates a uniform mesh of a 1D interval, starting at the start and ending at end. The resulting mesh uses n sub-intervals to discretize the geometry. The boundary markers are as follows:

  • 1 corresponds to \(x=start\)

  • 2 corresponds to \(x=end\)

Parameters:
  • n (int) – Number of elements for discretizing the interval, default is 10

  • start (float) – The start of the interval, default is 0.0

  • end (float) – The end of the interval, default is 1.0

  • partitions (Optional[List[float]]) – Points in the interval at which a partition in subdomains should be made. The resulting volume measure is sorted ascendingly according to the sub-intervals defined in partitions (starting at 1). Defaults to None.

  • comm (Optional[MPI.Comm]) – MPI communicator that is to be used for creating the mesh.

Returns:

A tuple (mesh, subdomains, boundaries, dx, ds, dS), where mesh is the imported FEM mesh, subdomains is a mesh function for the subdomains, boundaries is a mesh function for the boundaries, dx is a volume measure, ds is a surface measure, and dS is a measure for the interior facets.

Return type:

_typing.MeshTuple

cashocs.linear_solve(linear_form, u, bcs, ksp_options=None, preconditioner_form=None, A_tensor=None, b_tensor=None, linear_solver=None)[source]#

Solves a linear problem.

Parameters:
  • linear_form (ufl.Form) – The linear variational form of the problem, i.e., linear_form == 0

  • u (fenics.Function) – The function to be solved for

  • bcs (Union[fenics.DirichletBC, List[fenics.DirichletBC]]) – The boundary conditions for the problem

  • ksp_options (Optional[_typing.KspOption]) – The options for the PETSc KSP solver, optional. Default is None, where the linear solver MUMPS is used

  • A_tensor (Optional[fenics.PETScMatrix]) – A fenics.PETScMatrix for storing the left-hand side of the linear sub-problem.

  • b_tensor (Optional[fenics.PETScVector]) – A fenics.PETScVector for storing the right-hand side of the linear sub-problem.

  • linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver used to solve the (discretized) linear problem.

  • preconditioner_form (ufl.Form)

Returns:

The computed solution, this overwrites the input function u.

Return type:

fenics.Function

cashocs.load_config(path)[source]#

Loads a config object from a config file.

Loads the config from a .ini file via the configparser package.

Parameters:

path (str) – The path to the .ini file storing the configuration.

Returns:

The output config file, which includes the path to the .ini file.

Return type:

ConfigParser

cashocs.newton_solve(nonlinear_form, u, bcs, derivative=None, shift=None, rtol=1e-10, atol=1e-10, max_iter=50, convergence_type='combined', norm_type='l2', damped=True, inexact=True, verbose=True, ksp_options=None, A_tensor=None, b_tensor=None, is_linear=False, preconditioner_form=None, linear_solver=None)[source]#

Solves a nonlinear problem with Newton's method.

Parameters:
  • nonlinear_form (ufl.Form) – The variational form of the nonlinear problem to be solved by Newton’s method.

  • u (fenics.Function) – The sought solution / initial guess. It is not assumed that the initial guess satisfies the Dirichlet boundary conditions, they are applied automatically. The method overwrites / updates this Function.

  • bcs (Union[fenics.DirichletBC, List[fenics.DirichletBC]]) – A list of DirichletBCs for the nonlinear variational problem.

  • derivative (Optional[ufl.Form]) – The Jacobian of nonlinear_form, used for the Newton method. Default is None, and in this case the Jacobian is computed automatically with AD.

  • shift (Optional[ufl.Form]) – A shift term, if the right-hand side of the nonlinear problem is not zero, but shift.

  • rtol (float) – Relative tolerance of the solver if convergence_type is either 'combined' or 'rel' (default is rtol = 1e-10).

  • atol (float) – Absolute tolerance of the solver if convergence_type is either 'combined' or 'abs' (default is atol = 1e-10).

  • max_iter (int) – Maximum number of iterations carried out by the method (default is max_iter = 50).

  • convergence_type (Literal['combined', 'rel', 'abs']) – Determines the type of stopping criterion that is used.

  • norm_type (Literal['l2', 'linf']) – Determines which norm is used in the stopping criterion.

  • damped (bool) – If True, then a damping strategy is used. If False, the classical Newton-Raphson iteration (without damping) is used (default is True).

  • inexact (bool) – If True, an inexact Newton's method is used. Default is True.

  • verbose (bool) – If True, prints status of the iteration to the console (default is True).

  • ksp_options (Optional[_typing.KspOption]) – The list of options for the linear solver.

  • A_tensor (Optional[fenics.PETScMatrix]) – A fenics.PETScMatrix for storing the left-hand side of the linear sub-problem.

  • b_tensor (Optional[fenics.PETScVector]) – A fenics.PETScVector for storing the right-hand side of the linear sub-problem.

  • is_linear (bool) – A boolean flag, which indicates whether the problem is actually linear.

  • preconditioner_form (Optional[ufl.Form]) – A UFL form which defines the preconditioner matrix.

  • linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the linear systems arising from the discretized PDE.

Returns:

The solution of the nonlinear variational problem, if converged. This overwrites the input function u.

Return type:

fenics.Function

Examples

Consider the problem

\[\begin{split}\begin{alignedat}{2} - \Delta u + u^3 &= 1 \quad &&\text{ in } \Omega=(0,1)^2 \\ u &= 0 \quad &&\text{ on } \Gamma. \end{alignedat}\end{split}\]

This is solved with the code

from fenics import *
import cashocs

mesh, _, boundaries, dx, _, _ = cashocs.regular_mesh(25)
V = FunctionSpace(mesh, 'CG', 1)

u = Function(function_space)
v = TestFunction(function_space)
F = inner(grad(u), grad(v))*dx + pow(u,3)*v*dx - Constant(1)*v*dx
bcs = cashocs.create_dirichlet_bcs(V, Constant(0.0), boundaries, [1,2,3,4])
cashocs.newton_solve(F, u, bcs)
cashocs.picard_iteration(form_list, u_list, bcs_list, max_iter=50, rtol=1e-10, atol=1e-10, verbose=True, inner_damped=True, inner_inexact=True, inner_verbose=False, inner_max_iter=25, ksp_options=None, A_tensors=None, b_tensors=None, inner_is_linear=False, preconditioner_forms=None, linear_solver=None, newton_linearizations=None)[source]#

Solves a system of coupled PDEs via a Picard iteration.

Parameters:
  • form_list (Union[List[ufl.form], ufl.Form]) – List of the coupled PDEs.

  • u_list (Union[List[fenics.Function], fenics.Function]) – List of the state variables (to be solved for).

  • bcs_list (Union[List[fenics.DirichletBC], List[List[fenics.DirichletBC]]]) – List of boundary conditions for the PDEs.

  • max_iter (int) – The maximum number of iterations for the Picard iteration.

  • rtol (float) – The relative tolerance for the Picard iteration, default is 1e-10.

  • atol (float) – The absolute tolerance for the Picard iteration, default is 1e-10.

  • verbose (bool) – Boolean flag, if True, output is written to stdout, default is True.

  • inner_damped (bool) – Boolean flag, if True, the inner problems are solved with a damped Newton method, default is True

  • inner_inexact (bool) – Boolean flag, if True, the inner problems are solved with an inexact Newton method, default is True

  • inner_verbose (bool) – Boolean flag, if True, the inner problems write the history to stdout, default is False.

  • inner_max_iter (int) – Maximum number of iterations for the inner Newton solver; default is 25.

  • ksp_options (Optional[List[_typing.KspOption]]) – List of options for the KSP objects.

  • A_tensors (Optional[List[fenics.PETScMatrix]]) – List of matrices for the right-hand sides of the inner (linearized) equations.

  • b_tensors (Optional[List[fenics.PETScVector]]) – List of vectors for the left-hand sides of the inner (linearized) equations.

  • inner_is_linear (bool) – Boolean flag, if this is True, all problems are actually linear ones, and only a linear solver is used.

  • preconditioner_forms (Optional[Union[List[ufl.Form], ufl.Form]]) – The list of forms for the preconditioner. The default is None, so that the preconditioner matrix is the same as the system matrix.

  • linear_solver (Optional[_utils.linalg.LinearSolver]) – The linear solver (KSP) which is used to solve the linear systems arising from the discretized PDE.

  • newton_linearizations (Optional[List[ufl.Form]]) – A list of UFL forms describing which (alternative) linearizations should be used for the (nonlinear) equations when solving them (with Newton’s method). The default is None, so that the Jacobian of the supplied state forms is used.

Return type:

None

cashocs.regular_box_mesh(n=10, start_x=0.0, start_y=0.0, start_z=None, end_x=1.0, end_y=1.0, end_z=None, diagonal='right', comm=None)[source]#

Creates a mesh corresponding to a rectangle or cube.

This function creates a uniform mesh of either a rectangle or a cube, with specified start (S_) and end points (E_). The resulting mesh uses n elements along the shortest direction and accordingly many along the longer ones. The resulting domain is

\[\begin{split}\begin{alignedat}{2} &[start_x, end_x] \times [start_y, end_y] \quad &&\text{ in } 2D, \\ &[start_x, end_x] \times [start_y, end_y] \times [start_z, end_z] \quad &&\text{ in } 3D. \end{alignedat}\end{split}\]

The boundary markers are ordered as follows:

  • 1 corresponds to \(x=start_x\).

  • 2 corresponds to \(x=end_x\).

  • 3 corresponds to \(y=start_y\).

  • 4 corresponds to \(y=end_y\).

  • 5 corresponds to \(z=start_z\) (only in 3D).

  • 6 corresponds to \(z=end_z\) (only in 3D).

Parameters:
  • n (int) – Number of elements in the shortest coordinate direction.

  • start_x (float) – Start of the x-interval.

  • start_y (float) – Start of the y-interval.

  • start_z (Optional[float]) – Start of the z-interval, mesh is 2D if this is None (default is None).

  • end_x (float) – End of the x-interval.

  • end_y (float) – End of the y-interval.

  • end_z (Optional[float]) – End of the z-interval, mesh is 2D if this is None (default is None).

  • diagonal (Literal['right', 'left', 'left/right', 'right/left', 'crossed']) – This defines the type of diagonal used to create the box mesh in 2D. This can be one of "right", "left", "left/right", "right/left" or "crossed".

  • comm (Optional[MPI.Comm]) – MPI communicator that is to be used for creating the mesh.

Returns:

A tuple (mesh, subdomains, boundaries, dx, ds, dS), where mesh is the imported FEM mesh, subdomains is a mesh function for the subdomains, boundaries is a mesh function for the boundaries, dx is a volume measure, ds is a surface measure, and dS is a measure for the interior facets.

Return type:

_typing.MeshTuple

cashocs.regular_mesh(n=10, length_x=1.0, length_y=1.0, length_z=None, diagonal='right', comm=None)[source]#

Creates a mesh corresponding to a rectangle or cube.

This function creates a uniform mesh of either a rectangle or a cube, starting at the origin and having length specified in length_x, length_y, and length_z. The resulting mesh uses n elements along the shortest direction and accordingly many along the longer ones. The resulting domain is

\[\begin{split}\begin{alignedat}{2} &[0, length_x] \times [0, length_y] \quad &&\text{ in } 2D, \\ &[0, length_x] \times [0, length_y] \times [0, length_z] \quad &&\text{ in } 3D. \end{alignedat}\end{split}\]

The boundary markers are ordered as follows:

  • 1 corresponds to \(x=0\).

  • 2 corresponds to \(x=length_x\).

  • 3 corresponds to \(y=0\).

  • 4 corresponds to \(y=length_y\).

  • 5 corresponds to \(z=0\) (only in 3D).

  • 6 corresponds to \(z=length_z\) (only in 3D).

Parameters:
  • n (int) – Number of elements in the shortest coordinate direction.

  • length_x (float) – Length in x-direction.

  • length_y (float) – Length in y-direction.

  • length_z (Optional[float]) – Length in z-direction, if this is None, then the geometry will be two-dimensional (default is None).

  • diagonal (Literal['left', 'right', 'left/right', 'right/left', 'crossed']) – This defines the type of diagonal used to create the box mesh in 2D. This can be one of "right", "left", "left/right", "right/left" or "crossed".

  • comm (Optional[MPI.Comm]) – MPI communicator that is to be used for creating the mesh.

Returns:

A tuple (mesh, subdomains, boundaries, dx, ds, dS), where mesh is the imported FEM mesh, subdomains is a mesh function for the subdomains, boundaries is a mesh function for the boundaries, dx is a volume measure, ds is a surface measure, and dS is a measure for the interior facets.

Return type:

_typing.MeshTuple

cashocs.set_log_level(level)[source]#

Determines the log level of cashocs.

Can be used to show, e.g., info and warning messages or to hide them. There are a total of five different levels for the logs: DEBUG, INFO, WARNING, ERROR, and CRITICAL. The usage of this method is explained in the examples section.

Parameters:

level (int) – Should be one of cashocs.LogLevel.DEBUG, cashocs.LogLevel.INFO, cashocs.LogLevel.WARNING, cashocs.LogLevel.ERROR, cashocs.LogLevel.CRITICAL

Return type:

None

Notes

The log level setting is global, so if you use this interactively, you have to restart / reload your interactive console to return to the default settings.

Examples

To set the log level of cashocs, use this method as follows:

import cashocs

cashocs.set_log_level(cashocs.LogLevel.WARNING)

which only shows messages with a level of WARNING or higher. To use a different level, replace WARNING by DEBUG, INFO, ERROR, or CRITICAL.

Modules

cashocs.geometry

Mesh generation, quality, and management tools.

cashocs.io

Inputs and outputs.

cashocs.nonlinear_solvers

Custom solvers for nonlinear equations.

cashocs.space_mapping

Space mapping algorithms.