37#ifndef ITHACAcoeffsMass_H
38#define ITHACAcoeffsMass_H
42#include "freestreamFvPatchField.H"
45#pragma GCC diagnostic push
46#pragma GCC diagnostic ignored "-Wold-style-cast"
48#pragma GCC diagnostic pop
51#include "polyMeshTools.H"
53#include "mixedFvPatchFields.H"
54#include "fvMeshSubset.H"
55using namespace std::placeholders;
78template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
80 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes,
81 Eigen::MatrixXd& coeff_matrix, label Nmodes);
95template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
97 PtrList<GeometricField<Type, PatchField, GeoMesh >>
99 fields, label nModes = 0,
100 bool consider_volumes =
true);
114template<
class Type,
template<
class>
class PatchField,
class GeoMesh >
116 PtrList<GeometricField<Type, PatchField, GeoMesh >>
119 PtrList<GeometricField<Type, PatchField, GeoMesh >>&
120 fields2, label nModes = 0,
121 bool consider_volumes =
true);
137template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
139 PtrList<GeometricField<Type, PatchField, GeoMesh >>
141 PtrList<GeometricField<Type, PatchField, GeoMesh >>& fields2,
142 Eigen::VectorXd weights,
144 bool consider_volumes =
true);
158template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
160 GeometricField<Type, PatchField, GeoMesh>& snapshot);
177template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
178Eigen::VectorXd
getCoeffs(GeometricField<Type, PatchField, GeoMesh>&
180 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes, label Nmodes = 0,
181 bool consider_volumes =
true);
198template<
class Type,
template<
class>
class PatchField,
class GeoMesh >
199Eigen::MatrixXd
getCoeffs(PtrList<GeometricField<Type, PatchField, GeoMesh >>&
201 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes, label Nmodes = 0,
202 bool consider_volumes =
true);
243 acquiredSnapshotsTimes, Eigen::MatrixXd parameters);
Header file of the Foam2Eigen class.
Header file of the ITHACAfieldsOperations file.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Simple header and source file of the Color::Modifier class to change color to the output stream.
Namespace to implement some useful assign operation of OF fields.
Eigen::VectorXd getCoeffs(GeometricField< Type, PatchField, GeoMesh > &snapshot, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Projects a snapshot on a basis function and gets the coefficients of the projection.
PtrList< GeometricField< Type, PatchField, GeoMesh > > reconstructFromCoeff(PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, Eigen::MatrixXd &coeff_matrix, label Nmodes)
Exact reconstruction using a certain number of modes for a list of fields and the projection coeffici...
T project_to_POD_basis(T &field, PtrList< T > &modes)
Compute the projection of a field in a POD basis.
Eigen::VectorXd getMassMatrixFV(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Gets a vector containing the volumes of each cell of the mesh.
Eigen::MatrixXd parTimeCombMat(List< Eigen::VectorXd > acquiredSnapshotsTimes, Eigen::MatrixXd parameters)
A method to compute the time-parameter combined matrix whose any single element corresponds to a uniq...
Eigen::MatrixXd getMassMatrix(PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Gets the mass matrix using the eigen routine.