svZeroDSolver
Loading...
Searching...
No Matches
LevenbergMarquardtOptimizer.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Stanford University, The Regents of the
2// University of California, and others. SPDX-License-Identifier: BSD-3-Clause
3/**
4 * @file LevenbergMarquardtOptimizer.h
5 * @brief opt::LevenbergMarquardtOptimizer source file
6 */
7#ifndef SVZERODSOLVER_OPTIMIZE_LEVENBERGMARQUARDT_HPP_
8#define SVZERODSOLVER_OPTIMIZE_LEVENBERGMARQUARDT_HPP_
9
10#include <Eigen/Dense>
11#include <Eigen/Sparse>
12
13#include "Model.h"
14
15/**
16 * @brief Levenberg-Marquardt optimization class
17 *
18 * The 0D residual (assuming no time-dependency in parameters) is
19 *
20 * \f[
21 * \boldsymbol{r}(\boldsymbol{\alpha}, \boldsymbol{y}, \boldsymbol{\dot{y}}) =
22 * \boldsymbol{E}(\boldsymbol{\alpha}, \boldsymbol{y}) \cdot
23 * \dot{\boldsymbol{y}}+\boldsymbol{F}(\boldsymbol{\alpha}, \boldsymbol{y})
24 * \cdot \boldsymbol{y}+\boldsymbol{c}(\boldsymbol{\alpha}, \boldsymbol{y}) \f]
25 *
26 * with solution vector \f$\boldsymbol{y} \in \mathbb{R}^{N}\f$ (flow and
27 * pressure at nodes), LPN parameters \f$\boldsymbol{\alpha} \in
28 * \mathbb{R}^{P}\f$, system matrices \f$\boldsymbol{E},\boldsymbol{F} \in
29 * \mathbb{R}^{NxN}\f$, and system vector \f$\boldsymbol{c} \in
30 * \mathbb{R}^{N}\f$.
31 *
32 * The least squares problem can be formulated as
33 *
34 * \f[
35 * \min _\alpha S, \quad \mathrm { with } \quad S=\sum_i^D
36 * r_i^2\left(\boldsymbol{\alpha}, y_i, \dot{y}_i\right) \f]
37 *
38 * with given solution vectors \f$\boldsymbol{y}\f$, \f$\boldsymbol{\dot{y}}\f$
39 * at all datapoints \f$D\f$. The parameter vector is iteratively improved
40 * according to
41 *
42 * \f[
43 * \boldsymbol{\alpha}^{i+1}=\boldsymbol{\alpha}^{i}+\Delta
44 * \boldsymbol{\alpha}^{i+1} \f]
45 *
46 * wherein the increment \f$\Delta \boldsymbol{\alpha}^{i+1} \f$ is determined
47 * by solving the following system:
48 *
49 * \f[
50 * \left[\mathbf{J}^{\mathrm{T}} \mathbf{J}+\lambda
51 * \operatorname{diag}\left(\mathbf{J}^{\mathrm{T}} \mathbf{J}\right)\right]^{i}
52 * \cdot \Delta \boldsymbol{\alpha}^{i+1}=-\left[\mathbf{J}^{\mathrm{T}}
53 * \mathbf{r}\right]^{i}, \quad \lambda^{i}=\lambda^{i-1}
54 * \cdot\left\|\left[\mathbf{J}^{\mathrm{T}} \mathbf{r}\right]^{i}\right\|_2
55 * /\left\|\left[\mathbf{J}^{\mathrm{T}} \mathbf{r}\right]^{i-1}\right\|_2. \f]
56 *
57 * The algorithm terminates when the following tolerance thresholds are reached
58 *
59 * \f[
60 * \left\|\left[\mathbf{J}^{\mathrm{T}}
61 * \mathbf{r}\right]^{\mathrm{i}}\right\|_2<\operatorname{tol}_{\text {grad
62 * }}^\alpha \text { and }\left\|\Delta
63 * \boldsymbol{\alpha}^{\mathrm{i}+1}\right\|_2<\mathrm{tol}_{\text {inc
64 * }}^\alpha, \f]
65 *
66 * The Jacobian is derived from the residual as
67 *
68 * \f[
69 * J = \frac{\partial \boldsymbol{r}}{\partial \boldsymbol{\alpha}} =
70 * \frac{\partial \mathbf{E}}{\partial \boldsymbol{\alpha}} \cdot
71 * \dot{\mathbf{y}}+\frac{\partial \mathbf{F}}{\partial \boldsymbol{\alpha}}
72 * \cdot \mathbf{y}+\frac{\partial \mathbf{c}}{\partial \boldsymbol{\alpha}} \f]
73 *
74 *
75 */
77 public:
78 /**
79 * @brief Construct a new LevenbergMarquardtOptimizer object
80 *
81 * @param model The 0D model
82 * @param num_obs Number of observations in optimization
83 * @param num_params Number of parameters in optimization
84 * @param lambda0 Initial damping factor
85 * @param tol_grad Gradient tolerance
86 * @param tol_inc Parameter increment tolerance
87 * @param max_iter Maximum iterations
88 */
89 LevenbergMarquardtOptimizer(Model* model, int num_obs, int num_params,
90 double lambda0, double tol_grad, double tol_inc,
91 int max_iter);
92
93 /**
94 * @brief Run the optimization algorithm
95 *
96 * @param alpha Initial parameter vector alpha
97 * @param y_obs Matrix (num_obs x n) with all observations for y
98 * @param dy_obs Matrix (num_obs x n) with all observations for dy
99 * @return Eigen::Matrix<double, Eigen::Dynamic, 1> Optimized parameter vector
100 * alpha
101 */
102 Eigen::Matrix<double, Eigen::Dynamic, 1> run(
103 Eigen::Matrix<double, Eigen::Dynamic, 1> alpha,
104 std::vector<std::vector<double>>& y_obs,
105 std::vector<std::vector<double>>& dy_obs);
106
107 private:
108 Eigen::SparseMatrix<double> jacobian;
109 Eigen::Matrix<double, Eigen::Dynamic, 1> residual;
110 Eigen::Matrix<double, Eigen::Dynamic, 1> delta;
111 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> mat;
112 Eigen::Matrix<double, Eigen::Dynamic, 1> vec;
113 Model* model;
114 double lambda;
115
116 int num_obs;
117 int num_params;
118 int num_eqns;
119 int num_vars;
120 int num_dpoints;
121
122 double tol_grad;
123 double tol_inc;
124 int max_iter;
125
126 void update_gradient(Eigen::Matrix<double, Eigen::Dynamic, 1>& alpha,
127 std::vector<std::vector<double>>& y_obs,
128 std::vector<std::vector<double>>& dy_obs);
129
130 void update_delta(bool first_step);
131};
132
133#endif // SVZERODSOLVER_OPTIMIZE_LEVENBERGMARQUARDT_HPP_
model::Model source file
LevenbergMarquardtOptimizer(Model *model, int num_obs, int num_params, double lambda0, double tol_grad, double tol_inc, int max_iter)
Construct a new LevenbergMarquardtOptimizer object.
Definition LevenbergMarquardtOptimizer.cpp:7
Eigen::Matrix< double, Eigen::Dynamic, 1 > run(Eigen::Matrix< double, Eigen::Dynamic, 1 > alpha, std::vector< std::vector< double > > &y_obs, std::vector< std::vector< double > > &dy_obs)
Run the optimization algorithm.
Definition LevenbergMarquardtOptimizer.cpp:28
Model of 0D elements.
Definition Model.h:48