Struct varpro::solvers::levmar::LevMarProblemBuilder
source · pub struct LevMarProblemBuilder<Model, const MRHS: bool>where
Model::ScalarType: Scalar + ComplexField + Copy,
<Model::ScalarType as ComplexField>::RealField: Float,
Model: SeparableNonlinearModel,{ /* private fields */ }
Expand description
A builder structure to create a LevMarProblem, which can be used for fitting a separable model to data.
§Example
The following code shows how to create an unweighted least squares problem to fit the separable model
$\vec{f}(\vec{x},\vec{\alpha},\vec{c})$
given by model
to data $\vec{y}$
(given as y
) using the
independent variable $\vec{x}$
(given as x
). Furthermore we need an initial guess params
for the nonlinear model parameters $\vec{\alpha}$
let problem = LevMarProblemBuilder::new(model)
.observations(y)
.build()
.unwrap();
§Building a Model
A new builder is constructed with the new constructor. It must be filled with content using the methods described in the following. After all mandatory fields have been filled, the build method can be called. This returns a Result type that contains the finished model iff all mandatory fields have been set with valid values. Otherwise it contains an error variant.
§Multiple Right Hand Sides
We can also construct a problem with multiple right hand sides, using the mrhs constructor, see LevMarProblem for additional details.
Implementations§
source§impl<Model> LevMarProblemBuilder<Model, false>where
Model::ScalarType: Scalar + ComplexField + Zero + Copy,
<Model::ScalarType as ComplexField>::RealField: Float,
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Mul<Model::ScalarType, Output = Model::ScalarType> + Float,
Model: SeparableNonlinearModel,
impl<Model> LevMarProblemBuilder<Model, false>where
Model::ScalarType: Scalar + ComplexField + Zero + Copy,
<Model::ScalarType as ComplexField>::RealField: Float,
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Mul<Model::ScalarType, Output = Model::ScalarType> + Float,
Model: SeparableNonlinearModel,
sourcepub fn new(model: Model) -> Self
pub fn new(model: Model) -> Self
Create a new builder based on the given model for a problem with a single right hand side. This is the standard use case, where the data is a vector that is fitted by the model.
sourcepub fn observations(self, observed: OVector<Model::ScalarType, Dyn>) -> Self
pub fn observations(self, observed: OVector<Model::ScalarType, Dyn>) -> Self
Mandatory: Set the data which we want to fit: since this
is called on a model builder for problems with single right hand sides,
this is a column vector $\vec{y}=\vec{y}(\vec{x})$
containing the
values we want to fit with the model.
The length of $\vec{y}
and the output dimension of the model must be
the same.
source§impl<Model> LevMarProblemBuilder<Model, true>where
Model::ScalarType: Scalar + ComplexField + Zero + Copy,
<Model::ScalarType as ComplexField>::RealField: Float,
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Mul<Model::ScalarType, Output = Model::ScalarType> + Float,
Model: SeparableNonlinearModel,
impl<Model> LevMarProblemBuilder<Model, true>where
Model::ScalarType: Scalar + ComplexField + Zero + Copy,
<Model::ScalarType as ComplexField>::RealField: Float,
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Mul<Model::ScalarType, Output = Model::ScalarType> + Float,
Model: SeparableNonlinearModel,
sourcepub fn mrhs(model: Model) -> Self
pub fn mrhs(model: Model) -> Self
Create a new builder based on the given model for a problem with multiple right hand sides and perform a global fit.
That means the observations are expected to be a matrix, where the columns correspond to the individual observations. The nonlinear parameters will be optimized across all the observations, but the linear coefficients are calculated for each observation individually. Hence, they also become a matrix where the columns correspond to the linear coefficients of the observation in the same column.
sourcepub fn observations(
self,
observed: OMatrix<Model::ScalarType, Dyn, Dyn>
) -> Self
pub fn observations( self, observed: OMatrix<Model::ScalarType, Dyn, Dyn> ) -> Self
Mandatory: Set the data which we want to fit: This is either a single vector
$\vec{y}=\vec{y}(\vec{x})$
or a matrix $\boldsymbol{Y}$
of multiple
vectors. In the former case this corresponds to fitting a single right hand side,
in the latter case, this corresponds to global fitting of a problem with
multiple right hand sides.
The length of $\vec{x}$
and the number of rows in the data must
be the same.
source§impl<Model, const MRHS: bool> LevMarProblemBuilder<Model, MRHS>where
Model::ScalarType: Scalar + ComplexField + Zero + Copy,
<Model::ScalarType as ComplexField>::RealField: Float,
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Mul<Model::ScalarType, Output = Model::ScalarType> + Float,
Model: SeparableNonlinearModel,
impl<Model, const MRHS: bool> LevMarProblemBuilder<Model, MRHS>where
Model::ScalarType: Scalar + ComplexField + Zero + Copy,
<Model::ScalarType as ComplexField>::RealField: Float,
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Mul<Model::ScalarType, Output = Model::ScalarType> + Float,
Model: SeparableNonlinearModel,
sourcepub fn epsilon(
self,
eps: <Model::ScalarType as ComplexField>::RealField
) -> Self
pub fn epsilon( self, eps: <Model::ScalarType as ComplexField>::RealField ) -> Self
Optional This value is relevant for the solver, because it uses singular value decomposition
internally. This method sets a value \epsilon
for which smaller (i.e. absolute - wise) singular
values are considered zero. In essence this gives a truncation of the SVD. This might be
helpful if two basis functions become linear dependent when the nonlinear model parameters
align in an unfortunate way. In this case a higher epsilon might increase the robustness
of the fitting process.
If this value is not given, it will be set to machine epsilon.
The given epsilon is automatically converted to a non-negative number.
sourcepub fn weights(self, weights: OVector<Model::ScalarType, Dyn>) -> Self
pub fn weights(self, weights: OVector<Model::ScalarType, Dyn>) -> Self
Optional Add diagonal weights to the problem (meaning data points are statistically independent). If this is not given, the problem is unweighted, i.e. each data point has unit weight.
Note The weighted residual is calculated as $||W(\vec{y}-\vec{f}(\vec{\alpha})||^2$
, so
to make weights that have a statistical meaning, the diagonal elements of the weight matrix should be
set to w_{jj} = 1/\sigma_j
where $\sigma_j$
is the (estimated) standard deviation associated with
data point $y_j$
.
sourcepub fn build(self) -> Result<LevMarProblem<Model, MRHS>, LevMarBuilderError>
pub fn build(self) -> Result<LevMarProblem<Model, MRHS>, LevMarBuilderError>
build the least squares problem from the builder.
§Prerequisites
- All mandatory parameters have been set (see individual builder methods for details)
$\vec{x}$
and$\vec{y}$
have the same number of elements$\vec{x}$
and$\vec{y}$
have a nonzero number of elements- the length of the initial guesses vector is the same as the number of model parameters
§Returns
If all prerequisites are fulfilled, returns a LevMarProblem with the given content and the parameters set to the initial guess. Otherwise returns an error variant.
Trait Implementations§
source§impl<Model, const MRHS: bool> Clone for LevMarProblemBuilder<Model, MRHS>where
Model::ScalarType: Scalar + ComplexField + Copy + Clone,
<Model::ScalarType as ComplexField>::RealField: Float,
Model: SeparableNonlinearModel + Clone,
impl<Model, const MRHS: bool> Clone for LevMarProblemBuilder<Model, MRHS>where
Model::ScalarType: Scalar + ComplexField + Copy + Clone,
<Model::ScalarType as ComplexField>::RealField: Float,
Model: SeparableNonlinearModel + Clone,
source§fn clone(&self) -> LevMarProblemBuilder<Model, MRHS>
fn clone(&self) -> LevMarProblemBuilder<Model, MRHS>
1.0.0 · source§fn clone_from(&mut self, source: &Self)
fn clone_from(&mut self, source: &Self)
source
. Read moreAuto Trait Implementations§
impl<Model, const MRHS: bool> Freeze for LevMarProblemBuilder<Model, MRHS>where
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Sized + Freeze,
<Model as SeparableNonlinearModel>::ScalarType: Sized,
Model: Freeze,
impl<Model, const MRHS: bool> RefUnwindSafe for LevMarProblemBuilder<Model, MRHS>where
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Sized + RefUnwindSafe,
<Model as SeparableNonlinearModel>::ScalarType: Sized + RefUnwindSafe,
Model: RefUnwindSafe,
impl<Model, const MRHS: bool> Send for LevMarProblemBuilder<Model, MRHS>where
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Sized,
<Model as SeparableNonlinearModel>::ScalarType: Sized,
Model: Send,
impl<Model, const MRHS: bool> Sync for LevMarProblemBuilder<Model, MRHS>where
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Sized,
<Model as SeparableNonlinearModel>::ScalarType: Sized,
Model: Sync,
impl<Model, const MRHS: bool> Unpin for LevMarProblemBuilder<Model, MRHS>where
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Sized + Unpin,
<Model as SeparableNonlinearModel>::ScalarType: Sized + Unpin,
Model: Unpin,
impl<Model, const MRHS: bool> UnwindSafe for LevMarProblemBuilder<Model, MRHS>where
<<Model as SeparableNonlinearModel>::ScalarType as ComplexField>::RealField: Sized + UnwindSafe,
<Model as SeparableNonlinearModel>::ScalarType: Sized + UnwindSafe,
Model: UnwindSafe,
Blanket Implementations§
source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
source§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
source§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self
from the equivalent element of its
superset. Read moresource§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
self
is actually part of its subset T
(and can be converted to it).source§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset
but without any property checks. Always succeeds.source§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self
to the equivalent element of its superset.