Loading...
Searching...
No Matches
ITHACAPOD.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 ITHACAPOD
26Description
27 Implementation of a POD decomposition of a general field
28SourceFiles
29 ITHACAPOD.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ITHACAPOD_H
38#define ITHACAPOD_H
39
40#include "fvCFD.H"
41#include "ITHACAutilities.H"
42#include "ITHACAstream.H"
43#include "ITHACAparameters.H"
44#include "Foam2Eigen.H"
45#include "EigenFunctions.H"
46#pragma GCC diagnostic push
47#pragma GCC diagnostic ignored "-Wold-style-cast"
48#pragma GCC diagnostic ignored "-Wnon-virtual-dtor"
49#include <Spectra/GenEigsSolver.h>
50#include <Spectra/SymEigsSolver.h>
51#include <Eigen/Eigen>
52#include <unsupported/Eigen/SparseExtra>
53#pragma GCC diagnostic pop
54
55/*---------------------------------------------------------------------------*\
56 Class reductionProblem Declaration
57\*---------------------------------------------------------------------------*/
58
60namespace ITHACAPOD
61{
62
63// Public Members
64
65//------------------------------------------------------------------------------
84template<class Type, template<class> class PatchField, class GeoMesh>
85void getModes(
86 PtrList<GeometricField<Type, PatchField, GeoMesh >> & snapshots,
87 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes,
88 word fieldName, bool podex, bool supex = 0, bool sup = 0,
89 label nmodes = 0, bool correctBC = true);
90
91//------------------------------------------------------------------------------
111template<class Type, template<class> class PatchField, class GeoMesh>
112void getModesSVD(
113 PtrList<GeometricField<Type, PatchField, GeoMesh >> &
114 snapshots, PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes,
115 word fieldName, bool podex, bool supex = 0, bool sup = 0,
116 label nmodes = 0, bool correctBC = true);
117
118
119//------------------------------------------------------------------------------
136template<class Type, template<class> class PatchField, class GeoMesh>
138 PtrList<GeometricField<Type, PatchField, GeoMesh >> & snapshots,
139 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes,
140 word fieldName, label Npar,
141 label NnestedOut);
142
143
144//------------------------------------------------------------------------------
164template<class Type, template<class> class PatchField, class GeoMesh>
166 PtrList<GeometricField<Type, PatchField, GeoMesh >> & snapshots,
167 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes,
168 word fieldName, bool podex, bool supex = 0, bool sup = 0,
169 label nmodes = 0, bool correctBC = true);
170
171//------------------------------------------------------------------------------
176void GrammSchmidt(Eigen::MatrixXd& Matrix);
177
178//------------------------------------------------------------------------------
191template<class Field_type>
192Eigen::MatrixXd corMatrix(Field_type& snapshots);
193
194
195//------------------------------------------------------------------------------
213template<class Type, template<class> class PatchField, class GeoMesh>
214void exportBases(
215 PtrList<GeometricField<Type, PatchField, GeoMesh >> & snapshots,
216 PtrList<GeometricField<Type, PatchField, GeoMesh >>& bases,
217 word fieldName,
218 bool sup = 0);
219
220//------------------------------------------------------------------------------
229void exportEigenvalues(scalarField Eigenvalues, fileName name,
230 bool sup = 0);
231
232//------------------------------------------------------------------------------
241void exportcumEigenvalues(scalarField cumEigenvalues, fileName name,
242 bool sup = 0);
243
244//------------------------------------------------------------------------------
255template<class Type, template<class> class PatchField, class GeoMesh>
256PtrList<GeometricField<Type, PatchField, GeoMesh >> DEIMmodes(
257 PtrList<GeometricField<Type, PatchField, GeoMesh >>& SnapShotsMatrix,
258 label nmodes, word FunctionName, word FieldName);
259
260//------------------------------------------------------------------------------
275template<typename type_matrix>
276std::tuple<List<Eigen::SparseMatrix<double >>, List<Eigen::VectorXd >>
277DEIMmodes(PtrList<type_matrix> & MatrixList, label nmodesA, label nmodesB,
278 word MatrixName = "Matrix");
279
280//------------------------------------------------------------------------------
296std::tuple<List<Eigen::SparseMatrix<double >>, List<Eigen::VectorXd >>
297DEIMmodes(List<Eigen::SparseMatrix<double >> & A, List<Eigen::VectorXd>& b,
298 label nmodesA, label nmodesB, word MatrixName = "Matrix");
299
300
301//------------------------------------------------------------------------------
321template<class Type, template<class> class PatchField, class GeoMesh>
322void getModes(PtrList<GeometricField<Type, PatchField, GeoMesh >> & snapshots,
323 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes,
324 PtrList<volScalarField>& Volumes, word fieldName, bool podex, bool supex = 0,
325 bool sup = 0,
326 label nmodes = 0);
327
328//------------------------------------------------------------------------------
351template<class Field_type, class Field_type_2>
352void getModes(PtrList<Field_type>& snapshots, PtrList<Field_type>& modes,
353 PtrList<Field_type_2>& fields2, word fieldName, bool podex, bool supex = 0,
354 bool sup = 0,
355 label nmodes = 0);
356
357//------------------------------------------------------------------------------
377template<class Type, template<class> class PatchField, class GeoMesh>
379 GeometricField<Type, PatchField, GeoMesh>& templateField,
380 word snapshotsPath,
381 PtrList<GeometricField<Type, PatchField, GeoMesh >> & modes,
382 word fieldName,
383 bool podex = false,
384 bool supex = false,
385 bool sup = false,
386 label nmodes = 0,
387 bool correctBC = true,
388 autoPtr<GeometricField<Type, PatchField, GeoMesh >> meanField = NULL);
389
390//------------------------------------------------------------------------------
392template<class Type, template<class> class PatchField, class GeoMesh>
394 GeometricField<Type, PatchField, GeoMesh>& templateField,
395 word snapshotsPath,
396 autoPtr<GeometricField<Type, PatchField, GeoMesh >> & meanField,
397 bool meanex = false);
398
399//------------------------------------------------------------------------------
407template<class Type, template<class> class PatchField, class GeoMesh >
409 const GeometricField<Type, PatchField, GeoMesh>& field1,
410 const GeometricField<Type, PatchField, GeoMesh>& field2);
411
412//------------------------------------------------------------------------------
420template<class Type, template<class> class PatchField, class GeoMesh>
422 const GeometricField<Type, PatchField, GeoMesh>& field1,
423 const GeometricField<Type, PatchField, GeoMesh>& field2);
424
425};
426
427// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
428
429#endif
430
431
432
Header file of the EigenFunctions class.
Header file of the Foam2Eigen class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
volScalarField & A
namespace for the computation of the POD, it exploits bot the method of snapshots and SVD.
scalar computeFrobeniusInnerProduct(const GeometricField< scalar, fvPatchField, volMesh > &field1, const GeometricField< scalar, fvPatchField, volMesh > &field2)
Definition ITHACAPOD.C:2047
void exportcumEigenvalues(scalarField cumEigenvalues, fileName name, bool sup)
Export the cumulative eigenvalues as a txt file.
Definition ITHACAPOD.C:1107
void getModesMemoryEfficient(GeometricField< Type, PatchField, GeoMesh > &templateField, word snapshotsPath, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC, autoPtr< GeometricField< Type, PatchField, GeoMesh > > meanField)
Gets the modes in a memory-efficient manner.
Definition ITHACAPOD.C:317
void exportEigenvalues(scalarField Eigenvalues, fileName name, bool sup)
Exports the eigenvalues as a txt file.
Definition ITHACAPOD.C:1078
scalar computeInnerProduct(const GeometricField< scalar, fvPatchField, volMesh > &field1, const GeometricField< scalar, fvPatchField, volMesh > &field2)
Definition ITHACAPOD.C:2031
void getMeanMemoryEfficient(GeometricField< Type, PatchField, GeoMesh > &templateField, word snapshotsPath, autoPtr< GeometricField< Type, PatchField, GeoMesh > > &meanField, bool meanex)
Gets the mean field in a memory-efficient manner.
Definition ITHACAPOD.C:633
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:90
std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > DEIMmodes(List< Eigen::SparseMatrix< double > > &A, List< Eigen::VectorXd > &b, label nmodesA, label nmodesB, word MatrixName)
Get the DEIM modes for a generic non a parametrized matrix coming from a differential operator functi...
Definition ITHACAPOD.C:1138
void getNestedSnapshotMatrix(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &ModesGlobal, word fieldName, label Npar, label NnestedOut)
Nested-POD approach.
Definition ITHACAPOD.C:41
void exportBases(PtrList< GeometricField< Type, PatchField, GeoMesh > > &s, PtrList< GeometricField< Type, PatchField, GeoMesh > > &bases, word fieldName, bool sup)
Export the Bases.
Definition ITHACAPOD.C:1040
void getModesSVD(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Gets the bases for a scalar field using SVD instead of the method of snapshots.
Definition ITHACAPOD.C:826
Eigen::MatrixXd corMatrix(PtrList< volScalarField > &snapshots)
Construct the Correlation Matrix for Scalar Field.
Definition ITHACAPOD.C:923
void GrammSchmidt(Eigen::MatrixXd &Matrix)
Performs GrammSchmidt orthonormalization on an Eigen Matrix.
Definition ITHACAPOD.C:1324
void getWeightedModes(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 weighted bases (using the nested-pod approach) or read them for a vector field.
Definition ITHACAPOD.C:687