Loading...
Searching...
No Matches
09DEIM_ROM.C
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-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Description
25 Example of model reduction problem using the DEIM for a Heat Transfer Problem
26SourceFiles
27 09DEIM_ROM.C
28\*---------------------------------------------------------------------------*/
29
30#include "fvCFD.H"
31#include "IOmanip.H"
32#include "fvOptions.H"
33#include "simpleControl.H"
34#include "ITHACAutilities.H"
35#include "Foam2Eigen.H"
36#include <chrono>
37#include "fvMeshSubset.H"
38#include "ITHACAstream.H"
39#include "ITHACAPOD.H"
40#include "EigenFunctions.H"
41#include "DEIM.H"
42#include <chrono>
43#include <Eigen/SVD>
44#include <Eigen/SparseLU>
45#include "laplacianProblem.H"
46
47class DEIM_function : public DEIM<fvScalarMatrix>
48{
49 using DEIM::DEIM;
50 public:
51 static fvScalarMatrix evaluate_expression(volScalarField& T, Eigen::MatrixXd mu)
52 {
53 volScalarField yPos = T.mesh().C().component(vector::Y).ref();
54 volScalarField xPos = T.mesh().C().component(vector::X).ref();
55 volScalarField nu(T);
56
57 for (auto i = 0; i < nu.size(); i++)
58 {
59 nu[i] = std::exp(- 2 * std::pow(xPos[i] - mu(0) - 1,
60 2) - 2 * std::pow(yPos[i] - mu(1) - 0.5, 2)) + 1;
61 }
62
63 nu.correctBoundaryConditions();
64 fvScalarMatrix TiEqn22
65 (
66 fvm::laplacian(nu, T, "Gauss linear")
67 );
68 return TiEqn22;
69 }
70
71 Eigen::MatrixXd onlineCoeffsA(Eigen::MatrixXd mu)
72 {
73 Eigen::MatrixXd theta(magicPointsAcol().size(), 1);
74 fvScalarMatrix Aof = evaluate_expression(fieldA(), mu);
75 Eigen::SparseMatrix<double> Mr;
76 Eigen::VectorXd br;
77 Foam2Eigen::fvMatrix2Eigen(Aof, Mr, br);
78
79 for (int i = 0; i < magicPointsAcol().size(); i++)
80 {
81 int ind_row = localMagicPointsArow[i] + xyz_Arow()[i] *
82 fieldA().size();
83 int ind_col = localMagicPointsAcol[i] + xyz_Acol()[i] *
84 fieldA().size();
85 theta(i) = Mr.coeffRef(ind_row, ind_col);
86 }
87
88 return theta;
89 }
90
91 Eigen::MatrixXd onlineCoeffsB(Eigen::MatrixXd mu)
92 {
93 Eigen::MatrixXd theta(magicPointsB().size(), 1);
94 fvScalarMatrix Aof = evaluate_expression(fieldB(), mu);
95 Eigen::SparseMatrix<double> Mr;
96 Eigen::VectorXd br;
97 Foam2Eigen::fvMatrix2Eigen(Aof, Mr, br);
98
99 for (int i = 0; i < magicPointsB().size(); i++)
100 {
101 int ind_row = localMagicPointsB[i] + xyz_B()[i] * fieldB().size();
102 theta(i) = br(ind_row);
103 }
104
105 return theta;
106 }
107
108 PtrList<volScalarField> fieldsA;
109 autoPtr<volScalarField> fieldA;
110 autoPtr<volScalarField> fieldB;
111 PtrList<volScalarField> fieldsB;
112};
113
114class DEIMLaplacian: public laplacianProblem
115{
116 public:
117 explicit DEIMLaplacian(int argc, char* argv[])
118 :
119 laplacianProblem(argc, argv),
120 nu(_nu()),
121 S(_S()),
122 T(_T())
123 {
124 fvMesh& mesh = _mesh();
125 ITHACAdict = new IOdictionary
126 (
127 IOobject
128 (
129 "ITHACAdict",
130 "./system",
131 mesh,
132 IOobject::MUST_READ,
133 IOobject::NO_WRITE
134 )
135 );
136 NTmodes = readInt(ITHACAdict->lookup("N_modes_T"));
137 NmodesDEIMA = readInt(ITHACAdict->lookup("N_modes_DEIM_A"));
138 NmodesDEIMB = readInt(ITHACAdict->lookup("N_modes_DEIM_B"));
139 }
140
141 volScalarField& nu;
142 volScalarField& S;
143 volScalarField& T;
144
145
146 DEIM_function* DEIMmatrice;
147 PtrList<fvScalarMatrix> Mlist;
148 Eigen::MatrixXd ModesTEig;
149 std::vector<Eigen::MatrixXd> ReducedMatricesA;
150 std::vector<Eigen::MatrixXd> ReducedVectorsB;
151
152 int NTmodes;
153 int NmodesDEIMA;
154 int NmodesDEIMB;
155
156 double time_full;
157 double time_rom;
158
159 void OfflineSolve(Eigen::MatrixXd par, word Folder)
160 {
161 if (offline)
162 {
163 ITHACAstream::read_fields(Tfield, T, "./ITHACAoutput/Offline/");
164 }
165 else
166 {
167 for (int i = 0; i < par.rows(); i++)
168 {
169 fvScalarMatrix Teqn = DEIMmatrice->evaluate_expression(T, par.row(i));
170 Teqn.solve();
171 Mlist.append((Teqn).clone());
172 ITHACAstream::exportSolution(T, name(i + 1), "./ITHACAoutput/" + Folder);
173 Tfield.append((T).clone());
174 }
175 }
176 };
177
178 void OnlineSolveFull(Eigen::MatrixXd par, word Folder)
179 {
180 auto t1 = std::chrono::high_resolution_clock::now();
181 auto t2 = std::chrono::high_resolution_clock::now();
182 auto time_span = std::chrono::duration_cast<std::chrono::duration<double >>
183 (t2 - t1);
184 time_full = 0;
185
186 for (int i = 0; i < par.rows(); i++)
187 {
188 fvScalarMatrix Teqn = DEIMmatrice->evaluate_expression(T, par.row(i));
189 t1 = std::chrono::high_resolution_clock::now();
190 Teqn.solve();
191 t2 = std::chrono::high_resolution_clock::now();
192 time_span = std::chrono::duration_cast<std::chrono::duration<double >>
193 (t2 - t1);
194 time_full += time_span.count();
195 ITHACAstream::exportSolution(T, name(i + 1), "./ITHACAoutput/" + Folder);
196 Tfield.append((T).clone());
197 }
198 };
199
200 void PODDEIM()
201 {
202 PODDEIM(NTmodes, NmodesDEIMA, NmodesDEIMB);
203 }
204
205 void PODDEIM(int NmodesT, int NmodesDEIMA, int NmodesDEIMB)
206 {
207 DEIMmatrice = new DEIM_function(Mlist, NmodesDEIMA, NmodesDEIMB, "T_matrix");
208 fvMesh& mesh = const_cast<fvMesh&>(T.mesh());
209 // Differential Operator
210 DEIMmatrice->fieldA = autoPtr<volScalarField>(new volScalarField(
211 DEIMmatrice->generateSubmeshMatrix(2, mesh, T)));
212 DEIMmatrice->fieldB = autoPtr<volScalarField>(new volScalarField(
213 DEIMmatrice->generateSubmeshVector(2, mesh, T)));
214 // Source Terms
216 ModesTEig.conservativeResize(ModesTEig.rows(), NmodesT);
217 ReducedMatricesA.resize(NmodesDEIMA);
218 ReducedVectorsB.resize(NmodesDEIMB);
219
220 for (int i = 0; i < NmodesDEIMA; i++)
221 {
222 ReducedMatricesA[i] = ModesTEig.transpose() * DEIMmatrice->MatrixOnlineA[i] *
223 ModesTEig;
224 }
225
226 for (int i = 0; i < NmodesDEIMB; i++)
227 {
228 ReducedVectorsB[i] = ModesTEig.transpose() * DEIMmatrice->MatrixOnlineB;
229 }
230 };
231
232 void OnlineSolve(Eigen::MatrixXd par_new, word Folder)
233 {
234 auto t1 = std::chrono::high_resolution_clock::now();
235 auto t2 = std::chrono::high_resolution_clock::now();
236 auto time_span = std::chrono::duration_cast<std::chrono::duration<double >>
237 (t2 - t1);
238 time_rom = 0;
239
240 for (int i = 0; i < par_new.rows(); i++)
241 {
242 // solve
243 t1 = std::chrono::high_resolution_clock::now();
244 Eigen::MatrixXd thetaonA = DEIMmatrice->onlineCoeffsA(par_new.row(i));
245 Eigen::MatrixXd thetaonB = DEIMmatrice->onlineCoeffsB(par_new.row(i));
246 Eigen::MatrixXd A = EigenFunctions::MVproduct(ReducedMatricesA, thetaonA);
247 Eigen::VectorXd B = EigenFunctions::MVproduct(ReducedVectorsB, thetaonB);
248 Eigen::VectorXd x = A.fullPivLu().solve(B);
249 Eigen::VectorXd full = ModesTEig * x;
250 t2 = std::chrono::high_resolution_clock::now();
251 time_span = std::chrono::duration_cast<std::chrono::duration<double >>
252 (t2 - t1);
253 time_rom += time_span.count();
254 // Export
255 volScalarField Tred("Tred", T);
256 Tred = Foam2Eigen::Eigen2field(Tred, full);
257 ITHACAstream::exportSolution(Tred, name(i + 1), "./ITHACAoutput/" + Folder);
258 Tonline.append((Tred).clone());
259 }
260 }
261};
262
263int main(int argc, char* argv[])
264{
265 // Construct the case
266 DEIMLaplacian example(argc, argv);
267 // Read some parameters from file
268 ITHACAparameters* para = ITHACAparameters::getInstance(example._mesh(),
269 example._runTime());
270 // Create the offline parameters for the solve
271 example.mu = ITHACAutilities::rand(100, 2, -0.5, 0.5);
272 // Solve the offline problem to compute the snapshots for the projections
273 example.OfflineSolve(example.mu, "Offline");
274 // Compute the POD modes
275 ITHACAPOD::getModes(example.Tfield, example.Tmodes, example._T().name(),
276 example.podex, 0, 0, 20);
277 // Compute the offline part of the DEIM procedure
278 example.PODDEIM();
279 // Construct a new set of parameters
280 Eigen::MatrixXd par_new1 = ITHACAutilities::rand(100, 2, -0.5, 0.5);
281 // Solve the online problem with the new parameters
282 example.OnlineSolve(par_new1, "Online_red");
283 // Solve a new full problem with the new parameters (necessary to compute speed up and error)
284 DEIMLaplacian example_new(argc, argv);
285 example_new.OnlineSolveFull(par_new1, "Online_full");
286 // Output some infos
287 Info << endl << "The FOM Solve took: " << example_new.time_full <<
288 " seconds." << endl;
289 Info << endl << "The ROM Solve took: " << example.time_rom <<
290 " seconds." << endl;
291 Info << endl << "The Speed-up is: " << example_new.time_full /
292 example.time_rom << endl << endl;
293 Eigen::MatrixXd error = ITHACAutilities::errorL2Abs(example_new.Tfield,
294 example.Tonline);
295 Info << "The mean L2 error is: " << error.mean() << endl;
296 return 0;
297}
Header file of the EigenFunctions class.
Header file of the Foam2Eigen class.
Header file of the ITHACAPOD class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
Definition DEIM.H:41
Eigen::VectorXd theta
Definition DEIM.H:89
autoPtr< IOList< label > > xyz_Arow
Definition DEIM.H:152
List< label > localMagicPointsAcol
Definition DEIM.H:161
autoPtr< IOList< label > > magicPointsAcol
Definition DEIM.H:125
List< label > localMagicPointsB
Definition DEIM.H:162
autoPtr< IOList< label > > xyz_B
Definition DEIM.H:154
DEIM(PtrList< T > &SnapShotsMatrix, label MaxModes, word FunctionName, word FieldName)
Construct DEIM for non-linear function.
Definition DEIM.C:34
List< label > localMagicPointsArow
Definition DEIM.H:160
autoPtr< IOList< label > > xyz_Acol
Definition DEIM.H:153
autoPtr< IOList< label > > magicPointsB
Definition DEIM.H:127
static GeometricField< scalar, PatchField, GeoMesh > Eigen2field(GeometricField< scalar, PatchField, GeoMesh > &field, Eigen::VectorXd &eigen_vector, bool correctBC=true)
Convert a vector in Eigen format into an OpenFOAM scalar GeometricField.
Definition Foam2Eigen.C:533
static Eigen::MatrixXd PtrList2Eigen(PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert a PtrList of snapshots to Eigen matrix (only internal field).
Definition Foam2Eigen.C:665
static void fvMatrix2Eigen(fvMatrix< type_foam_matrix > foam_matrix, type_A &A, type_B &b)
Convert a FvMatrix OpenFOAM matrix (Linear System) into a Eigen Matrix A and a source vector b.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
PtrList< volScalarField > Tonline
List of snapshots for the solution.
autoPtr< volScalarField > _T
Temperature field.
autoPtr< volScalarField > _nu
Diffusivity.
autoPtr< fvMesh > _mesh
Mesh.
volScalarModes Tmodes
List of POD modes.
autoPtr< volScalarField > _S
Source Term.
PtrList< volScalarField > Tfield
List of snapshots for the solution.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
IOdictionary * ITHACAdict
dictionary to store input output infos
Header file of the laplacianProblem class.
Eigen::SparseMatrix< T > MVproduct(List< Eigen::SparseMatrix< T > > &A, Eigen::DenseBase< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > > &C)
Sparse Matrix-Vector product between a list of sparse matrices and a vector of coefficients.
void getModes(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Computes the bases or reads them for a field.
Definition ITHACAPOD.C:93
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
Eigen::MatrixXd rand(label rows, label cols, double min, double max)
Generates random matrix with random values in an interval.