pymor.models package¶
Submodules¶
basic module¶
-
class
pymor.models.basic.InstationaryModel(*args, **kwargs)[source]¶ Bases:
pymor.models.interface.ModelGeneric class for models of instationary problems.
This class describes instationary problems given by the equations:
M * ∂_t u(t, μ) + L(u(μ), t, μ) = F(t, μ) u(0, μ) = u_0(μ)for t in [0,T], where L is a (possibly non-linear) time-dependent
Operator, F is a time-dependent vector-likeOperator, and u_0 the initial data. The massOperatorM is assumed to be linear.Parameters
- T
The final time T.
- initial_data
The initial data
u_0. Either aVectorArrayof length 1 or (for theParameter-dependent case) a vector-likeOperator(i.e. a linearOperatorwithsource.dim == 1) which applied toNumpyVectorArray(np.array([1]))will yield the initial data for givenparameter values.- operator
The
OperatorL.- rhs
The right-hand side F.
- mass
The mass
OperatorM. IfNone, the identity is assumed.- time_stepper
The
time-stepperto be used bysolve.- num_values
The number of returned vectors of the solution trajectory. If
None, each intermediate vector that is calculated is returned.- output_functional
Operatormapping a given solution to the model output. In many applications, this will be aFunctional, i.e. anOperatormapping to scalars. This is not required, however.- products
A dict of product
Operatorsdefined on the discrete space the problem is posed on. For each product with key'x'a corresponding attributex_product, as well as a norm methodx_normis added to the model.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, m)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, m, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the model.
Methods
to_lti,with_time_steppercompute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
_compute_solution(mu=None, **kwargs)[source]¶ Compute the model’s solution for
parameter valuesmu.This method is called by the default implementation of
computeinpymor.models.interface.Model.Parameters
- mu
Parameter valuesfor which to compute the solution.- kwargs
Additional keyword arguments to customize how the solution is computed or to select additional data to be returned.
Returns
VectorArraywith the computed solution or a dict which at least must contain the key'solution'.
-
class
pymor.models.basic.StationaryModel(*args, **kwargs)[source]¶ Bases:
pymor.models.interface.ModelGeneric class for models of stationary problems.
This class describes discrete problems given by the equation:
L(u(μ), μ) = F(μ)
with a vector-like right-hand side F and a (possibly non-linear) operator L.
Note that even when solving a variational formulation where F is a functional and not a vector, F has to be specified as a vector-like
Operator(mapping scalars to vectors). This ensures that in the complex case both L and F are anti-linear in the test variable.Parameters
- operator
The
OperatorL.- rhs
The vector F. Either a
VectorArrayof length 1 or a vector-likeOperator.- output_functional
Operatormapping a given solution to the model output. In many applications, this will be aFunctional, i.e. anOperatormapping to scalars. This is not required, however.- products
A dict of inner product
Operatorsdefined on the discrete space the problem is posed on. For each product with key'x'a corresponding attributex_product, as well as a norm methodx_normis added to the model.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, m)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, m, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the model.
Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
_compute_output_d_mu(solution, mu, return_array=False, use_adjoint=None)[source]¶ compute the gradient of the output functional w.r.t. the parameters
Parameters
- solution
Internal model state for the given
Parameter value- mu
Parameter valuefor which to compute the gradient- return_array
if
True, return the output gradient as aNumPy array. Otherwise, return a dict of gradients for eachParameter.- use_adjoint
if
Noneuse standard approach, ifTrue, use the adjoint solution for a more efficient way of computing the gradient. See Section 1.6.2 in [HPUU09] for more details. So far, the adjoint approach is only valid for linear models.
Returns
The gradient as a
NumPy arrayor a dict ofNumPy arrays.
-
_compute_solution(mu=None, **kwargs)[source]¶ Compute the model’s solution for
parameter valuesmu.This method is called by the default implementation of
computeinpymor.models.interface.Model.Parameters
- mu
Parameter valuesfor which to compute the solution.- kwargs
Additional keyword arguments to customize how the solution is computed or to select additional data to be returned.
Returns
VectorArraywith the computed solution or a dict which at least must contain the key'solution'.
-
_compute_solution_d_mu_single_direction(parameter, index, solution, mu)[source]¶ Compute the partial derivative of the solution w.r.t. a parameter index
Parameters
- parameter
parameter for which to compute the sensitivity
- index
parameter index for which to compute the sensitivity
- solution
Internal model state for the given
Parameter value.- mu
Parameter valuefor which to solve
Returns
The sensitivity of the solution as a
VectorArray.
interface module¶
-
class
pymor.models.interface.Model(*args, **kwargs)[source]¶ Bases:
pymor.core.cache.CacheableObject,pymor.parameters.base.ParametricObjectInterface for model objects.
A model object defines a discrete problem via its
classand theOperatorsit contains. Furthermore, models can besolvedfor givenparameter valuesresulting in a solutionVectorArray.Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
solution_space¶ VectorSpaceof the solutionVectorArraysreturned bysolve.
-
linear¶ Trueif the model describes a linear problem.
-
products¶ Dict of inner product operators associated with the model.
-
_compute_output(solution, mu=None, **kwargs)[source]¶ Compute the model’s output for
parameter valuesmu.This method is called by the default implementation of
computeinpymor.models.interface.Model. The assumption is made that the output is a derived quantity from the model’s internal state as returned my_compute_solution. When this is not the case, the computation of the output should be implemented in_compute.Note
The default implementation applies the
Operatorgiven by theoutput_functionalattribute to the givensolutionVectorArray.Parameters
- solution
Internal model state for the given
parameter values.- mu
Parameter valuesfor which to compute the output.- kwargs
Additional keyword arguments to customize how the output is computed or to select additional data to be returned.
Returns
NumPy arraywith the computed output or a dict which at least must contain the key'output'.
-
_compute_output_d_mu(solution, mu=None, return_array=False, **kwargs)[source]¶ compute the gradient w.r.t. the parameter of the output functional
Parameters
- solution
Internal model state for the given
Parameter value.- mu
Parameter valuefor which to compute the gradient- return_array
if
True, return the output gradient as aNumPy array. Otherwise, return a dict of gradients for eachParameter.
Returns
The gradient as a
NumPy arrayor a dict ofNumPy arrays.
-
_compute_output_error_estimate(solution, mu=None, **kwargs)[source]¶ Compute an error estimate for the computed model output.
This method is called by the default implementation of
computeinpymor.models.interface.Model. The assumption is made that the error estimate is a derived quantity from the model’s internal state as returned my_compute_solution. When this is not the case, the computation of the error estimate should be implemented in_compute.Note
The default implementation calls the
estimate_output_errormethod of the object given by theerror_estimatorattribute, passingsolution,mu,selfand**kwargs.Parameters
- solution
Internal model state for the given
parameter values.- mu
Parameter valuesfor which to compute the error estimate.- kwargs
Additional keyword arguments to customize how the error estimate is computed or to select additional data to be returned.
Returns
The computed error estimate or a dict which at least must contain the key
'solution_error_estimate'.
-
_compute_solution(mu=None, **kwargs)[source]¶ Compute the model’s solution for
parameter valuesmu.This method is called by the default implementation of
computeinpymor.models.interface.Model.Parameters
- mu
Parameter valuesfor which to compute the solution.- kwargs
Additional keyword arguments to customize how the solution is computed or to select additional data to be returned.
Returns
VectorArraywith the computed solution or a dict which at least must contain the key'solution'.
-
_compute_solution_d_mu(solution, mu=None, **kwargs)[source]¶ Compute all partial derivative of the solution w.r.t. a parameter index
Parameters
- solution
Internal model state for the given
Parameter value.- mu
Parameter valuefor which to solve
Returns
A dict of all partial sensitivities of the solution.
-
_compute_solution_d_mu_single_direction(parameter, index, solution, mu=None, **kwargs)[source]¶ Compute the partial derivative of the solution w.r.t. a parameter index
Parameters
- parameter
parameter for which to compute the sensitivity
- index
parameter index for which to compute the sensitivity
- solution
Internal model state for the given
Parameter value.- mu
Parameter valuefor which to solve
Returns
The sensitivity of the solution as a
VectorArray.
-
_compute_solution_error_estimate(solution, mu=None, **kwargs)[source]¶ Compute an error estimate for the computed internal state.
This method is called by the default implementation of
computeinpymor.models.interface.Model. The assumption is made that the error estimate is a derived quantity from the model’s internal state as returned my_compute_solution. When this is not the case, the computation of the error estimate should be implemented in_compute.Note
The default implementation calls the
estimate_errormethod of the object given by theerror_estimatorattribute, passingsolution,mu,selfand**kwargs.Parameters
- solution
Internal model state for the given
parameter values.- mu
Parameter valuesfor which to compute the error estimate.- kwargs
Additional keyword arguments to customize how the error estimate is computed or to select additional data to be returned.
Returns
The computed error estimate or a dict which at least must contain the key
'solution_error_estimate'.
-
compute(solution=False, output=False, solution_d_mu=False, output_d_mu=False, solution_error_estimate=False, output_error_estimate=False, output_d_mu_return_array=False, *, mu=None, **kwargs)[source]¶ Compute the solution of the model and associated quantities.
This methods computes the output of the model it’s internal state and various associated quantities for given
parameter valuesmu.Note
The default implementation defers the actual computations to the methods
_compute_solution,_compute_output,_compute_solution_error_estimateand_compute_output_error_estimate. The call to_compute_solutioniscached. In addition,Modelimplementors may implement_computeto simultaneously compute multiple values in an optimized way. The corresponding_compute_XXXmethods will not be called for values already returned by_compute.Parameters
- solution
If
True, return the model’s internal state.- output
If
True, return the model output.- solution_d_mu
If not
False, eitherTrueto return the derivative of the model’s internal state w.r.t. all parameter components or a tuple(parameter, index)to return the derivative of a single parameter component.- output_d_mu
If
True, return the gradient of the model output w.r.t. theParameter.- solution_error_estimate
If
True, return an error estimate for the computed internal state.- output_error_estimate
If
True, return an error estimate for the computed output.- output_d_mu_return_array
if
True, return the output gradient as aNumPy array. Otherwise, return a dict of gradients for eachParameter.- mu
Parameter valuesfor which to compute the values.- kwargs
Further keyword arguments to select further quantities that sould be returned or to customize how the values are computed.
Returns
A dict with the computed values.
-
estimate_error(mu=None, **kwargs)[source]¶ Estimate the error for the computed internal state.
For given
parameter valuesmuthis method returns an error estimate for the computed internal model state as returned bysolve. It is a convenience wrapper aroundcompute.The model error could be the error w.r.t. the analytical solution of the given problem or the model reduction error w.r.t. a corresponding high-dimensional
Model.Parameters
- mu
Parameter valuesfor which to estimate the error.- kwargs
Additional keyword arguments passed to
computethat might affect how the error estimate (or the solution) is computed.
Returns
The estimated error.
-
estimate_output_error(mu=None, **kwargs)[source]¶ Estimate the error for the computed output.
For given
parameter valuesmuthis method returns an error estimate for the computed model output as returned byoutput. It is a convenience wrapper aroundcompute.The output error could be the error w.r.t. the analytical solution of the given problem or the model reduction error w.r.t. a corresponding high-dimensional
Model.Parameters
- mu
Parameter valuesfor which to estimate the error.- kwargs
Additional keyword arguments passed to
computethat might affect how the error estimate (or the output) is computed.
Returns
The estimated error.
-
output(mu=None, return_error_estimate=False, **kwargs)[source]¶ Return the model output for given
parameter valuesmu.This method is a convenience wrapper around
compute.Parameters
- mu
Parameter valuesfor which to compute the output.- return_error_estimate
If
True, also return an error estimate for the computed output.- kwargs
Additional keyword arguments passed to
computethat might affect how the solution is computed.
Returns
The computed model output as a 2D
NumPy array. The dimension of axis 1 is : attr:dim_output. (For stationary problems, axis 0 has dimension 1. For time-dependent problems, the dimension of axis 0 depends on the number of time steps.) Whenreturn_error_estimateisTrue, the estimate is returned as second value.
-
output_d_mu(mu=None, return_array=False, **kwargs)[source]¶ compute the gradient w.r.t. the parameter of the output functional
Parameters
- mu
Parameter valuefor which to compute the gradient- return_array
if
True, return the output gradient as aNumPy array. Otherwise, return a dict of gradients for eachParameter.
Returns
The gradient as a
NumPy arrayor a dict ofNumPy arrays.
-
solve(mu=None, return_error_estimate=False, **kwargs)[source]¶ Solve the discrete problem for the
parameter valuesmu.This method returns a
VectorArraywith a internal state representation of the model’s solution for givenparameter values. It is a convenience wrapper aroundcompute.The result may be
cachedin case caching has been activated for the given model.Parameters
- mu
Parameter valuesfor which to solve.- return_error_estimate
If
True, also return an error estimate for the computed solution.- kwargs
Additional keyword arguments passed to
computethat might affect how the solution is computed.
Returns
The solution
VectorArray. Whenreturn_error_estimateisTrue, the estimate is returned as second value.
-
solve_d_mu(parameter, index, mu=None, **kwargs)[source]¶ Solve for the partial derivative of the solution w.r.t. a parameter index
Parameters
- parameter
parameter for which to compute the sensitivity
- index
parameter index for which to compute the sensitivity
- mu
Parameter valuefor which to solve
Returns
The sensitivity of the solution as a
VectorArray.
-
visualize(U, **kwargs)[source]¶ Visualize a
VectorArrayU of the model’ssolution_space.Parameters
- U
The
VectorArrayfromsolution_spacethat shall be visualized.- kwargs
Additional keyword arguments to customize the visualization. See the docstring of
self.visualizer.visualize.
-
iosys module¶
-
class
pymor.models.iosys.BilinearModel(*args, **kwargs)[source]¶ Bases:
pymor.models.iosys.InputStateOutputModelClass for bilinear systems.
This class describes input-output systems given by
\[\begin{split}E x'(t) & = A x(t) + \sum_{i = 1}^m{N_i x(t) u_i(t)} + B u(t), \\ y(t) & = C x(t) + D u(t),\end{split}\]if continuous-time, or
\[\begin{split}E x(k + 1) & = A x(k) + \sum_{i = 1}^m{N_i x(k) u_i(k)} + B u(k), \\ y(k) & = C x(k) + D u(t),\end{split}\]if discrete-time, where \(E\), \(A\), \(N_i\), \(B\), \(C\), and \(D\) are linear operators and \(m\) is the number of inputs.
Parameters
- A
The
OperatorA.- N
The tuple of
OperatorsN_i.- B
The
OperatorB.- C
The
OperatorC.- D
The
OperatorD orNone(then D is assumed to be zero).- E
The
OperatorE orNone(then E is assumed to be identity).- cont_time
Trueif the system is continuous-time, otherwiseFalse.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
ordercache_region,input_dim,output_dimparameters,parameters_inherited,parameters_internal,parameters_own,parametric-
order¶ The order of the system (equal to A.source.dim).
-
dim_input¶ The number of inputs.
-
dim_output¶ The number of outputs.
-
class
pymor.models.iosys.InputOutputModel(*args, **kwargs)[source]¶ Bases:
pymor.models.interface.ModelBase class for input-output systems.
Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
cache_region,input_dim,output_dimparameters,parameters_inherited,parameters_internal,parameters_own,parametric-
bode(w, mu=None)[source]¶ Compute magnitudes and phases.
Parameters
- w
A sequence of angular frequencies at which to compute the transfer function.
- mu
Parameter valuesfor which to evaluate the transfer function.
Returns
- mag
Transfer function magnitudes at frequencies in
w,NumPy arrayof shape(len(w), self.dim_output, self.dim_input).- phase
Transfer function phases (in radians) at frequencies in
w,NumPy arrayof shape(len(w), self.dim_output, self.dim_input).
-
bode_plot(w, mu=None, ax=None, Hz=False, dB=False, deg=True, **mpl_kwargs)[source]¶ Draw the Bode plot for all input-output pairs.
Parameters
- w
A sequence of angular frequencies at which to compute the transfer function.
- mu
Parameterfor which to evaluate the transfer function.- ax
Axis of shape (2 *
self.dim_output,self.dim_input) to which to plot. If not given,matplotlib.pyplot.gcfis used to get the figure and create axis.- Hz
Should the frequency be in Hz on the plot.
- dB
Should the magnitude be in dB on the plot.
- deg
Should the phase be in degrees (otherwise in radians).
- mpl_kwargs
Keyword arguments used in the matplotlib plot function.
Returns
- artists
List of matplotlib artists added.
-
freq_resp(w, mu=None)[source]¶ Evaluate the transfer function on the imaginary axis.
Parameters
- w
A sequence of angular frequencies at which to compute the transfer function.
- mu
Parameter valuesfor which to evaluate the transfer function.
Returns
- tfw
Transfer function values at frequencies in
w,NumPy arrayof shape(len(w), self.dim_output, self.dim_input).
-
mag_plot(w, mu=None, ax=None, ord=None, Hz=False, dB=False, **mpl_kwargs)[source]¶ Draw the magnitude plot.
Parameters
- w
A sequence of angular frequencies at which to compute the transfer function.
- mu
Parameter valuesfor which to evaluate the transfer function.- ax
Axis to which to plot. If not given,
matplotlib.pyplot.gcais used.- ord
The order of the norm used to compute the magnitude (the default is the Frobenius norm).
- Hz
Should the frequency be in Hz on the plot.
- dB
Should the magnitude be in dB on the plot.
- mpl_kwargs
Keyword arguments used in the matplotlib plot function.
Returns
- out
List of matplotlib artists added.
-
-
class
pymor.models.iosys.InputStateOutputModel(*args, **kwargs)[source]¶ Bases:
pymor.models.iosys.InputOutputModelBase class for input-output systems with state space.
Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
ordercache_region,input_dim,output_dimparameters,parameters_inherited,parameters_internal,parameters_own,parametric
-
class
pymor.models.iosys.LTIModel(*args, **kwargs)[source]¶ Bases:
pymor.models.iosys.InputStateOutputModelClass for linear time-invariant systems.
This class describes input-state-output systems given by
\[\begin{split}E(\mu) \dot{x}(t, \mu) & = A(\mu) x(t, \mu) + B(\mu) u(t), \\ y(t, \mu) & = C(\mu) x(t, \mu) + D(\mu) u(t),\end{split}\]if continuous-time, or
\[\begin{split}E(\mu) x(k + 1, \mu) & = A(\mu) x(k, \mu) + B(\mu) u(k), \\ y(k, \mu) & = C(\mu) x(k, \mu) + D(\mu) u(k),\end{split}\]if discrete-time, where \(A\), \(B\), \(C\), \(D\), and \(E\) are linear operators.
Parameters
- A
The
OperatorA.- B
The
OperatorB.- C
The
OperatorC.- D
The
OperatorD orNone(then D is assumed to be zero).- E
The
OperatorE orNone(then E is assumed to be identity).- cont_time
Trueif the system is continuous-time, otherwiseFalse.- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Methods
eval_dtf,eval_tf,from_abcde_files,from_files,from_mat_file,from_matrices,gramian,h2_norm,hankel_norm,hinf_norm,hsv,poles,__add__,__mul__,__neg__,__sub__compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
ordercache_region,input_dim,output_dimparameters,parameters_inherited,parameters_internal,parameters_own,parametric-
order¶ The order of the system.
-
dim_input¶ The number of inputs.
-
dim_output¶ The number of outputs.
-
_hsv_U_V(mu=None)[source]¶ Compute Hankel singular values and vectors.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
- hsv
One-dimensional
NumPy arrayof singular values.- Uh
NumPy arrayof left singular vectors.- Vh
NumPy arrayof right singular vectors.
-
eval_dtf(s, mu=None)[source]¶ Evaluate the derivative of the transfer function.
The derivative of the transfer function at \(s\) is
\[-C(\mu) (s E(\mu) - A(\mu))^{-1} E(\mu) (s E(\mu) - A(\mu))^{-1} B(\mu).\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- dtfs
Derivative of transfer function evaluated at the complex number
s,NumPy arrayof shape(self.dim_output, self.dim_input).
-
eval_tf(s, mu=None)[source]¶ Evaluate the transfer function.
The transfer function at \(s\) is
\[C(\mu) (s E(\mu) - A(\mu))^{-1} B(\mu) + D(\mu).\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- tfs
Transfer function evaluated at the complex number
s,NumPy arrayof shape(self.dim_output, self.dim_input).
-
classmethod
from_abcde_files(files_basename, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶ Create
LTIModelfrom matrices stored in a .[ABCDE] files.Parameters
- files_basename
The basename of files containing A, B, C, and optionally D and E.
- cont_time
Trueif the system is continuous-time, otherwiseFalse.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Returns
- lti
The
LTIModelwith operators A, B, C, D, and E.
-
classmethod
from_files(A_file, B_file, C_file, D_file=None, E_file=None, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶ Create
LTIModelfrom matrices stored in separate files.Parameters
- A_file
The name of the file (with extension) containing A.
- B_file
The name of the file (with extension) containing B.
- C_file
The name of the file (with extension) containing C.
- D_file
Noneor the name of the file (with extension) containing D.- E_file
Noneor the name of the file (with extension) containing E.- cont_time
Trueif the system is continuous-time, otherwiseFalse.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Returns
- lti
The
LTIModelwith operators A, B, C, D, and E.
-
classmethod
from_mat_file(file_name, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶ Create
LTIModelfrom matrices stored in a .mat file.Parameters
- file_name
The name of the .mat file (extension .mat does not need to be included) containing A, B, C, and optionally D and E.
- cont_time
Trueif the system is continuous-time, otherwiseFalse.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Returns
- lti
The
LTIModelwith operators A, B, C, D, and E.
-
classmethod
from_matrices(A, B, C, D=None, E=None, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶ Create
LTIModelfrom matrices.Parameters
- A
The
NumPy arrayorSciPy spmatrixA.- B
The
NumPy arrayorSciPy spmatrixB.- C
The
NumPy arrayorSciPy spmatrixC.- D
The
NumPy arrayorSciPy spmatrixD orNone(then D is assumed to be zero).- E
The
NumPy arrayorSciPy spmatrixE orNone(then E is assumed to be identity).- cont_time
Trueif the system is continuous-time, otherwiseFalse.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Returns
- lti
The
LTIModelwith operators A, B, C, D, and E.
-
gramian(typ, mu=None)[source]¶ Compute a Gramian.
Parameters
- typ
The type of the Gramian:
'c_lrcf': low-rank Cholesky factor of the controllability Gramian,'o_lrcf': low-rank Cholesky factor of the observability Gramian,'c_dense': dense controllability Gramian,'o_dense': dense observability Gramian.
Note
For
'c_lrcf'and'o_lrcf'types, the method assumes the system is asymptotically stable. For'c_dense'and'o_dense'types, the method assumes there are no two system poles which add to zero.- mu
Returns
If typ is
'c_lrcf'or'o_lrcf', then the Gramian factor as aVectorArrayfromself.A.source. If typ is'c_dense'or'o_dense', then the Gramian as aNumPy array.
-
h2_norm(mu=None)[source]¶ Compute the H2-norm of the
LTIModel.Note
Assumes the system is asymptotically stable.
Parameters
Returns
- norm
H_2-norm.
-
hankel_norm(mu=None)[source]¶ Compute the Hankel-norm of the
LTIModel.Note
Assumes the system is asymptotically stable.
Parameters
Returns
- norm
Hankel-norm.
-
hinf_norm(mu=None, return_fpeak=False, ab13dd_equilibrate=False)[source]¶ Compute the H_infinity-norm of the
LTIModel.Note
Assumes the system is asymptotically stable.
Parameters
- mu
- return_fpeak
Whether to return the frequency at which the maximum is achieved.
- ab13dd_equilibrate
Whether
slycot.ab13ddshould use equilibration.
Returns
- norm
H_infinity-norm.
- fpeak
Frequency at which the maximum is achieved (if
return_fpeakisTrue).
-
hsv(mu=None)[source]¶ Hankel singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
- sv
One-dimensional
NumPy arrayof singular values.
-
poles(mu=None)[source]¶ Compute system poles.
Note
Assumes the systems is small enough to use a dense eigenvalue solver.
Parameters
- mu
Parameter valuesfor which to compute the systems poles.
Returns
One-dimensional
NumPy arrayof system poles.
-
class
pymor.models.iosys.LinearDelayModel(*args, **kwargs)[source]¶ Bases:
pymor.models.iosys.InputStateOutputModelClass for linear delay systems.
This class describes input-state-output systems given by
\[\begin{split}E x'(t) & = A x(t) + \sum_{i = 1}^q{A_i x(t - \tau_i)} + B u(t), \\ y(t) & = C x(t) + D u(t),\end{split}\]if continuous-time, or
\[\begin{split}E x(k + 1) & = A x(k) + \sum_{i = 1}^q{A_i x(k - \tau_i)} + B u(k), \\ y(k) & = C x(k) + D u(k),\end{split}\]if discrete-time, where \(E\), \(A\), \(A_i\), \(B\), \(C\), and \(D\) are linear operators.
Parameters
- A
The
OperatorA.- Ad
The tuple of
OperatorsA_i.- tau
The tuple of delay times (positive floats or ints).
- B
The
OperatorB.- C
The
OperatorC.- D
The
OperatorD orNone(then D is assumed to be zero).- E
The
OperatorE orNone(then E is assumed to be identity).- cont_time
Trueif the system is continuous-time, otherwiseFalse.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Methods
eval_dtf,eval_tf,__add__,__mul__,__neg__,__radd__,__rmul__,__rsub__,__sub__compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
ordercache_region,input_dim,output_dimparameters,parameters_inherited,parameters_internal,parameters_own,parametric-
order¶ The order of the system (equal to A.source.dim).
-
dim_input¶ The number of inputs.
-
dim_output¶ The number of outputs.
-
q¶ The number of delay terms.
-
tau¶ The tuple of delay times.
-
__add__(other)[source]¶ Add an
LTIModel,SecondOrderModelorLinearDelayModel.
-
__mul__(other)[source]¶ Postmultiply by a
SecondOrderModelor anLTIModel.
-
__neg__()[source]¶ Negate the
LinearDelayModel.
-
__radd__(other)[source]¶ Add to an
LTIModelor aSecondOrderModel.
-
__rmul__(other)[source]¶ Premultiply by an
LTIModelor aSecondOrderModel.
-
__rsub__(other)[source]¶ Subtract from an
LTIModelor aSecondOrderModel.
-
__sub__(other)[source]¶ Subtract an
LTIModel,SecondOrderModelorLinearDelayModel.
-
eval_dtf(s, mu=None)[source]¶ Evaluate the derivative of the transfer function.
The derivative of the transfer function at \(s\) is
\[-C \left(s E - A - \sum_{i = 1}^q{e^{-\tau_i s} A_i}\right)^{-1} \left(E + \sum_{i = 1}^q{\tau_i e^{-\tau_i s} A_i}\right) \left(s E - A - \sum_{i = 1}^q{e^{-\tau_i s} A_i}\right)^{-1} B.\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- dtfs
Derivative of transfer function evaluated at the complex number
s,NumPy arrayof shape(self.dim_output, self.dim_input).
-
eval_tf(s, mu=None)[source]¶ Evaluate the transfer function.
The transfer function at \(s\) is
\[C \left(s E - A - \sum_{i = 1}^q{e^{-\tau_i s} A_i}\right)^{-1} B + D.\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- tfs
Transfer function evaluated at the complex number
s,NumPy arrayof shape(self.dim_output, self.dim_input).
-
class
pymor.models.iosys.LinearStochasticModel(*args, **kwargs)[source]¶ Bases:
pymor.models.iosys.InputStateOutputModelClass for linear stochastic systems.
This class describes input-state-output systems given by
\[\begin{split}E \mathrm{d}x(t) & = A x(t) \mathrm{d}t + \sum_{i = 1}^q{A_i x(t) \mathrm{d}\omega_i(t)} + B u(t) \mathrm{d}t, \\ y(t) & = C x(t) + D u(t),\end{split}\]if continuous-time, or
\[\begin{split}E x(k + 1) & = A x(k) + \sum_{i = 1}^q{A_i x(k) \omega_i(k)} + B u(k), \\ y(k) & = C x(k) + D u(t),\end{split}\]if discrete-time, where \(E\), \(A\), \(A_i\), \(B\), \(C\), and \(D\) are linear operators and \(\omega_i\) are stochastic processes.
Parameters
- A
The
OperatorA.- As
The tuple of
OperatorsA_i.- B
The
OperatorB.- C
The
OperatorC.- D
The
OperatorD orNone(then D is assumed to be zero).- E
The
OperatorE orNone(then E is assumed to be identity).- cont_time
Trueif the system is continuous-time, otherwiseFalse.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
ordercache_region,input_dim,output_dimparameters,parameters_inherited,parameters_internal,parameters_own,parametric-
order¶ The order of the system (equal to A.source.dim).
-
dim_input¶ The number of inputs.
-
dim_output¶ The number of outputs.
-
q¶ The number of stochastic processes.
-
class
pymor.models.iosys.SecondOrderModel(*args, **kwargs)[source]¶ Bases:
pymor.models.iosys.InputStateOutputModelClass for linear second order systems.
This class describes input-output systems given by
\[\begin{split}M(\mu) \ddot{x}(t, \mu) + E(\mu) \dot{x}(t, \mu) + K(\mu) x(t, \mu) & = B(\mu) u(t), \\ y(t, \mu) & = C_p(\mu) x(t, \mu) + C_v(\mu) \dot{x}(t, \mu) + D(\mu) u(t),\end{split}\]if continuous-time, or
\[\begin{split}M(\mu) x(k + 2, \mu) + E(\mu) x(k + 1, \mu) + K(\mu) x(k, \mu) & = B(\mu) u(k), \\ y(k, \mu) & = C_p(\mu) x(k, \mu) + C_v(\mu) x(k + 1, \mu) + D(\mu) u(k),\end{split}\]if discrete-time, where \(M\), \(E\), \(K\), \(B\), \(C_p\), \(C_v\), and \(D\) are linear operators.
Parameters
- M
The
OperatorM.- E
The
OperatorE.- K
The
OperatorK.- B
The
OperatorB.- Cp
The
OperatorCp.- Cv
The
OperatorCv orNone(then Cv is assumed to be zero).- D
The
OperatorD orNone(then D is assumed to be zero).- cont_time
Trueif the system is continuous-time, otherwiseFalse.- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Methods
eval_dtf,eval_tf,from_files,from_matrices,gramian,h2_norm,hankel_norm,hinf_norm,poles,psv,pvsv,to_lti,vpsv,vsv,__add__,__mul__,__neg__,__radd__,__rmul__,__rsub__,__sub__compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
ordercache_region,input_dim,output_dimparameters,parameters_inherited,parameters_internal,parameters_own,parametric-
order¶ The order of the system (equal to M.source.dim).
-
dim_input¶ The number of inputs.
-
dim_output¶ The number of outputs.
-
__add__(other)[source]¶ Add a
SecondOrderModelor anLTIModel.
-
__mul__(other)[source]¶ Postmultiply by a
SecondOrderModelor anLTIModel.
-
__neg__()[source]¶ Negate the
SecondOrderModel.
-
__sub__(other)[source]¶ Subtract a
SecondOrderModelor anLTIModel.
-
eval_dtf(s, mu=None)[source]¶ Evaluate the derivative of the transfer function.
\[s C_v(\mu) (s^2 M(\mu) + s E(\mu) + K(\mu))^{-1} B(\mu) - (C_p(\mu) + s C_v(\mu)) (s^2 M(\mu) + s E(\mu) + K(\mu))^{-1} (2 s M(\mu) + E(\mu)) (s^2 M(\mu) + s E(\mu) + K(\mu))^{-1} B(\mu).\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- dtfs
Derivative of transfer function evaluated at the complex number
s,NumPy arrayof shape(self.dim_output, self.dim_input).
-
eval_tf(s, mu=None)[source]¶ Evaluate the transfer function.
The transfer function at \(s\) is
\[(C_p(\mu) + s C_v(\mu)) (s^2 M(\mu) + s E(\mu) + K(\mu))^{-1} B(\mu) + D(\mu).\]Note
Assumes that the number of inputs and outputs is much smaller than the order of the system.
Parameters
- s
Complex number.
- mu
Returns
- tfs
Transfer function evaluated at the complex number
s,NumPy arrayof shape(self.dim_output, self.dim_input).
-
classmethod
from_files(M_file, E_file, K_file, B_file, Cp_file, Cv_file=None, D_file=None, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶ Create
LTIModelfrom matrices stored in separate files.Parameters
- M_file
The name of the file (with extension) containing A.
- E_file
The name of the file (with extension) containing E.
- K_file
The name of the file (with extension) containing K.
- B_file
The name of the file (with extension) containing B.
- Cp_file
The name of the file (with extension) containing Cp.
- Cv_file
Noneor the name of the file (with extension) containing Cv.- D_file
Noneor the name of the file (with extension) containing D.- cont_time
Trueif the system is continuous-time, otherwiseFalse.- state_id
Id of the state space.
- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Returns
- som
The
SecondOrderModelwith operators M, E, K, B, Cp, Cv, and D.
-
classmethod
from_matrices(M, E, K, B, Cp, Cv=None, D=None, cont_time=True, state_id='STATE', solver_options=None, error_estimator=None, visualizer=None, name=None)[source]¶ Create a second order system from matrices.
Parameters
- M
The
NumPy arrayorSciPy spmatrixM.- E
The
NumPy arrayorSciPy spmatrixE.- K
The
NumPy arrayorSciPy spmatrixK.- B
The
NumPy arrayorSciPy spmatrixB.- Cp
The
NumPy arrayorSciPy spmatrixCp.- Cv
The
NumPy arrayorSciPy spmatrixCv orNone(then Cv is assumed to be zero).- D
The
NumPy arrayorSciPy spmatrixD orNone(then D is assumed to be zero).- cont_time
Trueif the system is continuous-time, otherwiseFalse.- solver_options
The solver options to use to solve the Lyapunov equations.
- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, model)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, model, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the system.
Returns
- lti
The SecondOrderModel with operators M, E, K, B, Cp, Cv, and D.
-
gramian(typ, mu=None)[source]¶ Compute a second-order Gramian.
Parameters
- typ
The type of the Gramian:
'pc_lrcf': low-rank Cholesky factor of the position controllability Gramian,'vc_lrcf': low-rank Cholesky factor of the velocity controllability Gramian,'po_lrcf': low-rank Cholesky factor of the position observability Gramian,'vo_lrcf': low-rank Cholesky factor of the velocity observability Gramian,'pc_dense': dense position controllability Gramian,'vc_dense': dense velocity controllability Gramian,'po_dense': dense position observability Gramian,'vo_dense': dense velocity observability Gramian.
Note
For
'*_lrcf'types, the method assumes the system is asymptotically stable. For'*_dense'types, the method assumes there are no two system poles which add to zero.- mu
Returns
If typ is
'pc_lrcf','vc_lrcf','po_lrcf'or'vo_lrcf', then the Gramian factor as aVectorArrayfromself.M.source. If typ is'pc_dense','vc_dense','po_dense'or'vo_dense', then the Gramian as aNumPy array.
-
h2_norm(mu=None)[source]¶ Compute the H2-norm.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
- norm
H_2-norm.
-
hankel_norm(mu=None)[source]¶ Compute the Hankel-norm.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
- norm
Hankel-norm.
-
hinf_norm(mu=None, return_fpeak=False, ab13dd_equilibrate=False)[source]¶ Compute the H_infinity-norm.
Note
Assumes the system is asymptotically stable.
Parameters
- mu
- return_fpeak
Should the frequency at which the maximum is achieved should be returned.
- ab13dd_equilibrate
Should
slycot.ab13dduse equilibration.
Returns
- norm
H_infinity-norm.
- fpeak
Frequency at which the maximum is achieved (if
return_fpeakisTrue).
-
poles(mu=None)[source]¶ Compute system poles.
Note
Assumes the systems is small enough to use a dense eigenvalue solver.
Parameters
Returns
One-dimensional
NumPy arrayof system poles.
-
psv(mu=None)[source]¶ Position singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
One-dimensional
NumPy arrayof singular values.
-
pvsv(mu=None)[source]¶ Position-velocity singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
One-dimensional
NumPy arrayof singular values.
-
to_lti()[source]¶ Return a first order representation.
The first order representation
\[\begin{split}\begin{bmatrix} I & 0 \\ 0 & M \end{bmatrix} \frac{\mathrm{d}}{\mathrm{d}t}\! \begin{bmatrix} x(t) \\ \dot{x}(t) \end{bmatrix} & = \begin{bmatrix} 0 & I \\ -K & -E \end{bmatrix} \begin{bmatrix} x(t) \\ \dot{x}(t) \end{bmatrix} + \begin{bmatrix} 0 \\ B \end{bmatrix} u(t), \\ y(t) & = \begin{bmatrix} C_p & C_v \end{bmatrix} \begin{bmatrix} x(t) \\ \dot{x}(t) \end{bmatrix} + D u(t)\end{split}\]is returned.
Returns
- lti
LTIModelequivalent to the second-order model.
-
vpsv(mu=None)[source]¶ Velocity-position singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
One-dimensional
NumPy arrayof singular values.
-
vsv(mu=None)[source]¶ Velocity singular values.
Note
Assumes the system is asymptotically stable.
Parameters
Returns
One-dimensional
NumPy arrayof singular values.
-
class
pymor.models.iosys.TransferFunction(*args, **kwargs)[source]¶ Bases:
pymor.models.iosys.InputOutputModelClass for systems represented by a transfer function.
This class describes input-output systems given by a transfer function \(H(s, mu)\).
Parameters
- dim_input
The number of inputs.
- dim_output
The number of outputs.
- tf
The transfer function defined at least on the open right complex half-plane.
tf(s, mu)is aNumPy arrayof shape(p, m).- dtf
The complex derivative of
Hwith respect tos.- cont_time
Trueif the system is continuous-time, otherwiseFalse.- name
Name of the system.
Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
cache_region,input_dim,output_dimparameters,parameters_inherited,parameters_internal,parameters_own,parametric-
dim_input¶ The number of inputs.
-
dim_output¶ The number of outputs.
-
tf¶ The transfer function.
-
dtf¶ The complex derivative of the transfer function.
-
h2_norm(return_norm_only=True, **quad_kwargs)[source]¶ Compute the H2-norm using quadrature.
This method uses
scipy.integrate.quadand makes no assumptions on the form of the transfer function.By default, the absolute error tolerance in
scipy.integrate.quadis set to zero (see its optional argumentepsabs). It can be changed by using theepsabskeyword argument.Parameters
- return_norm_only
Whether to only return the approximate H2-norm.
- quad_kwargs
Keyword arguments passed to
scipy.integrate.quad.
Returns
- norm
Computed H2-norm.
- norm_relerr
Relative error estimate (returned if
return_norm_onlyisFalse).- info
Quadrature info (returned if
return_norm_onlyisFalseandfull_outputisTrue). Seescipy.integrate.quaddocumentation for more details.
-
pymor.models.iosys.sparse_min_size(value=1000)[source]¶ Return minimal sparse problem size for which to warn about converting to dense.
Defaults
value (see
pymor.core.defaults)
mpi module¶
-
class
pymor.models.mpi.MPIModel(obj_id, *args, **kwargs)[source]¶ Bases:
objectWrapper class mixin for MPI distributed
Models.Given a single-rank implementation of a
Model, this wrapper class uses the event loop frompymor.tools.mpito allow an MPI distributed usage of theModel. The underlying implementation needs to be MPI aware. In particular, the model’ssolvemethod has to perform an MPI parallel solve of the model.Note that this class is not intended to be instantiated directly. Instead, you should use
mpi_wrap_model.Methods
visualize
-
class
pymor.models.mpi.MPIVisualizer(*args, **kwargs)[source]¶ Bases:
pymor.core.base.ImmutableObjectMethods
visualizeAttributes
-
class
pymor.models.mpi._OperatorToWrap(operator, mpi_range, mpi_source)¶ Bases:
tupleMethods
tuplecount,indexAttributes
-
__getnewargs__()¶ Return self as a plain tuple. Used by copy and pickle.
-
static
__new__(_cls, operator, mpi_range, mpi_source)¶ Create new instance of _OperatorToWrap(operator, mpi_range, mpi_source)
-
__repr__()¶ Return a nicely formatted representation string
-
_asdict()¶ Return a new OrderedDict which maps field names to their values.
-
classmethod
_make(iterable)¶ Make a new _OperatorToWrap object from a sequence or iterable
-
_replace(**kwds)¶ Return a new _OperatorToWrap object replacing specified fields with new values
-
property
mpi_range¶ Alias for field number 1
-
property
mpi_source¶ Alias for field number 2
-
property
operator¶ Alias for field number 0
-
-
pymor.models.mpi.mpi_wrap_model(local_models, mpi_spaces=('STATE', ), use_with=True, with_apply2=False, pickle_local_spaces=True, space_type=<class 'pymor.vectorarrays.mpi.MPIVectorSpace'>, base_type=None)[source]¶ Wrap MPI distributed local
Modelsto a globalModelon rank 0.Given MPI distributed local
Modelsreferred to by theObjectIdlocal_models, return a newModelwhich manages these distributed models from rank 0. This is done by first wrapping allOperatorsof theModelusingmpi_wrap_operator.Alternatively,
local_modelscan be a callable (with no arguments) which is then called on each rank to instantiate the localModels.When
use_withisFalse, anMPIModelis instantiated with the wrapped operators. A call tosolvewill then use an MPI parallel call to thesolvemethods of the wrapped localModelsto obtain the solution. This is usually what you want when the actual solve is performed by an implementation in the external solver.When
use_withisTrue,with_is called on the localModelon rank 0, to obtain a newModelwith the wrapped MPIOperators. This is mainly useful when the local models are genericModelsas inpymor.models.basicandsolveis implemented directly in pyMOR via operations on the containedOperators.Parameters
- local_models
ObjectIdof the localModelson each rank or a callable generating theModels.- mpi_spaces
List of types or ids of
VectorSpaceswhich are MPI distributed and need to be wrapped.- use_with
See above.
- with_apply2
See
MPIOperator.- pickle_local_spaces
See
MPIOperator.- space_type
See
MPIOperator.
neural_network module¶
-
class
pymor.models.neural_network.FullyConnectedNN(layers_sizes, activation_function=<built-in method tanh of type object>)[source]¶ Bases:
torch.nn.modules.module.Module,pymor.core.base.BasicObjectClass for neural networks with fully connected layers.
This class implements neural networks consisting of linear and fully connected layers. Furthermore, the same activation function is used between each layer, except for the last one where no activation function is applied.
Parameters
- layers_sizes
List of sizes (i.e. number of neurons) for the layers of the neural network.
- activation_function
Function to use as activation function between the single layers.
Methods
Moduleadd_module,apply,bfloat16,buffers,children,cpu,cuda,double,eval,extra_repr,float,half,load_state_dict,modules,named_buffers,named_children,named_modules,named_parameters,parameters,register_backward_hook,register_buffer,register_forward_hook,register_forward_pre_hook,register_parameter,requires_grad_,share_memory,state_dict,to,train,type,zero_gradAttributes
Moduledump_patches,T_destination
-
class
pymor.models.neural_network.NeuralNetworkInstationaryModel(*args, **kwargs)[source]¶ Bases:
pymor.models.interface.ModelClass for models of instationary problems that use artificial neural networks.
This class implements a
Modelthat uses a neural network for solving.Parameters
- T
The final time T.
- nt
The number of time steps.
- neural_network
The neural network that approximates the mapping from parameter space to solution space. Should be an instance of
FullyConnectedNNwith input size that matches the (total) number of parameters and output size equal to the dimension of the reduced space.- parameters
Parametersof the reduced order model (the same as used in the full-order model).- output_functional
Operatormapping a given solution to the model output. In many applications, this will be aFunctional, i.e. anOperatormapping to scalars. This is not required, however.- products
A dict of inner product
Operatorsdefined on the discrete space the problem is posed on. For each product with key'x'a corresponding attributex_product, as well as a norm methodx_normis added to the model.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, m)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, m, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the model.
Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
_compute_solution(mu=None, **kwargs)[source]¶ Compute the model’s solution for
parameter valuesmu.This method is called by the default implementation of
computeinpymor.models.interface.Model.Parameters
- mu
Parameter valuesfor which to compute the solution.- kwargs
Additional keyword arguments to customize how the solution is computed or to select additional data to be returned.
Returns
VectorArraywith the computed solution or a dict which at least must contain the key'solution'.
-
class
pymor.models.neural_network.NeuralNetworkModel(*args, **kwargs)[source]¶ Bases:
pymor.models.interface.ModelClass for models of stationary problems that use artificial neural networks.
This class implements a
Modelthat uses a neural network for solving.Parameters
- neural_network
The neural network that approximates the mapping from parameter space to solution space. Should be an instance of
FullyConnectedNNwith input size that matches the (total) number of parameters and output size equal to the dimension of the reduced space.- parameters
Parametersof the reduced order model (the same as used in the full-order model).- output_functional
Operatormapping a given solution to the model output. In many applications, this will be aFunctional, i.e. anOperatormapping to scalars. This is not required, however.- products
A dict of inner product
Operatorsdefined on the discrete space the problem is posed on. For each product with key'x'a corresponding attributex_product, as well as a norm methodx_normis added to the model.- error_estimator
An error estimator for the problem. This can be any object with an
estimate_error(U, mu, m)method. Iferror_estimatoris notNone, anestimate_error(U, mu)method is added to the model which will callerror_estimator.estimate_error(U, mu, self).- visualizer
A visualizer for the problem. This can be any object with a
visualize(U, m, ...)method. Ifvisualizeris notNone, avisualize(U, *args, **kwargs)method is added to the model which forwards its arguments to the visualizer’svisualizemethod.- name
Name of the model.
Methods
compute,estimate,estimate_error,estimate_output_error,output,output_d_mu,solve,solve_d_mu,visualizeAttributes
parameters,parameters_inherited,parameters_internal,parameters_own,parametric-
_compute_solution(mu=None, **kwargs)[source]¶ Compute the model’s solution for
parameter valuesmu.This method is called by the default implementation of
computeinpymor.models.interface.Model.Parameters
- mu
Parameter valuesfor which to compute the solution.- kwargs
Additional keyword arguments to customize how the solution is computed or to select additional data to be returned.
Returns
VectorArraywith the computed solution or a dict which at least must contain the key'solution'.