bice.continuation package
Submodules
bice.continuation.bifurcations module
Bifurcation tracking classes for the bice package.
- class BifurcationConstraint(phi: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]], free_parameter: tuple[object, str])
Bases:
EquationImplements a constraint equation for bifurcation tracking.
This class extends the original problem with an additional constraint that forces the system to remain at a bifurcation point (e.g., by ensuring the existence of a null vector).
Initialize the BifurcationConstraint.
- Parameters:
phi – The initial guess for the null-eigenvector.
free_parameter – A tuple (object, attribute_name) specifying the parameter to vary.
- actions_after_newton_solve() None
Perform actions after the Newton solver has finished.
Updates the free parameter from the stored unknowns.
- actions_before_evaluation(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) None
Perform actions before evaluating the system.
Updates the free parameter from the unknowns vector.
- Parameters:
u – The vector of unknowns.
- free_parameter
reference to the free parameter
- jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray | spmatrix | float
Calculate the Jacobian of the extended system.
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian matrix or 0 if disabled.
- Return type:
Union[Matrix, float]
- mass_matrix() float
Return the mass matrix contribution.
- Returns:
Always 0 as it couples to no time-derivatives.
- Return type:
float
- original_jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the original (unextended) Jacobian of the problem.
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian of the original unconstrained system.
- Return type:
Array
- plot(ax: Any) None
Plot the constraint state (no-op).
- Parameters:
ax – The matplotlib axes.
- rhs(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the right-hand side of the extended system.
- Parameters:
u – The vector of unknowns.
- Returns:
The residuals of the extended system.
- Return type:
Array
- Raises:
Exception – If the dimension of the problem mismatches the null-eigenvector.
bice.continuation.constraints module
Predefined constraint equations for continuation.
- class ConstraintEquation(shape: int | tuple[int, ...] = (1,))
Bases:
EquationAbstract base class for constraint type equations.
For simple implementation of PDE constraints and less redundant code.
Initialize the ConstraintEquation.
- Parameters:
shape – The shape of the constraint unknowns (Lagrange multipliers).
- mass_matrix() float
Return the mass matrix contribution.
- Returns:
Always 0 as constraints usually couple to no time-derivatives.
- Return type:
float
- plot(ax: Any) None
Plot the constraint state (no-op).
- Parameters:
ax – The matplotlib axes.
- class TranslationConstraint(reference_equation: PartialDifferentialEquation, variable: int | None = None, direction: int = 0)
Bases:
ConstraintEquationAssures that the center of mass does not move.
A translation constraint assures that the center of mass of some reference equation’s unknowns does not move when solving the system. The additional constraint equations (one per spatial dimension) come with Lagrange multipliers, that correspond to the velocities of a comoving frame (advection term).
Initialize the TranslationConstraint.
- Parameters:
reference_equation – The equation to constrain.
variable – The variable index to constrain.
direction – The spatial direction (index) to apply the constraint to.
- direction
which spatial direction (index of [x, y, …]) should the constraint apply to
- jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) csr_matrix
Calculate the analytical Jacobian of the translation constraint.
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian matrix.
- Return type:
sp.csr_matrix
- ref_eq
on which equation/unknowns should the constraint be imposed?
- rhs(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the residuals of the translation constraint.
- Parameters:
u – The vector of unknowns.
- Returns:
The residuals vector.
- Return type:
Array
- u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
the unknowns (velocity vector)
- variable
on which variable (index) of the equation should the constraint be imposed?
- class VolumeConstraint(reference_equation: Equation, variable: int | None = None)
Bases:
ConstraintEquationAssures the conservation of the integral of the unknowns.
A volume constraint (or mass constraint) assures the conservation of the integral of the unknowns of some given equation when solving the system. We may even prescribe the target volume (or mass) with a parameter, but we don’t have to. The constraint equation comes with an additional (unknown) Lagrange multiplier that can be interpreted as an influx into the system.
Initialize the VolumeConstraint.
- Parameters:
reference_equation – The equation to constrain.
variable – The index of the variable to constrain (if the equation has multiple).
- fixed_volume: float | None
This parameter allows for prescribing a fixed volume (unless it is None)
- jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) csr_matrix
Calculate the analytical Jacobian of the volume constraint.
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian matrix.
- Return type:
sp.csr_matrix
- ref_eq
on which equation/unknowns should the constraint be imposed?
- rhs(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the residuals of the volume constraint.
- Parameters:
u – The vector of unknowns.
- Returns:
The residuals vector.
- Return type:
Array
- u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
this equation brings a single extra degree of freedom (influx Lagrange multiplier)
- variable
on which variable (index) of the equation should the constraint be imposed?
bice.continuation.continuation_steppers module
Parameter continuation stepping strategies.
- class ContinuationStepper(ds: float = 0.001)
Bases:
objectAbstract base class for all parameter continuation-steppers.
Specifies attributes and methods that all continuation-steppers should have.
Initialize the ContinuationStepper.
- Parameters:
ds – The initial continuation step size.
- ds
continuation step size
- factory_reset() None
Reset the continuation-stepper parameters & storage to default.
Useful when starting off a new solution point, switching branches or switching the principal continuation parameter.
- class NaturalContinuation(ds: float = 0.001)
Bases:
ContinuationStepperNatural parameter continuation stepper.
The simplest continuation method, where the parameter is incremented by a fixed step size, and the new solution is found using Newton’s method.
Initialize the ContinuationStepper.
- Parameters:
ds – The initial continuation step size.
- class PseudoArclengthContinuation(ds: float = 0.001)
Bases:
ContinuationStepperPseudo-arclength parameter continuation stepper.
A robust continuation method that can follow solution branches through limit points (folds) by adding an arclength constraint to the system.
Initialize the PseudoArclengthContinuation.
- Parameters:
ds – The initial step size.
- adapt_stepsize
should the step size be adapted while stepping?
- convergence_tolerance
convergence tolerance for the newton solver in the continuation step
- ds_decrease_factor
ds decreases by this factor when less than desired_newton_steps are performed
- ds_increase_factor
ds increases by this factor when more than desired_newton_steps are performed
- ds_max
maximum step size
- ds_min
minimum step size
- fd_epsilon
finite-difference for calculating parameter derivatives
- max_newton_iterations
maximum number of newton iterations for solving
- ndesired_newton_steps
the desired number of newton iterations for solving, step size is adapted if we over/undershoot this number
- nnewton_iter_taken: int | None
the actual number of newton iterations taken in the last continuation step
- parameter_arc_length_proportion
Rescale the parameter constraint, for numerical stability. May be decreased, e.g. for very sharp folds.
- step(problem: Problem) None
Perform a pseudo-arclength continuation step on a problem.
Calculates the tangent predictor and corrects with a Newton-Raphson iteration on the extended system.
- Parameters:
problem – The problem to step.
- Raises:
np.linalg.LinAlgError – If the Newton solver fails to converge even after step size reduction.
bice.continuation.deflation module
Deflation operator for detecting disconnected branches.
- class DeflationOperator
Bases:
objectA deflation operator M for deflated continuation.
Adds singularities to the equation at given solutions u_i: 0 = F(u) –> 0 = M(u) * F(u) with M(u) = product_i <u_i - u, u_i - u>^-p + shift
- The parameters are:
p: some exponent to the norm <u, v> shift: some constant added shift parameter for numerical stability
Initialize the DeflationOperator.
- D_operator(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the Jacobian of deflation operator for given u.
- Parameters:
u – The vector of unknowns.
- Returns:
The gradient of the operator.
- Return type:
Array
- add_solution(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) None
Add a solution to the list of solutions used for deflation.
- Parameters:
u – The solution to deflate.
- clear_solutions() None
Clear the list of solutions used for deflation.
- deflated_jacobian(rhs: Callable[[ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]], ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]], jacobian: Callable[[ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]], ndarray | spmatrix]) Callable[[ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]], ndarray | spmatrix]
Generate Jacobian of deflated rhs of some equation or problem.
- Parameters:
rhs – The original right-hand side function.
jacobian – The original Jacobian function.
- Returns:
The deflated Jacobian function.
- Return type:
callable
- deflated_rhs(rhs: Callable[[ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]], ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]]) Callable[[ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]], ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]]
Deflate the rhs of some equation.
Returns a new function that represents M(u) * rhs(u).
- Parameters:
rhs – The original right-hand side function.
- Returns:
The deflated rhs function.
- Return type:
callable
- operator(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) float
Obtain the value of the deflation operator for given u.
- Parameters:
u – The vector of unknowns.
- Returns:
The value of the operator.
- Return type:
float
- p
the order of the norm that will be used for the deflation operator
- remove_solution(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) None
Remove a solution from the list of solutions used for deflation.
- Parameters:
u – The solution to remove.
- shift
small constant in the deflation operator, for numerical stability
- solutions: list[ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]]
list of solutions, that will be suppressed by the deflation operator
bice.continuation.timeperiodic module
Time-periodic orbit handling.
- class TimePeriodicOrbitHandler(reference_equation: Equation, T: float, Nt: int)
Bases:
EquationHandles equations with implicit time-stepping on a periodic time-mesh.
The TimePeriodicOrbitHandler is an Equation that has an implicit time-stepping on a periodic time-mesh of (unknown) period length T. In a Problem, it can be used for solving periodic orbits, e.g. in a Hopf continuation. The referenced equation should not simultaneously be part of the problem.
Initialize the TimePeriodicOrbitHandler.
- Parameters:
reference_equation – The equation to solve the orbit for.
T – Initial guess for the period.
Nt – Initial number of time steps.
- property Nt: int
Return the number of points in time.
- Returns:
Number of time steps.
- Return type:
int
- property T: float
Access the period length.
- Returns:
The period length.
- Return type:
float
- adapt() tuple[float, float]
Adapt the time mesh to the solution.
- Returns:
(min_error, max_error) estimates.
- Return type:
tuple
- build_ddt_matrix() csr_matrix
Build the time-derivative operator ddt, using periodic finite differences.
- Returns:
The finite difference matrix.
- Return type:
sp.csr_matrix
- floquet_multipliers(k: int = 20, use_cache: bool = True) ndarray[tuple[Any, ...], dtype[float64]]
Calculate the Floquet multipliers to obtain the stability of the orbit.
The Floquet multipliers are the eigenvalues of the monodromy matrix.
- Parameters:
k – Number of desired Floquet multipliers to be calculated by the iterative eigensolver.
use_cache – Whether to use cached Jacobians.
- Returns:
The Floquet multipliers (eigenvalues).
- Return type:
RealArray
- jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) csr_matrix
Calculate the Jacobian of rhs(u).
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian matrix.
- Return type:
sp.csr_matrix
- load(data: dict[str, Any]) None
Load the state of the equation, including the dt-values.
- Parameters:
data – The state dictionary.
- monodromy_matrix(use_cache: bool = True) csr_matrix
Calculate the monodromy matrix A.
Used to calculate the stability of the orbit using the Floquet multipliers (eigenvalues of A).
The matrix is computed as the product of the reference equations’s Jacobians for each point in time, i.e.: A = J[u(t_N), t_N] * J[u(t_{N-1}), t_{N-1}] * … * J[u(t_0), t_0]
- Parameters:
use_cache – Whether to use cached Jacobians from the last solving step.
- Returns:
The monodromy matrix.
- Return type:
sp.csr_matrix
- plot(ax: Any) None
Plot the solutions for different timesteps.
- Parameters:
ax – The matplotlib axes.
- rhs(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the rhs of the full system of equations.
- Parameters:
u – The vector of unknowns (including period T).
- Returns:
The residuals.
- Return type:
Array
- save() dict[str, Any]
Save the state of the equation, including the dt-values.
- Returns:
The state dictionary.
- Return type:
dict
- property t: ndarray[tuple[Any, ...], dtype[float64]]
Return the temporal domain vector.
- Returns:
The accumulated time steps.
- Return type:
RealArray
- u_orbit() ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Return the unknowns in separate arrays for each point in time.
- Returns:
The reshaped unknowns array of shape (Nt, …).
- Return type:
Array
Module contents
Continuation and bifurcation tracking functionality.
This package provides classes for path continuation (natural, pseudo-arclength) and bifurcation analysis (branch switching, tracking).
- class BifurcationConstraint(phi: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]], free_parameter: tuple[object, str])
Bases:
EquationImplements a constraint equation for bifurcation tracking.
This class extends the original problem with an additional constraint that forces the system to remain at a bifurcation point (e.g., by ensuring the existence of a null vector).
Initialize the BifurcationConstraint.
- Parameters:
phi – The initial guess for the null-eigenvector.
free_parameter – A tuple (object, attribute_name) specifying the parameter to vary.
- actions_after_newton_solve() None
Perform actions after the Newton solver has finished.
Updates the free parameter from the stored unknowns.
- actions_before_evaluation(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) None
Perform actions before evaluating the system.
Updates the free parameter from the unknowns vector.
- Parameters:
u – The vector of unknowns.
- free_parameter
reference to the free parameter
- jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray | spmatrix | float
Calculate the Jacobian of the extended system.
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian matrix or 0 if disabled.
- Return type:
Union[Matrix, float]
- mass_matrix() float
Return the mass matrix contribution.
- Returns:
Always 0 as it couples to no time-derivatives.
- Return type:
float
- original_jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the original (unextended) Jacobian of the problem.
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian of the original unconstrained system.
- Return type:
Array
- plot(ax: Any) None
Plot the constraint state (no-op).
- Parameters:
ax – The matplotlib axes.
- rhs(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the right-hand side of the extended system.
- Parameters:
u – The vector of unknowns.
- Returns:
The residuals of the extended system.
- Return type:
Array
- Raises:
Exception – If the dimension of the problem mismatches the null-eigenvector.
- class ConstraintEquation(shape: int | tuple[int, ...] = (1,))
Bases:
EquationAbstract base class for constraint type equations.
For simple implementation of PDE constraints and less redundant code.
Initialize the ConstraintEquation.
- Parameters:
shape – The shape of the constraint unknowns (Lagrange multipliers).
- mass_matrix() float
Return the mass matrix contribution.
- Returns:
Always 0 as constraints usually couple to no time-derivatives.
- Return type:
float
- plot(ax: Any) None
Plot the constraint state (no-op).
- Parameters:
ax – The matplotlib axes.
- class NaturalContinuation(ds: float = 0.001)
Bases:
ContinuationStepperNatural parameter continuation stepper.
The simplest continuation method, where the parameter is incremented by a fixed step size, and the new solution is found using Newton’s method.
Initialize the ContinuationStepper.
- Parameters:
ds – The initial continuation step size.
- class PseudoArclengthContinuation(ds: float = 0.001)
Bases:
ContinuationStepperPseudo-arclength parameter continuation stepper.
A robust continuation method that can follow solution branches through limit points (folds) by adding an arclength constraint to the system.
Initialize the PseudoArclengthContinuation.
- Parameters:
ds – The initial step size.
- adapt_stepsize
should the step size be adapted while stepping?
- convergence_tolerance
convergence tolerance for the newton solver in the continuation step
- ds_decrease_factor
ds decreases by this factor when less than desired_newton_steps are performed
- ds_increase_factor
ds increases by this factor when more than desired_newton_steps are performed
- ds_max
maximum step size
- ds_min
minimum step size
- fd_epsilon
finite-difference for calculating parameter derivatives
- max_newton_iterations
maximum number of newton iterations for solving
- ndesired_newton_steps
the desired number of newton iterations for solving, step size is adapted if we over/undershoot this number
- nnewton_iter_taken: int | None
the actual number of newton iterations taken in the last continuation step
- parameter_arc_length_proportion
Rescale the parameter constraint, for numerical stability. May be decreased, e.g. for very sharp folds.
- step(problem: Problem) None
Perform a pseudo-arclength continuation step on a problem.
Calculates the tangent predictor and corrects with a Newton-Raphson iteration on the extended system.
- Parameters:
problem – The problem to step.
- Raises:
np.linalg.LinAlgError – If the Newton solver fails to converge even after step size reduction.
- class TimePeriodicOrbitHandler(reference_equation: Equation, T: float, Nt: int)
Bases:
EquationHandles equations with implicit time-stepping on a periodic time-mesh.
The TimePeriodicOrbitHandler is an Equation that has an implicit time-stepping on a periodic time-mesh of (unknown) period length T. In a Problem, it can be used for solving periodic orbits, e.g. in a Hopf continuation. The referenced equation should not simultaneously be part of the problem.
Initialize the TimePeriodicOrbitHandler.
- Parameters:
reference_equation – The equation to solve the orbit for.
T – Initial guess for the period.
Nt – Initial number of time steps.
- property Nt: int
Return the number of points in time.
- Returns:
Number of time steps.
- Return type:
int
- property T: float
Access the period length.
- Returns:
The period length.
- Return type:
float
- adapt() tuple[float, float]
Adapt the time mesh to the solution.
- Returns:
(min_error, max_error) estimates.
- Return type:
tuple
- build_ddt_matrix() csr_matrix
Build the time-derivative operator ddt, using periodic finite differences.
- Returns:
The finite difference matrix.
- Return type:
sp.csr_matrix
- floquet_multipliers(k: int = 20, use_cache: bool = True) ndarray[tuple[Any, ...], dtype[float64]]
Calculate the Floquet multipliers to obtain the stability of the orbit.
The Floquet multipliers are the eigenvalues of the monodromy matrix.
- Parameters:
k – Number of desired Floquet multipliers to be calculated by the iterative eigensolver.
use_cache – Whether to use cached Jacobians.
- Returns:
The Floquet multipliers (eigenvalues).
- Return type:
RealArray
- jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) csr_matrix
Calculate the Jacobian of rhs(u).
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian matrix.
- Return type:
sp.csr_matrix
- load(data: dict[str, Any]) None
Load the state of the equation, including the dt-values.
- Parameters:
data – The state dictionary.
- monodromy_matrix(use_cache: bool = True) csr_matrix
Calculate the monodromy matrix A.
Used to calculate the stability of the orbit using the Floquet multipliers (eigenvalues of A).
The matrix is computed as the product of the reference equations’s Jacobians for each point in time, i.e.: A = J[u(t_N), t_N] * J[u(t_{N-1}), t_{N-1}] * … * J[u(t_0), t_0]
- Parameters:
use_cache – Whether to use cached Jacobians from the last solving step.
- Returns:
The monodromy matrix.
- Return type:
sp.csr_matrix
- plot(ax: Any) None
Plot the solutions for different timesteps.
- Parameters:
ax – The matplotlib axes.
- rhs(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the rhs of the full system of equations.
- Parameters:
u – The vector of unknowns (including period T).
- Returns:
The residuals.
- Return type:
Array
- save() dict[str, Any]
Save the state of the equation, including the dt-values.
- Returns:
The state dictionary.
- Return type:
dict
- property t: ndarray[tuple[Any, ...], dtype[float64]]
Return the temporal domain vector.
- Returns:
The accumulated time steps.
- Return type:
RealArray
- u_orbit() ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Return the unknowns in separate arrays for each point in time.
- Returns:
The reshaped unknowns array of shape (Nt, …).
- Return type:
Array
- class TranslationConstraint(reference_equation: PartialDifferentialEquation, variable: int | None = None, direction: int = 0)
Bases:
ConstraintEquationAssures that the center of mass does not move.
A translation constraint assures that the center of mass of some reference equation’s unknowns does not move when solving the system. The additional constraint equations (one per spatial dimension) come with Lagrange multipliers, that correspond to the velocities of a comoving frame (advection term).
Initialize the TranslationConstraint.
- Parameters:
reference_equation – The equation to constrain.
variable – The variable index to constrain.
direction – The spatial direction (index) to apply the constraint to.
- direction
which spatial direction (index of [x, y, …]) should the constraint apply to
- jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) csr_matrix
Calculate the analytical Jacobian of the translation constraint.
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian matrix.
- Return type:
sp.csr_matrix
- ref_eq
on which equation/unknowns should the constraint be imposed?
- rhs(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the residuals of the translation constraint.
- Parameters:
u – The vector of unknowns.
- Returns:
The residuals vector.
- Return type:
Array
- u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
the unknowns (velocity vector)
- variable
on which variable (index) of the equation should the constraint be imposed?
- class VolumeConstraint(reference_equation: Equation, variable: int | None = None)
Bases:
ConstraintEquationAssures the conservation of the integral of the unknowns.
A volume constraint (or mass constraint) assures the conservation of the integral of the unknowns of some given equation when solving the system. We may even prescribe the target volume (or mass) with a parameter, but we don’t have to. The constraint equation comes with an additional (unknown) Lagrange multiplier that can be interpreted as an influx into the system.
Initialize the VolumeConstraint.
- Parameters:
reference_equation – The equation to constrain.
variable – The index of the variable to constrain (if the equation has multiple).
- fixed_volume: float | None
This parameter allows for prescribing a fixed volume (unless it is None)
- jacobian(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) csr_matrix
Calculate the analytical Jacobian of the volume constraint.
- Parameters:
u – The vector of unknowns.
- Returns:
The Jacobian matrix.
- Return type:
sp.csr_matrix
- ref_eq
on which equation/unknowns should the constraint be imposed?
- rhs(u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]) ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
Calculate the residuals of the volume constraint.
- Parameters:
u – The vector of unknowns.
- Returns:
The residuals vector.
- Return type:
Array
- u: ndarray[tuple[Any, ...], dtype[float64 | complexfloating]]
this equation brings a single extra degree of freedom (influx Lagrange multiplier)
- variable
on which variable (index) of the equation should the constraint be imposed?