Loading...
Searching...
No Matches
hyperReduction.H
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-------------------------------------------------------------------------------
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/>.
24Class
25 HyperReduction
26Description
27 Implementation of the hyper-reduction methods.
28SourceFiles
29 hyperReduction.C
30\*---------------------------------------------------------------------------*/
31#ifndef hyperReduction_H
32#define hyperReduction_H
33
34#include "fvCFD.H"
35#include "ITHACAPOD.H"
36#include "Foam2Eigen.H"
37#include "EigenFunctions.H"
38#include "ITHACAutilities.H"
39#include "fvMeshSubset.H"
40#include <set>
41#include "redsvd"
42
43
44
45// TODO uniform style for separating offline and onlie stages
46template <typename... SnapshotsLists>
48{
49 public:
51
52 using SnapshotsListTuple = std::tuple<std::decay_t<SnapshotsLists>...>;
54 std::tuple<typename std::decay_t<SnapshotsLists>::value_type ...>;
55 static constexpr auto n_fields = sizeof...(SnapshotsLists);
56
57 // TODO get vectorial dim from class template parameters
58 unsigned int vectorial_dim;
59 unsigned int sumFieldsDim;
60
61 template <std::size_t N>
62 using NthFieldListType = typename
63 std::tuple_element<N, SnapshotsListTuple>::type;
64
65 template <std::size_t N>
67
68 //----------------------------------------------------------------------
79 label n_nodes,
80 Eigen::VectorXi initialSeeds,
81 word problemName,
82 SnapshotsLists &&...snapshotsLists);
83
84 //----------------------------------------------------------------------
95 label n_nodes,
96 unsigned int vectorial_dim,
97 label n_cells,
98 Eigen::VectorXi initialSeeds,
99 word problemName);
100
101 // ITHACA parameters dict
103
104 //TODO
106
107 //TODO
109
110 //TODO
112
114 label n_modes;
115
117 label n_nodes;
118
121
122 // Initial nodes
123 Eigen::VectorXi initialSeeds;
124
125 // nodes
126 Eigen::VectorXi nodes;
127
130
133
135 List<word> fieldNames;
136
138 List<unsigned int> fieldDims;
139
141 label n_cells;
142
145
147 autoPtr<IOList<label >> nodePoints;
148
150 autoPtr<IOList<labelList >> totalNodePoints;
151
153 autoPtr<IOList<label >> uniqueNodePoints;
154
157
160
162 List<label> localNodePoints;
163
164 // Outputs
166 Eigen::MatrixXd basisMatrix;
167
168 // Normalizing weights of SVD decomposition
169 Eigen::VectorXd normalizingWeights;
170
172 Eigen::MatrixXd renormalizedBasisMatrix;
173
174 // TODO
175 Eigen::MatrixXd pinvPU;
176
177 // TODO
178 Eigen::MatrixXd wPU;
179
180 // TODO
181 Eigen::VectorXd eigenValueseig;
182
184 Eigen::SparseMatrix<double> P;
185
186 // Mask field to submesh
187 Eigen::SparseMatrix<double> field2submesh;
188
189 // Mask submesh to nodes
190 Eigen::SparseMatrix<double> submesh2nodes;
191
192 // Indices submesh to nodes
193 Eigen::VectorXi submesh2nodesMask;
194
196 Eigen::VectorXd quadratureWeights;
197
199 autoPtr<fvMeshSubset> submesh;
200
202 autoPtr<volVectorField> submesh_field;
203
205 Eigen::MatrixXd MatrixOnline;
206
207 //----------------------------------------------------------------------
210 void offlineGappyDEIM(Eigen::MatrixXd& snapshotsModes,
211 Eigen::VectorXd& normalizingWeights)
212 {
213 methodName = "GappyDEIM";
214 word name = "ITHACAoutput/" + problemName + "/" + methodName;
215 offlineGappyDEIM(snapshotsModes, normalizingWeights, name);
216 }
217
218 void offlineGappyDEIM(Eigen::MatrixXd& snapshotsModes,
219 Eigen::VectorXd& normalizingWeights, word folderMethodName);
220
221 //----------------------------------------------------------------------
224 void offlineECP(Eigen::MatrixXd& snapshotsModes,
225 Eigen::VectorXd& normalizingWeights)
226 {
227 methodName = "ECP";
228 word name = "ITHACAoutput/" + problemName + "/" + methodName;
229 offlineECP(snapshotsModes, normalizingWeights, name);
230 }
231
232 void offlineECP(Eigen::MatrixXd& snapshotsModes,
233 Eigen::VectorXd& normalizingWeights, word folderMethodName);
234
235 //----------------------------------------------------------------------
238 void initSeeds(Eigen::VectorXd mp_not_mask, std::set<label> nodePointsSet);
239
240 //----------------------------------------------------------------------
243 void computeLS(Eigen::MatrixXd& J, Eigen::MatrixXd& JWhole, Eigen::VectorXd& b,
244 Eigen::VectorXd& q);
245
246 //----------------------------------------------------------------------
253 Eigen::MatrixXd& modesSVD, Eigen::VectorXd& fieldWeights,
254 bool saveModesFlag = false);
255
256 //----------------------------------------------------------------------
263 Eigen::MatrixXd& modesSVD, Eigen::VectorXd& fieldWeights,
264 Eigen::MatrixXd& modesSVDBoundary, Eigen::VectorXd& fieldWeightsBoundary,
265 bool saveModesFlag);
266
267 //----------------------------------------------------------------------
273 void updateNodes(Eigen::SparseMatrix<double>& P, label& ind_max,
274 Eigen::VectorXd& mp_not_mask);
275
276 //----------------------------------------------------------------------
282 template <typename SnapshotsList>
283 void stackSnapshots(SnapshotsList sList, Eigen::MatrixXd& snapshotsMatrix,
284 Eigen::VectorXd& fieldWeights);
285
286 //----------------------------------------------------------------------
292 template <typename SnapshotsList>
293 void stackSnapshotsBoundary(SnapshotsList sList,
294 List<Eigen::MatrixXd>& snapshotsMatrixBoundary,
295 List<Eigen::VectorXd>& fieldWeightsBoundary);
296
297 //----------------------------------------------------------------------
303 template <typename SnapshotsList>
304 void saveModes(SnapshotsList sList, Eigen::MatrixXd& snapshotsMatrix,
305 unsigned int& rowIndex, unsigned int& modeIndex, word folder);
306
307 //----------------------------------------------------------------------
313 template <typename SnapshotsList>
314 void saveModes(SnapshotsList sList, Eigen::MatrixXd& snapshotsMatrix,
315 Eigen::MatrixXd& snapshotsMatrixBoundary, unsigned int& rowIndex,
316 unsigned int& rowIndexBoundary, unsigned int& modeIndex, word folder);
317
318 //----------------------------------------------------------------------
324 template <typename SnapshotsList>
325 void stackNames(SnapshotsList sList)
326 {
327 fieldNames.append(sList[0].name());
328 }
329
330 //----------------------------------------------------------------------
336 template <typename SnapshotsList>
337 inline void stackDimensions(SnapshotsList sList)
338 {
340 }
341
342 //----------------------------------------------------------------------
348 template <typename SnapshotsList>
349 inline void sumDimensions(double sum, SnapshotsList sList)
350 {
352 }
353
354 //----------------------------------------------------------------------
360 void evaluatePinv(Eigen::SparseMatrix<double>& Projector,
361 Eigen::MatrixXd& Modes, Eigen::VectorXd& fieldWeights);
362
363 //----------------------------------------------------------------------
369 void evaluateWPU(Eigen::SparseMatrix<double>& Projector,
370 Eigen::MatrixXd& Modes, Eigen::VectorXd& fieldWeights,
371 Eigen::VectorXd& quadratureWeights);
372
373 //----------------------------------------------------------------------
379 void generateSubmesh(label layers, const fvMesh& mesh);
380
381 //----------------------------------------------------------------------
389 List<label> global2local(List<label>& points, fvMeshSubset& submesh);
390
396 template <typename FieldType>
397 inline autoPtr<FieldType> interpolateField(const FieldType& field)
398 {
399 return autoPtr<FieldType>(new FieldType(submesh->interpolate(field)));
400 }
401
407 // TODO expand with sphericalTensor, symmTensor, tensor
408 template <typename Field>
409 inline constexpr unsigned int get_field_dim()
410 {
411 return std::is_same<Field, volScalarField>::value ? 1 : 3;
412 }
413
419 template <typename LastList>
420 inline constexpr unsigned int compute_vectorial_dim(LastList x)
421 {
423 }
424
430 template <typename List, typename... RemainingLists>
431 inline constexpr unsigned int compute_vectorial_dim(List &&head,
432 RemainingLists &&...tail)
433 {
436 }
437
443 static inline double s_optimality(Eigen::MatrixXd& A)
444 {
445 return std::pow(
446 std::sqrt((A.transpose() * A).determinant()) / \
447 A.colwise().lpNorm<2>().prod(),
448 1.0 / A.cols());
449 }
450
456 template <typename... FieldsArgs>
457 void eigen2fields(Eigen::VectorXd& eFields, FieldsArgs&&... oFields)
458 {
459 unsigned int cumulativeSize{0}, ith_field{0};
460 ([ & ](auto & oFieldsArg)
461 {
462 double ith_fieldSize = fieldDims[ith_field] * oFieldsArg.size();
463 Eigen::VectorXd vec = eFields.segment(cumulativeSize, ith_fieldSize);
464 oFieldsArg = Foam2Eigen::Eigen2field(oFieldsArg, vec);
465 cumulativeSize += ith_fieldSize;
466 ith_field++;
467 }
468
469 (oFields), ...);
470 }
471
477 void initReshapeMat(Eigen::SparseMatrix<double>& reshapeMat);
478
484 void createMasks(bool offlineStage = true);
485
491 void getSnapMatrix(Eigen::MatrixXd& snapMatrix, Eigen::VectorXd& fieldWeights);
492
498 void getSnapMatrix(Eigen::MatrixXd& snapMatrix, Eigen::VectorXd& fieldWeights,
499 List<Eigen::MatrixXd>& snapMatrixBoundary,
500 List<Eigen::VectorXd>& fieldWeightsBoundary);
501};
502
503#endif
Header file of the EigenFunctions class.
Header file of the Foam2Eigen class.
Header file of the ITHACAPOD class.
Foam::fvMesh & mesh
Definition createMesh.H:47
Header file of the ITHACAutilities namespace.
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:514
Eigen::VectorXi submesh2nodesMask
static constexpr auto n_fields
typename NthFieldListType< N >::value_type NthFieldType
void stackSnapshotsBoundary(SnapshotsList sList, List< Eigen::MatrixXd > &snapshotsMatrixBoundary, List< Eigen::VectorXd > &fieldWeightsBoundary)
TODO.
SnapshotsListTuple snapshotsListTuple
The snapshots matrix containing the nonlinear function or operator.
label n_modes
The maximum number of modes to be considered.
void evaluateWPU(Eigen::SparseMatrix< double > &Projector, Eigen::MatrixXd &Modes, Eigen::VectorXd &fieldWeights, Eigen::VectorXd &quadratureWeights)
Compute the pseudo-inverse of the matrix M restricted with the projector P.
label n_nodes
The maximum number of modes to be considered.
word folderMethod
Folder for the selected HR method.
void stackNames(SnapshotsList sList)
TODO.
List< word > fieldNames
Names of the fields.
void generateSubmesh(label layers, const fvMesh &mesh)
Compute the submesh common to all fields in SnapshotsLists.
Eigen::VectorXi nodes
autoPtr< FieldType > interpolateField(const FieldType &field)
TODO.
autoPtr< fvMeshSubset > submesh
Submesh of the HyperReduction method.
Eigen::MatrixXd renormalizedBasisMatrix
Renormalized basis of HR.
List< label > n_boundary_cells_list
void offlineGappyDEIM(Eigen::MatrixXd &snapshotsModes, Eigen::VectorXd &normalizingWeights)
Methods implemented: 'GappyDEIM' from "DEIM, Chaturantabut, Saifon, and Danny C. Sorensen....
void updateNodes(Eigen::SparseMatrix< double > &P, label &ind_max, Eigen::VectorXd &mp_not_mask)
TODO.
unsigned int sumFieldsDim
void evaluatePinv(Eigen::SparseMatrix< double > &Projector, Eigen::MatrixXd &Modes, Eigen::VectorXd &fieldWeights)
Compute the pseudo-inverse of the matrix M restricted with the projector P.
Eigen::VectorXd quadratureWeights
Quadrature weights. Ordered in the same order of matrix P.
autoPtr< IOList< labelList > > totalNodePoints
List of label lists of the nodes and corresponding surrounding nodes.
void sumDimensions(double sum, SnapshotsList sList)
TODO.
Eigen::MatrixXd MatrixOnline
Online Matrix.
List< label > global2local(List< label > &points, fvMeshSubset &submesh)
Get local indices in the submesh from indices in the global ones.
Eigen::SparseMatrix< double > submesh2nodes
void stackDimensions(SnapshotsList sList)
TODO.
constexpr unsigned int compute_vectorial_dim(LastList x)
TODO.
void getSnapMatrix(Eigen::MatrixXd &snapMatrix, Eigen::VectorXd &fieldWeights)
TODO.
Eigen::VectorXi initialSeeds
Eigen::MatrixXd pinvPU
List< label > localNodePoints
Indices of the local node points in the subMesh.
constexpr unsigned int get_field_dim()
TODO.
Eigen::VectorXd normalizingWeights
void offlineECP(Eigen::MatrixXd &snapshotsModes, Eigen::VectorXd &normalizingWeights)
Methods implemented: 'ECP' from "ECP, Hernandez, Joaquin Alberto, Manuel Alejandro Caicedo,...
word problemName
The name of the non-linear function e.g. HR_method/residual.
constexpr unsigned int compute_vectorial_dim(List &&head, RemainingLists &&...tail)
TODO.
void eigen2fields(Eigen::VectorXd &eFields, FieldsArgs &&... oFields)
TODO.
HyperReduction(label n_modes, label n_nodes, Eigen::VectorXi initialSeeds, word problemName, SnapshotsLists &&...snapshotsLists)
Construct HyperReduction class, interpolation-based.
ITHACAparameters * para
std::tuple< std::decay_t< SnapshotsLists >... > SnapshotsListTuple
Eigen::VectorXd eigenValueseig
Eigen::SparseMatrix< double > P
The P matrix of the HyperReduction method. The nodes are ordered in the order of insertion during the...
typename std::tuple_element< N, SnapshotsListTuple >::type NthFieldListType
unsigned int vectorial_dim
label n_cells
Int Number of Cells;.
void initReshapeMat(Eigen::SparseMatrix< double > &reshapeMat)
TODO.
void initSeeds(Eigen::VectorXd mp_not_mask, std::set< label > nodePointsSet)
TODO.
autoPtr< IOList< label > > nodePoints
Nodes in the case of the a nonlinear function.
static double s_optimality(Eigen::MatrixXd &A)
TODO.
void getModesSVD(SnapshotsListTuple &SnapshotsListTuple, Eigen::MatrixXd &modesSVD, Eigen::VectorXd &fieldWeights, bool saveModesFlag=false)
TODO.
autoPtr< volVectorField > submesh_field
Submeshes.
label n_snapshots
The length of the snapshots lists.
autoPtr< IOList< label > > uniqueNodePoints
List of the unique indices of the nodes that define the submesh.
Eigen::SparseMatrix< double > field2submesh
List< unsigned int > fieldDims
Dimensions of the fields.
void saveModes(SnapshotsList sList, Eigen::MatrixXd &snapshotsMatrix, unsigned int &rowIndex, unsigned int &modeIndex, word folder)
TODO.
word folderProblem
Folder for the HR problem.
std::tuple< typename std::decay_t< SnapshotsLists >::value_type ... > FieldsTuple
void stackSnapshots(SnapshotsList sList, Eigen::MatrixXd &snapshotsMatrix, Eigen::VectorXd &fieldWeights)
TODO.
Eigen::MatrixXd basisMatrix
Orthonormal basis of HR.
label n_cellsSubfields
Int Number of Cells in submeshes;.
void createMasks(bool offlineStage=true)
TODO.
void computeLS(Eigen::MatrixXd &J, Eigen::MatrixXd &JWhole, Eigen::VectorXd &b, Eigen::VectorXd &q)
TODO.
Eigen::MatrixXd wPU
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
Implementation of a container class derived from PtrList.
Definition Modes.H:69
volScalarField & A