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;
77template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
79 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes,
80 Eigen::MatrixXd& coeff_matrix, label Nmodes);
94template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
96 PtrList<GeometricField<Type, PatchField, GeoMesh >>
98 fields, label nModes = 0,
99 bool consider_volumes =
true);
113template<
class Type,
template<
class>
class PatchField,
class GeoMesh >
115 PtrList<GeometricField<Type, PatchField, GeoMesh >>
118 PtrList<GeometricField<Type, PatchField, GeoMesh >>&
119 fields2, label nModes = 0,
120 bool consider_volumes =
true);
136template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
138 PtrList<GeometricField<Type, PatchField, GeoMesh >>
140 PtrList<GeometricField<Type, PatchField, GeoMesh >>& fields2,
141 Eigen::VectorXd weights,
143 bool consider_volumes =
true);
157template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
159 GeometricField<Type, PatchField, GeoMesh>& snapshot);
176template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
177Eigen::VectorXd
getCoeffs(GeometricField<Type, PatchField, GeoMesh>&
179 PtrList<GeometricField<Type, PatchField, GeoMesh >> & modes, label Nmodes = 0,
180 bool consider_volumes =
true);
197template<
class Type,
template<
class>
class PatchField,
class GeoMesh >
198Eigen::MatrixXd
getCoeffs(PtrList<GeometricField<Type, PatchField, GeoMesh >> &
200 PtrList<GeometricField<Type, PatchField, GeoMesh >>& modes, label Nmodes = 0,
201 bool consider_volumes =
true);
216 acquiredSnapshotsTimes, Eigen::MatrixXd parameters);
Header file of the Foam2Eigen class.
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...
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.