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>
52#include <unsupported/Eigen/SparseExtra>
53#pragma GCC diagnostic pop
84template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
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);
111template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
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);
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,
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);
191template<
class Field_type>
213template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
215 PtrList<GeometricField<Type, PatchField, GeoMesh >> & snapshots,
216 PtrList<GeometricField<Type, PatchField, GeoMesh >>& bases,
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);
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");
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");
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,
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,
377template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
379 GeometricField<Type, PatchField, GeoMesh>& templateField,
381 PtrList<GeometricField<Type, PatchField, GeoMesh >> & modes,
387 bool correctBC =
true,
388 autoPtr<GeometricField<Type, PatchField, GeoMesh >> meanField = NULL);
392template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
394 GeometricField<Type, PatchField, GeoMesh>& templateField,
396 autoPtr<GeometricField<Type, PatchField, GeoMesh >> & meanField,
397 bool meanex =
false);
407template<
class Type,
template<
class>
class PatchField,
class GeoMesh >
409 const GeometricField<Type, PatchField, GeoMesh>& field1,
410 const GeometricField<Type, PatchField, GeoMesh>& field2);
420template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
422 const GeometricField<Type, PatchField, GeoMesh>& field1,
423 const GeometricField<Type, PatchField, GeoMesh>& field2);
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.
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)
void exportcumEigenvalues(scalarField cumEigenvalues, fileName name, bool sup)
Export the cumulative eigenvalues as a txt file.
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.
void exportEigenvalues(scalarField Eigenvalues, fileName name, bool sup)
Exports the eigenvalues as a txt file.
scalar computeInnerProduct(const GeometricField< scalar, fvPatchField, volMesh > &field1, const GeometricField< scalar, fvPatchField, volMesh > &field2)
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.
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.
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...
void getNestedSnapshotMatrix(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &ModesGlobal, word fieldName, label Npar, label NnestedOut)
Nested-POD approach.
void exportBases(PtrList< GeometricField< Type, PatchField, GeoMesh > > &s, PtrList< GeometricField< Type, PatchField, GeoMesh > > &bases, word fieldName, bool sup)
Export the Bases.
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.
Eigen::MatrixXd corMatrix(PtrList< volScalarField > &snapshots)
Construct the Correlation Matrix for Scalar Field.
void GrammSchmidt(Eigen::MatrixXd &Matrix)
Performs GrammSchmidt orthonormalization on an Eigen Matrix.
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.