37template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
39 GeometricField<Type, PatchField, GeoMesh>& snapshot,
40 double _tol, word _PODnorm)
42 Info <<
"WARNING: the projection of the BC has not been implemented yet!" <<
47 PODnorm ==
"Frobenius",
"The PODnorm can be only L2 or Frobenius");
51template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
53 GeometricField<Type, PatchField, GeoMesh>& snapshot)
55 Info <<
"Initializing the incremental POD" << endl;
56 M_Assert(
tolleranceSVD > 0,
"Set up the tollerance before initialization");
57 M_Assert(
rank == 0,
"POD already initialized");
69 GeometricField<Type, PatchField, GeoMesh> tmp(snapshot);
71 this->append(tmp.clone());
82 this->append(snapshot / snapNorm);
92 Info <<
"Initialization ended" << endl;
95template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
97 GeometricField<Type, PatchField, GeoMesh>& snapshot)
99 Info <<
"********************************************************************"
101 Info <<
"Adding a snapshot" << endl;
102 Info <<
"Initial rank = " <<
rank << endl;
104 double projectionError;
105 double relProjectionError;
112 Eigen::VectorXd redCoeff =
project(snapshot);
113 GeometricField<Type, PatchField, GeoMesh> projectSnap =
projectSnapshot(
115 GeometricField<Type, PatchField, GeoMesh> orthoSnap(snapshot - projectSnap);
123 else if (
PODnorm ==
"Frobenius")
125 projectionError = snapshotEig.transpose() * snapshotEig -
126 (redCoeff.transpose() * redCoeff)(0, 0);
127 projectionError = Foam::sqrt(std::abs(projectionError));
128 relProjectionError = projectionError / (snapshotEig.transpose() * snapshotEig);
131 Info <<
"Relative projection error = " << relProjectionError << endl;
132 Eigen::MatrixXd Q = Eigen::MatrixXd::Zero(
rank + 1,
rank + 1);
134 Q.topRightCorner(
rank, 1) = redCoeff;
139 Info <<
"Projection error is smaller than the tollerance" << endl;
143 Eigen::JacobiSVD<Eigen::MatrixXd> svd(Q,
144 Eigen::ComputeThinU | Eigen::ComputeThinV);
145 Eigen::MatrixXd newModes = svd.matrixU();
146 Eigen::VectorXd newSingVal = svd.singularValues();
148 if (relProjectionError < tolleranceSVD || rank >= snapshotEig.size())
150 Info <<
"Projection error is small or matrix is full rank." <<
151 endl <<
"SVD rank held constant" << endl;
157 Info <<
"New snapshot is not in the span of the POD space." << endl
158 <<
"SVD rank increases" << endl;
159 orthogonalSnap = orthogonalSnap / projectionError;
160 Eigen::MatrixXd temp(this->
EigenModes[0].rows(),
163 this->EigenModes[0] = temp * newModes;
168 Info <<
"New POD rank = " <<
rank << endl;
169 double orthogonalPar;
173 orthogonalPar = std::abs(
175 massVector.asDiagonal() * this->EigenModes[0].col(0));
177 else if (
PODnorm ==
"Frobenius")
179 orthogonalPar = std::abs(this->
EigenModes[0].col(0).transpose() *
183 double EPS = 2.2204e-16;
184 Info <<
"Orthogonality = " << orthogonalPar << endl;
188 Info <<
"Orthogonalization required" << endl;
194 else if (
PODnorm ==
"Frobenius")
196 Eigen::HouseholderQR<Eigen::MatrixXd> qr(this->
EigenModes[0]);
197 Eigen::MatrixXd thinQ(Eigen::MatrixXd::Identity(this->
EigenModes[0].rows(),
199 this->
EigenModes[0] = qr.householderQ() * thinQ;
204 Info <<
"********************************************************************"
208template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
213 for (label i = 0; i <
rank; i++)
215 GeometricField<Type, PatchField, GeoMesh> tmp(this->first().name(),
217 Eigen::VectorXd vec = this->
EigenModes[0].col(i);
219 this->set(i, tmp.clone());
223template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
225 GeometricField<Type, PatchField, GeoMesh>& inputField,
229 Eigen::VectorXd projField;
231 if (numberOfModes == 0)
235 projField = this->
EigenModes[0].transpose() * fieldEig;
237 else if (
PODnorm ==
"Frobenius")
239 projField = this->
EigenModes[0].transpose() * fieldEig;
244 M_Assert(numberOfModes <= this->
EigenModes[0].cols(),
245 "Number of required modes for projection is higher then the number of available ones");
249 projField = (this->
EigenModes[0].leftCols(numberOfModes)).transpose() *
252 else if (
PODnorm ==
"Frobenius")
254 projField = (this->
EigenModes[0].leftCols(numberOfModes)).transpose() *
262template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
263GeometricField<Type, PatchField, GeoMesh>
265 GeometricField<Type, PatchField, GeoMesh>& inputField,
266 Eigen::MatrixXd Coeff,
274 label Nmodes = Coeff.rows();
275 Eigen::VectorXd InField = this->
EigenModes[0].leftCols(Nmodes) * Coeff;
277 inputField.rename(Name);
278 Info <<
"WARNING: Boundary conditions are not reconstructed using incremental POD"
283template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
284GeometricField<Type, PatchField, GeoMesh>
286 GeometricField<Type, PatchField, GeoMesh>& snapshot,
289 Eigen::MatrixXd
Modes;
290 Eigen::MatrixXd M_vol;
292 Eigen::MatrixXd projSnapI;
293 Eigen::MatrixXd projSnapCoeff;
295 if (numberOfModes == 0)
304 GeometricField<Type, PatchField, GeoMesh> Fr = snapshot;
312 else if (
PODnorm ==
"Frobenius")
315 projSnapI =
Modes.transpose() * F_eigen;
319 Foam::Info <<
"Inner product not defined" << Foam::endl;
323 projSnapCoeff = M.fullPivLu().solve(projSnapI);
328template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
330 PtrList<GeometricField<Type, PatchField, GeoMesh >> snapshots,
331 PtrList<GeometricField<Type, PatchField, GeoMesh >>& projSnapshots,
334 M_Assert(numberOfModes <= this->size(),
335 "The number of Modes used for the projection cannot be bigger than the number of available modes");
336 projSnapshots.resize(snapshots.size());
337 Eigen::MatrixXd
Modes;
339 if (numberOfModes == 0)
348 Eigen::MatrixXd M_vol;
350 Eigen::MatrixXd projSnapI;
351 Eigen::MatrixXd projSnapCoeff;
353 for (label i = 0; i < snapshots.size(); i++)
355 GeometricField<Type, PatchField, GeoMesh> Fr = snapshots[0];
363 else if (
PODnorm ==
"Frobenius")
366 projSnapI =
Modes.transpose() * F_eigen;
370 Foam::Info <<
"Inner product not defined" << Foam::endl;
374 projSnapCoeff = M.fullPivLu().solve(projSnapI);
376 projSnapshots.set(i, Fr.clone());
382template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
Header file of the EigenFunctions class.
static GeometricField< scalar, PatchField, GeoMesh > Eigen2field(GeometricField< scalar, PatchField, GeoMesh > &field, Eigen::VectorXd &eigen_vector, bool correctBC=true)
Convert a vector in Eigen format into an OpenFOAM scalar GeometricField.
static Eigen::MatrixXd field2Eigen(GeometricField< Type, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
Implementation of a container class derived from PtrList.
List< Eigen::MatrixXd > EigenModes
List of Matrices that contains the internalField and the additional matrices for the boundary patches...
List< Eigen::MatrixXd > toEigen()
Method that convert a PtrList of modes into Eigen matrices filling the EigenModes object.
PtrList< GeometricField< Type, PatchField, GeoMesh > > & toPtrList()
Function that returns the Modes object as a standard PtrList.
Implementation of a incremental POD algorithm according to Oxberry et al.
void addSnapshot(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Add a snapshot to the POD space.
void initialize(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Initialize the incremental POD algorithm.
void fillPtrList()
Fill the POD modes prtList from the Eigen matrix.
Eigen::VectorXd singularValues
Vector of the singular values.
void writeModes()
Write to modes to file.
void projectSnapshots(PtrList< GeometricField< Type, PatchField, GeoMesh > > snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &projSnapshots, label numberOfModes=0)
Project a list of snapshots onto the reduced space.
Eigen::VectorXd project(GeometricField< Type, PatchField, GeoMesh > &inputField, label numberOfModes=0)
Project a field into the reduced space.
GeometricField< Type, PatchField, GeoMesh > reconstruct(GeometricField< Type, PatchField, GeoMesh > &inputField, Eigen::MatrixXd Coeff, word Name)
Reconstruct a field from the reduced coefficients.
word outputFolder
Folder where the modes and singular values are saved.
incrementalPOD()
Constructors.
double tolleranceSVD
Tollerance for the incremental SVD.
Eigen::VectorXd massVector
Vector of the cells volumes.
int rank
Rank of the reduced space.
GeometricField< Type, PatchField, GeoMesh > projectSnapshot(GeometricField< Type, PatchField, GeoMesh > &snapshot, label numberOfModes=0)
Project a snapshot onto the POD space.
word PODnorm
POD norm used for projection. It can be L2 or Frobenius.
Header file of the incrementalPOD class.
void weightedGramSchmidt(Eigen::MatrixXd &matrix, Eigen::VectorXd &weights)
Weighted Gram-Schmidt orthogonalization.
void exportFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &field, word folder, word fieldname)
Function to export a scalar of vector field.
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
Eigen::VectorXd getMassMatrixFV(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Gets a vector containing the volumes of each cell of the mesh.
double L2Norm(const T v)
Compute the L2 norm of v.