37template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
42 for (label
i = 0;
i < (this->first()).boundaryFieldRef().size();
i++)
44 if ((this->first()).boundaryFieldRef()[
i].type() !=
"processor")
50 EigenModes.resize(NBC + 1);
54 for (label
i = 0;
i < NBC;
i++)
56 EigenModes[
i + 1] = BC[
i];
62template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
64 fvMatrix<Type>& Af, label numberOfModes,
67 M_Assert(projType ==
"G" || projType ==
"PG",
68 "Projection type can be G for Galerkin or PG for Petrov-Galerkin");
69 List<Eigen::MatrixXd> LinSys;
72 if (EigenModes.size() == 0)
77 Eigen::SparseMatrix<double> Ae;
81 if (numberOfModes == 0)
85 LinSys[0] = EigenModes[0].transpose() * Ae * EigenModes[0];
86 LinSys[1] = EigenModes[0].transpose() * be;
91 LinSys[0] = (Ae * EigenModes[0]).transpose() * Ae * EigenModes[0];
92 LinSys[1] = (Ae * EigenModes[0]).transpose() * be;
97 M_Assert(numberOfModes <= EigenModes[0].cols(),
98 "Number of required modes for projection is higher then the number of available ones");
102 LinSys[0] = ((EigenModes[0]).leftCols(numberOfModes)).transpose() * Ae *
103 (EigenModes[0]).leftCols(numberOfModes);
104 LinSys[1] = ((EigenModes[0]).leftCols(numberOfModes)).transpose() * be;
107 if (projType ==
"PG")
109 LinSys[0] = (Ae * ((EigenModes[0]).leftCols(numberOfModes))).transpose() * Ae *
110 (EigenModes[0]).leftCols(numberOfModes);
111 LinSys[1] = (Ae * ((EigenModes[0]).leftCols(numberOfModes))).transpose() * be;
118template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
120 GeometricField<Type, PatchField, GeoMesh>&
121 field, label numberOfModes, word projType, fvMatrix<Type>* Af)
123 M_Assert(projType ==
"F" || projType ==
"G" || projType ==
"PG",
124 "Projection type can be F for Frobenius, G for Galerkin or PG for Petrov-Galerkin");
127 Eigen::MatrixXd projField;
131 if (EigenModes.size() == 0)
136 if (numberOfModes == 0)
140 vol = Eigen::VectorXd::Ones(vol.size());
141 b = EigenModes[0].transpose() * vol.asDiagonal() * fieldEig;
142 M = EigenModes[0].transpose() * EigenModes[0].transpose();
143 projField = M.fullPivLu().solve(b);
145 else if (projType ==
"G")
147 projField = EigenModes[0].transpose() * vol.asDiagonal() * fieldEig;
149 else if (projType ==
"PG")
152 "Using a Petrov-Galerkin projection you have to provide also the system matrix");
153 Eigen::SparseMatrix<double> Ae;
156 projField = (Ae * EigenModes[0]).transpose() * vol.asDiagonal() * fieldEig;
161 M_Assert(numberOfModes <= EigenModes[0].cols(),
162 "Number of required modes for projection is higher then the number of available ones");
166 vol = Eigen::VectorXd::Ones(vol.size());
167 M = EigenModes[0].leftCols(numberOfModes).transpose() * EigenModes[0].leftCols(
169 b = ((EigenModes[0]).leftCols(numberOfModes)).transpose() *
170 vol.asDiagonal() * fieldEig;
171 projField = M.fullPivLu().solve(b);
173 else if (projType ==
"G")
175 projField = ((EigenModes[0]).leftCols(numberOfModes)).transpose() *
176 vol.asDiagonal() * fieldEig;
178 else if (projType ==
"PG")
181 "Using a Petrov-Galerkin projection you have to provide also the system matrix");
182 Eigen::SparseMatrix<double> Ae;
185 projField = (Ae * ((EigenModes[0]).leftCols(numberOfModes))).transpose() *
186 vol.asDiagonal() * fieldEig;
193template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
194GeometricField<Type, PatchField, GeoMesh>
196 GeometricField<Type, PatchField, GeoMesh>&
197 field, label numberOfModes, word projType, fvMatrix<Type>* Af)
199 Eigen::MatrixXd proj = project(field, numberOfModes, projType, Af);
200 GeometricField<Type, PatchField, GeoMesh> projSnap = field;
201 reconstruct(projSnap, proj, projSnap.name());
205template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
207 PtrList<GeometricField<Type, PatchField, GeoMesh>>&
209 label numberOfModes, word projType, fvMatrix<Type>* Af)
212 "Projection type can be G for Galerking or PG for Petrov-Galerkin");
215 Eigen::MatrixXd projField;
217 if (EigenModes.size() == 0)
222 if (numberOfModes == 0)
226 vol = Eigen::VectorXd::Ones(vol.size());
227 projField = EigenModes[0].transpose() * vol.asDiagonal() * fieldEig;
229 else if (projType ==
"G")
231 projField = EigenModes[0].transpose() * vol.asDiagonal() * fieldEig;
233 else if (projType ==
"PG")
236 "Using a Petrov-Galerkin projection you have to provide also the system matrix");
237 Eigen::SparseMatrix<double> Ae;
240 projField = (Ae * EigenModes[0]).transpose() * vol.asDiagonal() * fieldEig;
246 "Number of required modes for projection is higher then the number of available ones");
250 vol = Eigen::VectorXd::Ones(vol.size());
251 projField = ((EigenModes[0]).leftCols(numberOfModes)).transpose() *
252 vol.asDiagonal() * fieldEig;
254 else if (projType ==
"G")
256 projField = ((EigenModes[0]).leftCols(numberOfModes)).transpose() *
257 vol.asDiagonal() * fieldEig;
259 else if (projType ==
"PG")
262 "Using a Petrov-Galerkin projection you have to provide also the system matrix");
263 Eigen::SparseMatrix<double> Ae;
266 projField = (Ae * ((EigenModes[0]).leftCols(numberOfModes))).transpose() *
267 vol.asDiagonal() * fieldEig;
273template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
274GeometricField<Type, PatchField, GeoMesh>
276 GeometricField<Type, PatchField, GeoMesh>& inputField,
277 Eigen::MatrixXd Coeff,
280 if (EigenModes.size() == 0)
285 label Nmodes = Coeff.rows();
286 Eigen::VectorXd InField = EigenModes[0].leftCols(Nmodes) * Coeff;
288 if (inputField.name() ==
"nut")
290 InField = (InField.array() < 0).select(0, InField);
294 inputField.rename(Name);
296 for (label
i = 0;
i < NBC;
i++)
298 Eigen::VectorXd BF = EigenModes[
i + 1].leftCols(Nmodes) * Coeff;
300 if (inputField.name() ==
"nut")
302 BF = (BF.array() < 0).select(0, BF);
310template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
311PtrList<GeometricField<Type, PatchField, GeoMesh>>
313 GeometricField<Type, PatchField, GeoMesh>& inputField,
314 List < Eigen::MatrixXd> Coeff,
317 PtrList<GeometricField<Type, PatchField, GeoMesh>> inputFields;
318 inputFields.resize(0);
320 for (label
i = 0;
i < Coeff.size();
i++)
322 inputField = reconstruct(inputField, Coeff[
i], Name);
323 inputFields.append(inputField.clone());
330template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
332 PtrList<GeometricField<Type, PatchField, GeoMesh>> snapshots,
333 PtrList<GeometricField<Type, PatchField, GeoMesh>>& projSnapshots,
334 PtrList<volScalarField> Volumes,
338 if (EigenModes.size() == 0)
343 M_Assert(snapshots.size() == Volumes.size(),
344 "The number of snapshots and the number of volumes vectors must be equal");
345 M_Assert(numberOfModes <= this->size(),
346 "The number of Modes used for the projection cannot be bigger than the number of available modes");
347 M_Assert(innerProduct ==
"L2" || innerProduct ==
"Frobenius",
348 "The chosen inner product is not implemented");
349 projSnapshots.resize(snapshots.size());
350 label dim = std::nearbyint(EigenModes[0].rows() /
352 Eigen::MatrixXd totVolumes(Volumes[0].size()*dim, Volumes.size());
354 for (label
i = 0;
i < Volumes.size();
i++)
359 totVolumes.replicate(dim, 1);
360 Eigen::MatrixXd
Modes;
362 if (numberOfModes == 0)
364 Modes = EigenModes[0];
368 Modes = EigenModes[0].leftCols(numberOfModes);
372 Eigen::MatrixXd projSnapI;
373 Eigen::MatrixXd projSnapCoeff;
375 for (label
i = 0;
i < snapshots.size();
i++)
377 GeometricField<Type, PatchField, GeoMesh> Fr = snapshots[0];
380 if (innerProduct ==
"L2")
382 M =
Modes.transpose() * (totVolumes.col(
i)).asDiagonal() *
Modes;
383 projSnapI =
Modes.transpose() * (totVolumes.col(
i)).asDiagonal() * F_eigen;
388 projSnapI =
Modes.transpose() * F_eigen;
391 projSnapCoeff = M.fullPivLu().solve(projSnapI);
392 reconstruct(Fr, projSnapCoeff,
"projSnap");
393 projSnapshots.set(
i, Fr.clone());
396template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
398 PtrList<GeometricField<Type, PatchField, GeoMesh>> snapshots,
399 PtrList<GeometricField<Type, PatchField, GeoMesh>>& projSnapshots,
400 PtrList<volScalarField> Volumes, word innerProduct)
402 label numberOfModes = 0;
403 projectSnapshots(snapshots, projSnapshots, Volumes, numberOfModes,
407template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
409 PtrList<GeometricField<Type, PatchField, GeoMesh>> snapshots,
410 PtrList<GeometricField<Type, PatchField, GeoMesh>>& projSnapshots,
414 if (EigenModes.size() == 0)
419 M_Assert(numberOfModes <= this->size(),
420 "The number of Modes used for the projection cannot be bigger than the number of available modes");
421 M_Assert(innerProduct ==
"L2" || innerProduct ==
"Frobenius",
422 "The chosen inner product is not implemented");
423 projSnapshots.resize(snapshots.size());
424 Eigen::MatrixXd
Modes;
426 if (numberOfModes == 0)
428 Modes = EigenModes[0];
432 Modes = EigenModes[0].leftCols(numberOfModes);
435 Eigen::MatrixXd M_vol;
437 Eigen::MatrixXd projSnapI;
438 Eigen::MatrixXd projSnapCoeff;
440 for (label
i = 0;
i < snapshots.size();
i++)
442 GeometricField<Type, PatchField, GeoMesh> Fr = snapshots[0];
445 if (innerProduct ==
"L2")
449 else if (innerProduct ==
"Frobenius")
451 M_vol = Eigen::VectorXd::Identity(F_eigen.rows(), 1);
455 std::cout <<
"Inner product not defined" << endl;
459 M =
Modes.transpose() * M_vol.asDiagonal() *
Modes;
460 projSnapI =
Modes.transpose() * M_vol.asDiagonal() * F_eigen;
461 projSnapCoeff = M.fullPivLu().solve(projSnapI);
462 reconstruct(Fr, projSnapCoeff,
"projSnap");
463 projSnapshots.set(
i, Fr.clone());
466template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
468 PtrList<GeometricField<Type, PatchField, GeoMesh>> snapshots,
469 PtrList<GeometricField<Type, PatchField, GeoMesh>>& projSnapshots,
472 label numberOfModes = 0;
473 projectSnapshots(snapshots, projSnapshots, numberOfModes, innerProduct);
476template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
478 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes)
480 this->
resize(modes.size());
482 for (label
i = 0;
i < modes.size();
i++)
484 (*this).set(
i, modes[
i].clone());
#define M_Assert(Expr, Msg)
Header file of the Modes 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::VectorXd field2Eigen(GeometricField< tensor, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
static List< Eigen::MatrixXd > PtrList2EigenBC(PtrList< GeometricField< scalar, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert an OpenFOAM scalar field to a List of Eigen Vectors, one for each boundary.
static Eigen::MatrixXd PtrList2Eigen(PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert a PtrList of snapshots to Eigen matrix (only internal field)
static void fvMatrix2Eigen(fvMatrix< type_foam_matrix > foam_matrix, type_A &A, type_B &b)
Convert a FvMatrix OpenFOAM matrix (Linear System) into a Eigen Matrix A and a source vector b.
Implementation of a container class derived from PtrList.
GeometricField< Type, PatchField, GeoMesh > reconstruct(GeometricField< Type, PatchField, GeoMesh > &inputField, Eigen::MatrixXd Coeff, word Name)
Function to reconstruct the solution starting from the coefficients, in this case the field is passed...
List< Eigen::MatrixXd > toEigen()
Method that convert a PtrList of modes into Eigen matrices filling the EigenModes object.
List< Eigen::MatrixXd > project(fvMatrix< Type > &Af, label numberOfModes=0, word projType="G")
A function that project an FvMatrix (OpenFoam linear System) on the modes.
void projectSnapshots(PtrList< GeometricField< Type, PatchField, GeoMesh > > snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &projSnapshots, label numberOfModes=0, word innerProduct="L2")
Function to project a list of fields into the modes manifold.
GeometricField< Type, PatchField, GeoMesh > projectSnapshot(GeometricField< Type, PatchField, GeoMesh > &field, label numberOfModes=0, word projType="G", fvMatrix< Type > *Af=NULL)
A function that project and reconstruct a snapshot on the modes.
void operator=(const PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes)
Eigen::VectorXd getMassMatrixFV(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Gets a vector containing the volumes of each cell of the mesh.
void assignBC(GeometricField< scalar, fvPatchField, volMesh > &s, label BC_ind, double value)
Assign uniform Boundary Condition to a volScalarField.