Loading...
Searching...
No Matches
ReducedProblem.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
33
34#include "ReducedProblem.H"
35
36// ******************** //
37// class reducedProblem //
38// ******************** //
39
40// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
41
42// Constructor
43
47
49 :
50 problem(&problem)
51{
52}
53
54// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55
57{
58 Info << "The method reducedProblem::solveOnline in reducedProblem.C is a virtual method"
59 << endl;
60 Info << "It must be overridden, exiting the code" << endl;
61 exit(0);
62}
63
64Eigen::MatrixXd reducedProblem::solveLinearSys(List<Eigen::MatrixXd> LinSys,
65 Eigen::MatrixXd x, Eigen::VectorXd& residual, const Eigen::MatrixXd& bc,
66 const std::string solverType)
67
68{
69 M_Assert(solverType == "fullPivLu" || solverType == "partialPivLu"
70 || solverType == "householderQR" || solverType == "colPivHouseholderQR"
71 || solverType == "fullPivHouseholderQR"
72 || solverType == "completeOrthogonalDecomposition" || solverType == "llt"
73 || solverType == "ldlt" || solverType == "bdcSvd"
74 || solverType == "jacobiSvd", "solver not defined");
75 Eigen::MatrixXd A;
76 Eigen::MatrixXd b;
77
78 if (LinSys[0].rows() < LinSys[0].cols())
79 {
80 FatalErrorInFunction << "The system is undetermined, it has more unknowns: " <<
81 LinSys[0].rows() << ", then equations: " << LinSys[0].rows() << abort(
82 FatalError);
83 }
84
85 Eigen::MatrixXd y;
86 residual = LinSys[0] * x - LinSys[1];
87 A = LinSys[0];
88 b = LinSys[1];
89
90 if (A.rows() > A.cols())
91 {
92 if (solverType != "completeOrthogonalDecomposition" && solverType != "bdcSvd"
93 && solverType != "jacobiSvd")
94 {
95 A.resize(A.cols(), A.cols());
96 b.resize(A.cols(), 1);
97 A = LinSys[0].transpose() * LinSys[0];
98 b = LinSys[0].transpose() * LinSys[1];
99 WarningInFunction <<
100 "Using normal equation, results might be inaccurate, better to rely on completeOrthogonalDecomposition, bdcSvd or jacobiSvd"
101 << endl;
102 }
103 }
104
105 for (int i = 0; i < bc.size(); i++)
106 {
107 A.row(i) *= 0;
108 A(i, i) = 1;
109 b(i, 0) = bc(i);
110 }
111
112 if (solverType == "fullPivLu")
113 {
114 y = A.fullPivLu().solve(b);
115 }
116 else if (solverType == "partialPivLu")
117 {
118 y = A.partialPivLu().solve(b);
119 }
120 else if (solverType == "householderQr")
121 {
122 y = A.householderQr().solve(b);
123 }
124 else if (solverType == "colPivHouseholderQr")
125 {
126 y = A.colPivHouseholderQr().solve(b);
127 }
128 else if (solverType == "fullPivHouseholderQr")
129 {
130 y = A.fullPivHouseholderQr().solve(b);
131 }
132 else if (solverType == "completeOrthogonalDecomposition")
133 {
134 y = A.completeOrthogonalDecomposition().solve(b);
135 }
136 else if (solverType == "llt")
137 {
138 y = A.llt().solve(b);
139 }
140 else if (solverType == "ldlt")
141 {
142 y = A.ldlt().solve(b);
143 }
144 else if (solverType == "bdcSvd")
145 {
146 y = A.bdcSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
147 }
148 else if (solverType == "jacobiSvd")
149 {
150 if (A.rows() > A.cols())
151 {
152 y = A.jacobiSvd(Eigen::ComputeThinU | Eigen::ComputeThinV).solve(b);
153 }
154 else
155 {
156 y = A.jacobiSvd().solve(b);
157 }
158 }
159
160 return y;
161}
162
163Eigen::MatrixXd reducedProblem::solveLinearSys(List<Eigen::MatrixXd> LinSys,
164 Eigen::MatrixXd x, Eigen::VectorXd& residual, const std::string solverType)
165{
166 const Eigen::MatrixXd& bc = Eigen::MatrixXd::Zero(0, 0);
167 Eigen::MatrixXd y = reducedProblem::solveLinearSys(LinSys, x, residual, bc,
168 solverType);
169 return y;
170}
171
172// ****************** //
173// class onlineInterp //
174// ****************** //
175
176// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
177
178// Constructor
179
181{
183 ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt").cols();
184}
185
186Eigen::MatrixXd onlineInterp::getInterpCoeffRBF(std::vector<SPLINTER::RBFSpline>
187 rbfVec, Eigen::MatrixXd mu_interp)
188{
189 M_Assert(mu_interp.cols() == Nmu_samples,
190 "Matrix 'mu_interp' must have same number and order of columns (i.e. parameters) as the matrix 'mu_samples'.");
191 int Nmodes = rbfVec.size();
192 // Nsamples (snapshots) to evaluate the interpolator at
193 int Nsamples = mu_interp.rows();
194 Eigen::MatrixXd muInterpRefined(mu_interp.rows(), 0);
195
196 // Discard columns of mu_interp that have constant value for all elements, and use a refined matrix muInterpRefined
197 for (int i = 0; i < mu_interp.cols(); i++)
198 {
199 double firstVal = mu_interp(0, i);
200 bool isConst = (mu_interp.col(i).array() == firstVal).all();
201
202 // Consider only columns with non-constant values
203 if (isConst == 0)
204 {
205 muInterpRefined.conservativeResize(muInterpRefined.rows(),
206 muInterpRefined.cols() + 1);
207 muInterpRefined.col(muInterpRefined.cols() - 1) = mu_interp.col(i);
208 }
209 }
210
211 int Nmu = muInterpRefined.cols();
212 Eigen::MatrixXd coeff_interp(Nmodes, Nsamples);
213
214 for (int i = 0; i < Nmodes; i++)
215 {
216 for (int j = 0; j < Nsamples; j++)
217 {
218 SPLINTER::DenseVector x(Nmu);
219
220 for (int k = 0; k < Nmu; k++)
221 {
222 x(k) = muInterpRefined(j, k);
223 }
224
225 coeff_interp(i, j) = rbfVec[i].eval(x);
226 }
227 }
228
229 return coeff_interp;
230}
231
232
233Eigen::MatrixXd onlineInterp::getInterpCoeffSPL(std::vector<SPLINTER::BSpline>
234 splVec, Eigen::MatrixXd mu_interp)
235{
236 M_Assert(mu_interp.cols() == Nmu_samples,
237 "Matrix 'mu_interp' must have same number and order of columns (i.e. parametes) as the matrix 'mu_samples'.");
238 int Nmodes = splVec.size();
239 // Nsamples (snapshots) to evaluate the interpolator at
240 int Nsamples = mu_interp.rows();
241 Eigen::MatrixXd muInterpRefined(mu_interp.rows(), 0);
242
243 // Discard columns of mu_interp that have constant value for all elements, and use a refined matrix muInterpRefined
244 for (int i = 0; i < mu_interp.cols(); i++)
245 {
246 double firstVal = mu_interp(0, i);
247 bool isConst = (mu_interp.col(i).array() == firstVal).all();
248
249 // Consider only columns with non-constant values
250 if (isConst == 0)
251 {
252 muInterpRefined.conservativeResize(muInterpRefined.rows(),
253 muInterpRefined.cols() + 1);
254 muInterpRefined.col(muInterpRefined.cols() - 1) = mu_interp.col(i);
255 }
256 }
257
258 int Nmu = muInterpRefined.cols();
259 Eigen::MatrixXd coeff_interp(Nmodes, Nsamples);
260
261 for (int i = 0; i < Nmodes; i++)
262 {
263 for (int j = 0; j < Nsamples; j++)
264 {
265 SPLINTER::DenseVector x(Nmu);
266
267 for (int k = 0; k < Nmu; k++)
268 {
269 x(k) = muInterpRefined(j, k);
270 }
271
272 coeff_interp(i, j) = splVec[i].eval(x);
273 }
274 }
275
276 return coeff_interp;
277}
278
279// ************************************************************************* //
280
281
#define M_Assert(Expr, Msg)
scalar residual
Header file of the reducedProblem class.
int Nmu_samples
Number of parameters in the mu_samples matrix in the FOM.
Eigen::MatrixXd getInterpCoeffRBF(std::vector< SPLINTER::RBFSpline > rbfVec, Eigen::MatrixXd mu_interp)
Get interpolated coefficients evaluated at points from matrix "mu_interp" using the constructed param...
onlineInterp()
Construct Null.
Eigen::MatrixXd getInterpCoeffSPL(std::vector< SPLINTER::BSpline > splVec, Eigen::MatrixXd mu_interp)
Get interpolated coefficients evaluated at points from matrix "mu_interp" using the constructed param...
virtual void solveOnline()
Virtual Method to perform and online Solve.
static Eigen::MatrixXd solveLinearSys(List< Eigen::MatrixXd > LinSys, Eigen::MatrixXd x, Eigen::VectorXd &residual, const Eigen::MatrixXd &bc=Eigen::MatrixXd::Zero(0, 0), const std::string solverType="fullPivLu")
Linear system solver for the online problem.
reducedProblem()
Construct Null.
A general class for the implementation of a full order parametrized problem.
volScalarField & A
List< Eigen::MatrixXd > readMatrix(word folder, word mat_name)
Read a three dimensional matrix from a txt file in Eigen format.
label i
Definition pEqn.H:46