Header file of the ITHACAcoeffsMass file. More...
#include "fvCFD.H"
#include "IOmanip.H"
#include "freestreamFvPatchField.H"
#include <sys/stat.h>
#include <unistd.h>
#include <Eigen/Eigen>
#include <functional>
#include "./colormod.H"
#include "polyMeshTools.H"
#include <chrono>
#include "mixedFvPatchFields.H"
#include "fvMeshSubset.H"
#include "Foam2Eigen.H"
#include "ITHACAstream.H"
Go to the source code of this file.
Namespaces | |
namespace | ITHACAutilities |
Namespace to implement some useful assign operation of OF fields. | |
Functions | |
template<class Type , template< class > class PatchField, class GeoMesh > | |
PtrList< GeometricField< Type, PatchField, GeoMesh > > | ITHACAutilities::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 coefficients. | |
template<class Type , template< class > class PatchField, class GeoMesh > | |
Eigen::MatrixXd | ITHACAutilities::getMassMatrix (PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label nModes=0, bool consider_volumes=true) |
Gets the mass matrix using the eigen routine. | |
template<class Type , template< class > class PatchField, class GeoMesh > | |
Eigen::MatrixXd | ITHACAutilities::getMassMatrix (PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields2, label nModes=0, bool consider_volumes=true) |
Gets the cross mass matrix using the eigen routine. | |
template<class Type , template< class > class PatchField, class GeoMesh > | |
Eigen::MatrixXd | ITHACAutilities::getMassMatrix (PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields2, Eigen::VectorXd weights, label Nmodes=0, bool consider_volumes=true) |
Gets the cross mass matrix using the eigen routine. | |
template<class Type , template< class > class PatchField, class GeoMesh > | |
Eigen::VectorXd | ITHACAutilities::getMassMatrixFV (GeometricField< Type, PatchField, GeoMesh > &snapshot) |
Gets a vector containing the volumes of each cell of the mesh. | |
template<class Type , template< class > class PatchField, class GeoMesh > | |
Eigen::VectorXd | ITHACAutilities::getCoeffs (GeometricField< Type, PatchField, GeoMesh > &snapshot, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes=0, bool consider_volumes=true) |
Projects a snapshot on a basis function and gets the coefficients of the projection. | |
template<class Type , template< class > class PatchField, class GeoMesh > | |
Eigen::MatrixXd | ITHACAutilities::getCoeffs (PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshot, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes=0, bool consider_volumes=true) |
Projects snapshots on a basis function and gets the coefficients of the projection. | |
Eigen::MatrixXd | ITHACAutilities::parTimeCombMat (List< Eigen::VectorXd > acquiredSnapshotsTimes, Eigen::MatrixXd parameters) |
A method to compute the time-parameter combined matrix whose any single element corresponds to a unique snapshot in the snapshots acquired for the offline stage. | |
Header file of the ITHACAcoeffsMass file.
Definition in file ITHACAcoeffsMass.H.