40template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
42 PtrList<GeometricField<Type, PatchField, GeoMesh>>& snapshots,
43 PtrList<GeometricField<Type, PatchField, GeoMesh>>& ModesGlobal,
45 label Npar, label NnestedOut)
48 List<PtrList<GeometricField<Type, PatchField, GeoMesh>>>
50 label Nt = snapshots.size() / Npar;
51 SnapMatrixNested.setSize(Nt);
53 for (label
i = 0;
i < Npar;
i++)
55 SnapMatrixNested[
i].resize(Nt);
57 for (label j = 0; j < Nt; j++)
59 SnapMatrixNested[
i].set(j, snapshots[j + Nt *
i].clone());
63 List<PtrList<GeometricField<Type, PatchField, GeoMesh>>> UModesNested;
64 PtrList<GeometricField<Type, PatchField, GeoMesh>> y;
65 UModesNested.setSize(Nt);
67 for (label
i = 0;
i < Npar;
i++)
73 for (label
i = 0;
i < Npar;
i++)
77 for (label j = 0; j < y.size(); j++)
79 ModesGlobal.append(y[j].clone());
85 PtrList<volScalarField>& snapshots, PtrList<volScalarField>& ModesGlobal,
86 word fieldName, label Npar, label NnestedOut);
89 PtrList<volVectorField>& snapshots, PtrList<volVectorField>& ModesGlobal,
90 word fieldName, label Npar, label NnestedOut);
92template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
94 PtrList<GeometricField<Type, PatchField, GeoMesh>>& snapshots,
95 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes,
96 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
100 word PODkey =
"POD_" + fieldName;
101 word PODnorm = para->
ITHACAdict->lookupOrDefault<word>(PODkey,
"L2");
103 PODnorm ==
"Frobenius",
"The PODnorm can be only L2 or Frobenius");
104 Info <<
"Performing POD for " << fieldName <<
" using the " << PODnorm <<
107 if ((podex == 0 && sup == 0) || (supex == 0 && sup == 1))
113 nmodes = snapshots.size() - 2;
116 M_Assert(nmodes <= snapshots.size() - 2,
117 "The number of requested modes cannot be bigger than the number of Snapshots - 2");
123 nmodes = snapshots.size();
126 M_Assert(nmodes <= snapshots.size(),
127 "The number of requested modes cannot be bigger than the number of Snapshots");
132 label NBC = snapshots[0].boundaryField().size();
133 Eigen::MatrixXd _corMatrix;
139 else if (PODnorm ==
"Frobenius")
144 if (Pstream::parRun())
146 reduce(_corMatrix, sumOp<Eigen::MatrixXd>());
149 Eigen::VectorXd eigenValueseig;
150 Eigen::MatrixXd eigenVectoreig;
151 modes.resize(nmodes);
152 Info <<
"####### Performing the POD using EigenDecomposition " <<
153 fieldName <<
" #######" << endl;
154 label ncv = snapshots.size();
155 Spectra::DenseSymMatProd<double> op(_corMatrix);
156 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> esEg;
160 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra::DenseSymMatProd<double>>
161 es(&op, nmodes, ncv);
162 std::cout <<
"Using Spectra EigenSolver " << std::endl;
164 es.compute(1000, 1e-10, Spectra::LARGEST_ALGE);
165 M_Assert(es.info() == Spectra::SUCCESSFUL,
166 "The Eigenvalue Decomposition did not succeed");
167 eigenVectoreig = es.eigenvectors().real();
168 eigenValueseig = es.eigenvalues().real();
172 std::cout <<
"Using Eigen EigenSolver " << std::endl;
173 esEg.compute(_corMatrix);
174 M_Assert(esEg.info() == Eigen::Success,
175 "The Eigenvalue Decomposition did not succeed");
176 eigenVectoreig = esEg.eigenvectors().real().rowwise().reverse().leftCols(
178 eigenValueseig = esEg.eigenvalues().real().array().reverse();
181 if (eigenValueseig.array().minCoeff() < 0)
183 eigenValueseig = eigenValueseig.array() + 2 * abs(
184 eigenValueseig.array().minCoeff());
187 Info <<
"####### End of the POD for " << snapshots[0].name() <<
" #######" <<
193 Eigen::MatrixXd modesEig = (SnapMatrix * eigenVectoreig);
196 Eigen::MatrixXd normFact(nmodes, 1);
198 for (label
i = 0;
i < nmodes;
i++)
202 normFact(
i, 0) = std::sqrt((modesEig.col(
i).transpose() * V.asDiagonal() *
203 modesEig.col(
i))(0, 0));
205 if (Pstream::parRun())
207 normFact(
i, 0) = (modesEig.col(
i).transpose() * V.asDiagonal() *
208 modesEig.col(
i))(0, 0);
211 else if (PODnorm ==
"Frobenius")
213 normFact(
i, 0) = std::sqrt((modesEig.col(
i).transpose() * modesEig.col(
i))(0,
216 if (Pstream::parRun())
218 normFact(
i, 0) = (modesEig.col(
i).transpose() * modesEig.col(
i))(0, 0);
223 if (Pstream::parRun())
225 reduce(normFact, sumOp<Eigen::MatrixXd>());
228 if (Pstream::parRun())
230 normFact = normFact.cwiseSqrt();
233 List<Eigen::MatrixXd> modesEigBC;
234 modesEigBC.resize(NBC);
236 for (label
i = 0;
i < NBC;
i++)
238 modesEigBC[
i] = (SnapMatrixBC[
i] * eigenVectoreig);
241 std::cout << normFact << std::endl;
243 for (label
i = 0;
i < nmodes;
i++)
245 modesEig.col(
i) = modesEig.col(
i).array() / normFact(
i, 0);
247 for (label j = 0; j < NBC; j++)
249 modesEigBC[j].col(
i) = modesEigBC[j].col(
i).array() / normFact(
i, 0);
253 for (label
i = 0;
i < modes.size();
i++)
255 GeometricField<Type, PatchField, GeoMesh> tmp2(snapshots[0].name(),
257 Eigen::VectorXd vec = modesEig.col(
i);
260 for (label k = 0; k < NBC; k++)
265 modes.set(
i, tmp2.clone());
268 eigenValueseig = eigenValueseig / eigenValueseig.sum();
269 Eigen::VectorXd cumEigenValues(eigenValueseig);
271 for (label j = 1; j < cumEigenValues.size(); ++j)
273 cumEigenValues(j) += cumEigenValues(j - 1);
276 Info <<
"####### Saving the POD bases for " << snapshots[0].name() <<
282 snapshots[0].name());
290 "./ITHACAoutput/POD/Eigenvalues_" + snapshots[0].name(), para->
precision,
293 "./ITHACAoutput/POD/CumEigenvalues_" + snapshots[0].name(), para->
precision,
298 Info <<
"Reading the existing modes" << endl;
303 "./ITHACAoutput/supremizer/");
313 PtrList<volVectorField>& snapshots, PtrList<volVectorField>& modes,
314 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
318 PtrList<volScalarField>& snapshots, PtrList<volScalarField>& modes,
319 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
323 PtrList<surfaceScalarField>& snapshots, PtrList<surfaceScalarField>& modes,
324 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
327template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
329 PtrList<GeometricField<Type, PatchField, GeoMesh>>& snapshots,
330 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes,
331 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
338 nmodes = snapshots.size() - 2;
341 M_Assert(nmodes <= snapshots.size() - 2,
342 "The number of requested modes cannot be bigger than the number of Snapshots - 2");
344 if ((podex == 0 && sup == 0) || (supex == 0 && sup == 1))
348 label NBC = snapshots[0].boundaryField().size();
350 Eigen::VectorXd eigenValueseig;
351 Eigen::MatrixXd eigenVectoreig;
352 modes.resize(nmodes);
353 Info <<
"####### Performing the POD using EigenDecomposition " <<
354 snapshots[0].name() <<
" #######" << endl;
355 label ncv = snapshots.size();
356 Spectra::DenseSymMatProd<double> op(_corMatrix);
357 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> esEg;
358 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra::DenseSymMatProd<double>>
359 es(&op, nmodes, ncv);
363 std::cout <<
"Using Spectra EigenSolver " << std::endl;
365 es.compute(1000, 1e-10, Spectra::LARGEST_ALGE);
366 M_Assert(es.info() == Spectra::SUCCESSFUL,
367 "The Eigenvalue Decomposition did not succeed");
368 eigenVectoreig = es.eigenvectors().real();
369 eigenValueseig = es.eigenvalues().real();
373 std::cout <<
"Using Eigen EigenSolver " << std::endl;
374 esEg.compute(_corMatrix);
375 M_Assert(esEg.info() == Eigen::Success,
376 "The Eigenvalue Decomposition did not succeed");
377 eigenVectoreig = esEg.eigenvectors().real().rowwise().reverse().leftCols(
379 eigenValueseig = esEg.eigenvalues().real().reverse();
382 Info <<
"####### End of the POD for " << snapshots[0].name() <<
" #######" <<
384 Eigen::VectorXd eigenValueseigLam =
385 eigenValueseig.real().array().cwiseInverse().sqrt() ;
386 Eigen::VectorXd eigenValueseigWeigted = eigenValueseig.head(
387 nmodes).real().array() ;
388 Eigen::MatrixXd modesEig = (SnapMatrix * eigenVectoreig) *
389 eigenValueseigLam.head(nmodes).asDiagonal() *
390 eigenValueseigWeigted.asDiagonal();
391 List<Eigen::MatrixXd> modesEigBC;
392 modesEigBC.resize(NBC);
394 for (label
i = 0;
i < NBC;
i++)
396 modesEigBC[
i] = (SnapMatrixBC[
i] * eigenVectoreig) *
397 eigenValueseigLam.asDiagonal() * eigenValueseigWeigted.asDiagonal();
400 for (label
i = 0;
i < modes.size();
i++)
402 GeometricField<Type, PatchField, GeoMesh> tmp2(snapshots[0].name(),
404 Eigen::VectorXd vec = modesEig.col(
i);
407 for (label k = 0; k < tmp2.boundaryField().size(); k++)
412 modes.set(
i, tmp2.clone());
415 eigenValueseig = eigenValueseig / eigenValueseig.sum();
416 Eigen::VectorXd cumEigenValues(eigenValueseig);
418 for (label j = 1; j < cumEigenValues.size(); ++j)
420 cumEigenValues(j) += cumEigenValues(j - 1);
423 Info <<
"####### Saving the POD bases for " << snapshots[0].name() <<
430 snapshots[0].name());
438 "./ITHACAoutput/POD/Eigenvalues_" + snapshots[0].name(), para->
precision,
441 "./ITHACAoutput/POD/CumEigenvalues_" + snapshots[0].name(), para->
precision,
446 Info <<
"Reading the existing modes" << endl;
451 "./ITHACAoutput/supremizer/");
461 PtrList<volScalarField>& snapshots, PtrList<volScalarField>& modes,
462 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
466 PtrList<volVectorField>& snapshots, PtrList<volVectorField>& modes,
467 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
470template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
472 PtrList<GeometricField<Type, PatchField, GeoMesh>>& snapshots,
473 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes,
474 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
479 if ((podex == 0 && sup == 0) || (supex == 0 && sup == 1))
481 PtrList<volVectorField> Bases;
482 modes.resize(nmodes);
483 Info <<
"####### Performing POD using Singular Value Decomposition for " <<
484 snapshots[0].name() <<
" #######" << endl;
487 Eigen::VectorXd V3dSqrt = V.array().sqrt();
488 Eigen::VectorXd V3dInv = V3dSqrt.array().cwiseInverse();
489 auto VMsqr = V3dSqrt.asDiagonal();
490 auto VMsqrInv = V3dInv.asDiagonal();
491 Eigen::MatrixXd SnapMatrix2 = VMsqr * SnapMatrix;
492 Eigen::JacobiSVD<Eigen::MatrixXd>
svd(SnapMatrix2,
493 Eigen::ComputeThinU | Eigen::ComputeThinV);
494 Info <<
"####### End of the POD for " << snapshots[0].name() <<
" #######" <<
496 Eigen::VectorXd eigenValueseig;
497 Eigen::MatrixXd eigenVectoreig;
498 eigenValueseig =
svd.singularValues().real();
499 eigenVectoreig =
svd.matrixU().real();
500 Eigen::MatrixXd modesEig = VMsqrInv * eigenVectoreig;
501 GeometricField<Type, PatchField, GeoMesh> tmb_bu(snapshots[0].name(),
504 for (label
i = 0;
i < nmodes;
i++)
506 Eigen::VectorXd vec = modesEig.col(
i);
508 modes.set(
i, tmb_bu.clone());
511 eigenValueseig = eigenValueseig / eigenValueseig.sum();
512 Eigen::VectorXd cumEigenValues(eigenValueseig);
514 for (label j = 1; j < cumEigenValues.size(); ++j)
516 cumEigenValues(j) += cumEigenValues(j - 1);
519 Info <<
"####### Saving the POD bases for " << snapshots[0].name() <<
525 snapshots[0].name());
534 "./ITHACAoutput/POD/Eigenvalues_" + snapshots[0].name(), para->
precision,
537 "./ITHACAoutput/POD/CumEigenvalues_" + snapshots[0].name(), para->
precision,
542 Info <<
"Reading the existing modes" << endl;
547 "./ITHACAoutput/supremizer/");
557 PtrList<volScalarField>& snapshots, PtrList<volScalarField>& modes,
558 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
562 PtrList<volVectorField>& snapshots, PtrList<volVectorField>& modes,
563 word fieldName,
bool podex,
bool supex,
bool sup, label nmodes,
568Eigen::MatrixXd
corMatrix(PtrList<volScalarField>& snapshots)
570 Info <<
"########## Filling the correlation matrix for " << snapshots[0].name()
571 <<
"##########" << endl;
572 Eigen::MatrixXd matrix( snapshots.size(), snapshots.size());
574 for (label
i = 0;
i < snapshots.size();
i++)
576 for (label j = 0; j <=
i; j++)
578 matrix(
i, j) = fvc::domainIntegrate(snapshots[
i] * snapshots[j]).value();
582 for (label
i = 1;
i < snapshots.size();
i++)
584 for (label j = 0; j <
i; j++)
586 matrix(j,
i) = matrix(
i, j);
596Eigen::MatrixXd
corMatrix(PtrList<volVectorField>& snapshots)
598 Info <<
"########## Filling the correlation matrix for " << snapshots[0].name()
599 <<
"##########" << endl;
600 Eigen::MatrixXd matrix( snapshots.size(), snapshots.size());
602 for (label
i = 0;
i < snapshots.size();
i++)
604 for (label j = 0; j <=
i; j++)
606 matrix(
i, j) = fvc::domainIntegrate(snapshots[
i] & snapshots[j]).value();
610 for (label
i = 1;
i < snapshots.size();
i++)
612 for (label j = 0; j <
i; j++)
614 matrix(j,
i) = matrix(
i, j);
623Eigen::MatrixXd
corMatrix(List<Eigen::SparseMatrix<double>>&
626 Info <<
"########## Filling the correlation matrix for the matrix list ##########"
628 Eigen::MatrixXd matrix( snapshots.size(), snapshots.size());
630 for (label
i = 0;
i < snapshots.size();
i++)
632 for (label j = 0; j <=
i; j++)
636 for (label k = 0; k < snapshots[
i].cols(); k++)
638 res += snapshots[
i].col(k).dot(snapshots[j].
col(k));
645 for (label
i = 1;
i < snapshots.size();
i++)
647 for (label j = 0; j <
i; j++)
649 matrix(j,
i) = matrix(
i, j);
658Eigen::MatrixXd
corMatrix(List<Eigen::VectorXd>& snapshots)
660 Info <<
"########## Filling the correlation matrix for the matrix list ##########"
662 Eigen::MatrixXd matrix( snapshots.size(), snapshots.size());
664 for (label
i = 0;
i < snapshots.size();
i++)
666 for (label j = 0; j <=
i; j++)
668 matrix(
i, j) = (snapshots[
i].transpose() * snapshots[j]).trace();
672 for (label
i = 1;
i < snapshots.size();
i++)
674 for (label j = 0; j <
i; j++)
676 matrix(j,
i) = matrix(
i, j);
686template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
687void exportBases(PtrList<GeometricField<Type, PatchField, GeoMesh>>& s,
688 PtrList<GeometricField<Type, PatchField, GeoMesh>>& bases,
689 word fieldName,
bool sup)
695 for (label
i = 0;
i < s.size();
i++)
697 mkDir(
"./ITHACAoutput/supremizer/" + name(
i + 1));
698 fieldname =
"./ITHACAoutput/supremizer/" + name(
i + 1) +
"/" +
700 OFstream os(fieldname);
701 bases[
i].writeHeader(os);
709 for (label
i = 0;
i < s.size();
i++)
711 mkDir(
"./ITHACAoutput/POD/" + name(
i + 1));
712 fieldname =
"./ITHACAoutput/POD/" + name(
i + 1) +
"/" + bases[
i].name();
713 OFstream os(fieldname);
714 bases[
i].writeHeader(os);
720 PtrList<volVectorField>& bases, word fieldName,
bool sup);
722 PtrList<volScalarField>& bases, word fieldName,
bool sup);
730 mkDir(
"./ITHACAoutput/supremizer/");
731 fieldname =
"./ITHACAoutput/supremizer/Eigenvalues_" + name;
732 OFstream os(fieldname);
734 for (label
i = 0;
i < Eigenvalues.size();
i++)
736 os << Eigenvalues[
i] << endl;
742 mkDir(
"./ITHACAoutput/POD/");
743 fieldname =
"./ITHACAoutput/POD/Eigenvalues_" + name;
744 OFstream os(fieldname);
746 for (label
i = 0;
i < Eigenvalues.size();
i++)
748 os << Eigenvalues[
i] << endl;
759 mkDir(
"./ITHACAoutput/supremizer/");
760 fieldname =
"./ITHACAoutput/supremizer/cumEigenvalues_" + name;
761 OFstream os(fieldname);
763 for (label
i = 0;
i < cumEigenvalues.size();
i++)
765 os << cumEigenvalues[
i] << endl;
771 mkDir(
"./ITHACAoutput/POD/");
772 fieldname =
"./ITHACAoutput/POD/cumEigenvalues_" + name;
773 OFstream os(fieldname);
775 for (label
i = 0;
i < cumEigenvalues.size();
i++)
777 os << cumEigenvalues[
i] << endl;
783std::tuple<List<Eigen::SparseMatrix<double>>, List<Eigen::VectorXd>>
785 List<Eigen::VectorXd>& b, label nmodesA, label nmodesB, word MatrixName)
788 List<Eigen::SparseMatrix<double>> ModesA(nmodesA);
789 List<Eigen::VectorXd> ModesB(nmodesB);
793 if (nmodesA >
A.size() - 2 || nmodesB >
A.size() - 2 )
796 "The number of requested modes cannot be bigger than the number of Snapshots - 2"
801 scalarField eigenValuesA(nmodesA);
802 scalarField cumEigenValuesA(nmodesA);
803 scalarField eigenValuesB(nmodesB);
804 scalarField cumEigenValuesB(nmodesB);
805 List<scalarField> eigenVectorA(nmodesA);
806 List<scalarField> eigenVectorB(nmodesB);
808 for (label
i = 0;
i < nmodesA;
i++)
810 eigenVectorA[
i].setSize(
A.size());
813 for (label
i = 0;
i < nmodesB;
i++)
815 eigenVectorB[
i].setSize(
A.size());
819 Eigen::MatrixXd corMatrixB =
corMatrix(b);
820 Info <<
"####### Performing the POD for the Matrix List #######" << endl;
821 Spectra::DenseSymMatProd<double> opA(corMatrixA);
822 Spectra::DenseSymMatProd<double> opB(corMatrixB);
823 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra:: DenseSymMatProd<double>>
824 esA(&opA, nmodesA, nmodesA + 10);
825 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra:: DenseSymMatProd<double>>
826 esB(&opB, nmodesB, nmodesB + 10);
836 Info <<
"####### End of the POD for the Matrix List #######" << endl;
837 Eigen::VectorXd eigenValueseigA;
838 Eigen::MatrixXd eigenVectorseigA;
839 Eigen::VectorXd eigenValueseigB;
840 Eigen::MatrixXd eigenVectorseigB;
842 if (esA.info() == Spectra::SUCCESSFUL)
844 eigenValueseigA = esA.eigenvalues().real();
845 eigenVectorseigA = esA.eigenvectors().real();
847 if (esB.info() == Spectra::SUCCESSFUL && nmodesB != 1)
849 eigenValueseigB = esB.eigenvalues().real();
850 eigenVectorseigB = esB.eigenvectors().real();
852 else if (nmodesB == 1)
854 eigenValueseigB.resize(1);
855 eigenVectorseigB.resize(
A.size(), nmodesB);
856 eigenValueseigB(0) = 1;
857 eigenVectorseigB = eigenVectorseigB * 0;
858 eigenVectorseigB(0, 0) = 1;
862 Info <<
"The Eigenvalue solver in ITHACAPOD.H did not converge, exiting the code"
869 Info <<
"The Eigenvalue solver in ITHACAPOD.H did not converge, exiting the code"
874 for (label
i = 0;
i < nmodesA;
i++)
876 eigenValuesA[
i] = eigenValueseigA(
i) / eigenValueseigA.sum();
879 for (label
i = 0;
i < nmodesB;
i++)
881 eigenValuesB[
i] = eigenValueseigB(
i) / eigenValueseigB.sum();
884 for (label
i = 0;
i < nmodesA;
i++)
886 for (label k = 0; k <
A.size(); k++)
888 eigenVectorA[
i][k] = eigenVectorseigA(k,
i);
892 for (label
i = 0;
i < nmodesB;
i++)
894 for (label k = 0; k <
A.size(); k++)
896 eigenVectorB[
i][k] = eigenVectorseigB(k,
i);
900 cumEigenValuesA[0] = eigenValuesA[0];
901 cumEigenValuesB[0] = eigenValuesB[0];
903 for (label
i = 1;
i < nmodesA;
i++)
905 cumEigenValuesA[
i] = cumEigenValuesA[
i - 1] + eigenValuesA[
i];
908 for (label
i = 1;
i < nmodesB;
i++)
910 cumEigenValuesB[
i] = cumEigenValuesB[
i - 1] + eigenValuesB[
i];
913 Eigen::SparseMatrix<double> tmp_A;
914 Eigen::VectorXd tmp_B;
916 for (label
i = 0;
i < nmodesA;
i++)
918 tmp_A = eigenVectorA[
i][0] *
A[0];
920 for (label k = 1; k <
A.size(); k++)
922 tmp_A += eigenVectorA[
i][k] *
A[k];
928 for (label
i = 0;
i < nmodesB;
i++)
930 tmp_B = eigenVectorB[
i][0] * b[0];
932 for (label k = 1; k <
A.size(); k++)
934 tmp_B += eigenVectorB[
i][k] * b[k];
941 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"eigenValuesA_" + MatrixName);
943 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"cumEigenValuesA_" + MatrixName);
945 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"eigenValuesB_" + MatrixName);
947 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"cumEigenValuesB_" + MatrixName);
949 for (label
i = 0;
i < ModesA.size();
i++)
952 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"A_" + MatrixName + name(
i));
955 for (label
i = 0;
i < ModesB.size();
i++)
958 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"B_" + MatrixName + name(
i));
963 for (label
i = 0;
i < nmodesA;
i++)
966 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"A_" + MatrixName + name(
i));
969 for (label
i = 0;
i < nmodesB;
i++)
972 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"B_" + MatrixName + name(
i));
976 std::tuple <List<Eigen::SparseMatrix<double>>, List<Eigen::VectorXd>> tupla;
977 tupla = std::make_tuple(ModesA, ModesB);
983 Eigen::MatrixXd Ortho = Matrix;
986 for (label
i = 0;
i < Matrix.cols();
i++)
988 for (label k = 0; k <
i; k++)
990 double num = Ortho.col(k).transpose() * Matrix.col(
i);
991 double den = (Ortho.col(k).transpose() * Ortho.col(k));
992 double fact = num / den;
993 Ortho.col(
i) -= fact * Ortho.col(k) ;
996 Ortho.col(
i).normalize();
1001template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1003 PtrList<GeometricField<Type, PatchField, GeoMesh>>& snapshots,
1004 PtrList<GeometricField<Type, PatchField, GeoMesh>>& modes,
1005 PtrList<volScalarField>& Volumes, word fieldName,
bool podex,
bool supex,
1006 bool sup, label nmodes,
bool correctBC)
1010 if (nmodes == 0 && para->
eigensolver ==
"spectra")
1012 nmodes = snapshots.size() - 2;
1017 nmodes = snapshots.size();
1022 M_Assert(nmodes <= snapshots.size() - 2,
1023 "The number of requested modes cannot be bigger than the number of Snapshots - 2");
1026 if ((podex == 0 && sup == 0) || (supex == 0 && sup == 1))
1030 label NBC = snapshots[0].boundaryField().size();
1032 Eigen::MatrixXd _corMatrix(snapshots.size(), snapshots.size());
1033 Info <<
"Filling the correlation matrix for field " << snapshots[0].name() <<
1036 for (label
i = 0;
i < snapshots.size();
i++)
1038 for (label j = 0; j <=
i; j++)
1040 Eigen::VectorXd Mij = (V.col(
i).array() * V.col(j).array());
1041 Mij = Mij.array().abs().sqrt();
1042 _corMatrix(
i, j) = SnapMatrix.col(
i).transpose() * Mij.asDiagonal() *
1047 std::cout << std::endl;
1049 for (label
i = 1;
i < snapshots.size();
i++)
1051 for (label j = 0; j <
i; j++)
1053 _corMatrix(j,
i) = _corMatrix(
i, j);
1057 Eigen::VectorXd eigenValueseig;
1058 Eigen::MatrixXd eigenVectoreig;
1059 modes.resize(nmodes);
1060 Info <<
"####### Performing the POD using EigenDecomposition for " <<
1061 snapshots[0].name() <<
" #######" << endl;
1062 label ncv = snapshots.size();
1063 Spectra::DenseSymMatProd<double> op(_corMatrix);
1064 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> esEg;
1068 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra::DenseSymMatProd<double>>
1069 es(&op, nmodes, ncv);
1070 std::cout <<
"Using Spectra EigenSolver " << std::endl;
1072 es.compute(1000, 1e-10, Spectra::LARGEST_ALGE);
1073 M_Assert(es.info() == Spectra::SUCCESSFUL,
1074 "The Eigenvalue Decomposition did not succeed");
1075 eigenVectoreig = es.eigenvectors().real();
1076 eigenValueseig = es.eigenvalues().real();
1080 std::cout <<
"Using Eigen EigenSolver " << std::endl;
1081 esEg.compute(_corMatrix);
1082 M_Assert(esEg.info() == Eigen::Success,
1083 "The Eigenvalue Decomposition did not succeed");
1084 eigenVectoreig = esEg.eigenvectors().real().rowwise().reverse().leftCols(
1086 eigenValueseig = esEg.eigenvalues().real().reverse().head(nmodes);
1089 Info <<
"####### End of the POD for " << snapshots[0].name() <<
" #######" <<
1091 Eigen::VectorXd eigenValueseigLam =
1092 eigenValueseig.real().array().cwiseInverse().sqrt() ;
1093 Eigen::MatrixXd modesEig = (SnapMatrix * eigenVectoreig) *
1094 eigenValueseigLam.asDiagonal();
1095 List<Eigen::MatrixXd> modesEigBC;
1096 modesEigBC.resize(NBC);
1098 for (label
i = 0;
i < NBC;
i++)
1100 modesEigBC[
i] = (SnapMatrixBC[
i] * eigenVectoreig) *
1101 eigenValueseigLam.asDiagonal();
1104 for (label
i = 0;
i < modes.size();
i++)
1106 GeometricField<Type, PatchField, GeoMesh> tmp2(snapshots[0].name(),
1108 Eigen::VectorXd vec = modesEig.col(
i);
1112 for (label k = 0; k < tmp2.boundaryField().size(); k++)
1117 modes.set(
i, tmp2.clone());
1120 eigenValueseig = eigenValueseig / eigenValueseig.sum();
1121 Eigen::VectorXd cumEigenValues(eigenValueseig);
1123 for (label j = 1; j < cumEigenValues.size(); ++j)
1125 cumEigenValues(j) += cumEigenValues(j - 1);
1128 Info <<
"####### Saving the POD bases for " << snapshots[0].name() <<
" #######"
1132 "./ITHACAoutput/POD/Eigenvalues_" + snapshots[0].name(), para->
precision,
1135 "./ITHACAoutput/POD/CumEigenvalues_" + snapshots[0].name(), para->
precision,
1140 Info <<
"Reading the existing modes" << endl;
1153 PtrList<volScalarField>& snapshots, PtrList<volScalarField>& modes,
1154 PtrList<volScalarField>& Volumes, word fieldName,
bool podex,
bool supex,
1155 bool sup, label nmodes,
bool correctBC);
1158 PtrList<volVectorField>& snapshots, PtrList<volVectorField>& modes,
1159 PtrList<volScalarField>& Volumes, word fieldName,
bool podex,
bool supex,
1160 bool sup, label nmodes,
bool correctBC);
1162template<
typename type_matrix>
1163std::tuple<List<Eigen::SparseMatrix<double>>, List<Eigen::VectorXd>>
1164 DEIMmodes(PtrList<type_matrix>& MatrixList, label nmodesA, label nmodesB,
1168 List<Eigen::SparseMatrix<double>> ModesA(nmodesA);
1169 List<Eigen::VectorXd> ModesB(nmodesB);
1173 M_Assert(nmodesA <= MatrixList.size() - 2
1174 && nmodesB <= MatrixList.size() - 2,
1175 "The number of requested modes cannot be bigger than the number of Snapshots - 2");
1176 std::tuple<List<Eigen::SparseMatrix<double>>, List<Eigen::VectorXd>> snapshots =
1178 Eigen::MatrixXd corMatrixA =
corMatrix(std::get<0>(snapshots));
1179 Eigen::MatrixXd corMatrixB =
corMatrix(std::get<1>(snapshots));
1180 Eigen::VectorXd eigenValueseigA;
1181 Eigen::MatrixXd eigenVectorseigA;
1182 Eigen::VectorXd eigenValueseigB;
1183 Eigen::MatrixXd eigenVectorseigB;
1187 Info <<
"####### Performing the POD decomposition for the Matrix List using Spectra #######"
1189 Spectra::DenseSymMatProd<double> opA(corMatrixA);
1190 Spectra::DenseSymMatProd<double> opB(corMatrixB);
1191 label ncvA = MatrixList.size();
1192 label ncvB = MatrixList.size();
1193 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra:: DenseSymMatProd<double>>
1194 esA(&opA, nmodesA, ncvA);
1195 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra:: DenseSymMatProd<double>>
1196 esB(&opB, nmodesB, ncvB);
1199 esA.compute(1000, 1e-10, Spectra::LARGEST_ALGE);
1200 M_Assert(esA.info() == Spectra::SUCCESSFUL,
1201 "The Eigenvalue Decomposition did not succeed");
1205 esB.compute(1000, 1e-10, Spectra::LARGEST_ALGE);
1206 M_Assert(esB.info() == Spectra::SUCCESSFUL,
1207 "The Eigenvalue Decomposition did not succeed");
1210 Info <<
"####### End of the POD decomposition for the Matrix List #######" <<
1212 eigenValueseigA = esA.eigenvalues().real();
1213 eigenVectorseigA = esA.eigenvectors().real();
1217 eigenValueseigB = esB.eigenvalues().real();
1218 eigenVectorseigB = esB.eigenvectors().real();
1222 eigenValueseigB.resize(1);
1223 eigenVectorseigB.resize(MatrixList.size(), nmodesB);
1224 eigenValueseigB(0) = 1;
1225 eigenVectorseigB = eigenVectorseigB * 0;
1226 eigenVectorseigB(0, 0) = 1;
1231 Info <<
"####### Performing the POD decomposition for the Matrix List using Eigen #######"
1233 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> esEgA;
1234 esEgA.compute(corMatrixA);
1235 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> esEgB;
1236 esEgB.compute(corMatrixB);
1237 M_Assert(esEgA.info() == Eigen::Success,
1238 "The Eigenvalue Decomposition did not succeed");
1239 M_Assert(esEgB.info() == Eigen::Success,
1240 "The Eigenvalue Decomposition did not succeed");
1241 eigenVectorseigA = esEgA.eigenvectors().real().rowwise().reverse().leftCols(
1243 eigenValueseigA = esEgA.eigenvalues().real().reverse().head(nmodesA);
1244 eigenVectorseigB = esEgB.eigenvectors().real().rowwise().reverse().leftCols(
1246 eigenValueseigB = esEgB.eigenvalues().real().reverse().head(nmodesB);
1249 Eigen::SparseMatrix<double> tmp_A;
1250 Eigen::VectorXd tmp_B;
1252 for (label
i = 0;
i < nmodesA;
i++)
1254 tmp_A = eigenVectorseigA(0,
i) * std::get<0>(snapshots)[0];
1256 for (label k = 1; k < MatrixList.size(); k++)
1258 tmp_A += eigenVectorseigA(k,
i) * std::get<0>(snapshots)[k];
1264 for (label
i = 0;
i < nmodesB;
i++)
1266 tmp_B = eigenVectorseigB(0,
i) * std::get<1>(snapshots)[0];
1268 for (label k = 1; k < MatrixList.size(); k++)
1270 tmp_B += eigenVectorseigB(k,
i) * std::get<1>(snapshots)[k];
1276 eigenValueseigA = eigenValueseigA / eigenValueseigA.sum();
1277 eigenValueseigB = eigenValueseigB / eigenValueseigB.sum();
1278 Eigen::VectorXd cumEigenValuesA(eigenValueseigA);
1279 Eigen::VectorXd cumEigenValuesB(eigenValueseigB);
1281 for (label j = 1; j < cumEigenValuesA.size(); ++j)
1283 cumEigenValuesA(j) += cumEigenValuesA(j - 1);
1286 for (label j = 1; j < cumEigenValuesB.size(); ++j)
1288 cumEigenValuesB(j) += cumEigenValuesB(j - 1);
1291 for (label
i = 0;
i < ModesA.size();
i++)
1294 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"A_" + MatrixName + name(
i));
1297 for (label
i = 0;
i < ModesB.size();
i++)
1300 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"B_" + MatrixName + name(
i));
1304 "./ITHACAoutput/DEIM/" + MatrixName +
"/eigenValuesA", para->
precision,
1307 "./ITHACAoutput/DEIM/" + MatrixName +
"/eigenValuesB", para->
precision,
1310 "./ITHACAoutput/DEIM/" + MatrixName +
"/cumEigenValuesA", para->
precision,
1313 "./ITHACAoutput/DEIM/" + MatrixName +
"/cumEigenValuesB", para->
precision,
1318 for (label
i = 0;
i < nmodesA;
i++)
1321 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"A_" + MatrixName + name(
i));
1324 for (label
i = 0;
i < nmodesB;
i++)
1327 "./ITHACAoutput/DEIM/" + MatrixName +
"/",
"B_" + MatrixName + name(
i));
1331 std::tuple <List<Eigen::SparseMatrix<double>>, List<Eigen::VectorXd>> tupla;
1332 tupla = std::make_tuple(ModesA, ModesB);
1336template std::tuple<List<Eigen::SparseMatrix<double>>, List<Eigen::VectorXd>>
1337DEIMmodes(PtrList<fvScalarMatrix>& MatrixList, label nmodesA,
1341template std::tuple<List<Eigen::SparseMatrix<double>>, List<Eigen::VectorXd>>
1342DEIMmodes(PtrList<fvVectorMatrix>& MatrixList, label nmodesA,
1346template<
class Type,
template<
class>
class PatchField,
class GeoMesh>
1348 PtrList<GeometricField<Type, PatchField, GeoMesh>>& snapshots, label nmodes,
1349 word FunctionName, word fieldName)
1352 word PODkey =
"POD_" + fieldName;
1353 word PODnorm = para->
ITHACAdict->lookupOrDefault<word>(PODkey,
"L2");
1355 PODnorm ==
"Frobenius",
"The PODnorm can be only L2 or Frobenius");
1356 Info <<
"Performing POD for " << fieldName <<
" using the " << PODnorm <<
1358 PtrList<GeometricField<Type, fvPatchField, volMesh>> modes;
1359 bool correctBC =
true;
1361 if (nmodes == 0 && para->
eigensolver ==
"spectra")
1363 nmodes = snapshots.size() - 2;
1368 nmodes = snapshots.size();
1373 M_Assert(nmodes <= snapshots.size() - 2,
1374 "The number of requested modes cannot be bigger than the number of Snapshots - 2");
1381 label NBC = snapshots[0].boundaryField().size();
1382 Eigen::MatrixXd _corMatrix;
1384 if (PODnorm ==
"L2")
1388 else if (PODnorm ==
"Frobenius")
1393 if (Pstream::parRun())
1395 reduce(_corMatrix, sumOp<Eigen::MatrixXd>());
1398 Eigen::VectorXd eigenValueseig;
1399 Eigen::MatrixXd eigenVectoreig;
1400 modes.resize(nmodes);
1401 Info <<
"####### Performing the POD using EigenDecomposition " <<
1402 fieldName <<
" #######" << endl;
1403 label ncv = snapshots.size();
1404 Spectra::DenseSymMatProd<double> op(_corMatrix);
1405 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> esEg;
1409 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra::DenseSymMatProd<double>>
1410 es(&op, nmodes, ncv);
1411 std::cout <<
"Using Spectra EigenSolver " << std::endl;
1413 es.compute(1000, 1e-10, Spectra::LARGEST_ALGE);
1414 M_Assert(es.info() == Spectra::SUCCESSFUL,
1415 "The Eigenvalue Decomposition did not succeed");
1416 eigenVectoreig = es.eigenvectors().real();
1417 eigenValueseig = es.eigenvalues().real();
1421 std::cout <<
"Using Eigen EigenSolver " << std::endl;
1422 esEg.compute(_corMatrix);
1423 M_Assert(esEg.info() == Eigen::Success,
1424 "The Eigenvalue Decomposition did not succeed");
1425 eigenVectoreig = esEg.eigenvectors().real().rowwise().reverse().leftCols(
1427 eigenValueseig = esEg.eigenvalues().real().array().reverse();
1430 if (eigenValueseig.array().minCoeff() < 0)
1432 eigenValueseig = eigenValueseig.array() + 2 * abs(
1433 eigenValueseig.array().minCoeff());
1436 Info <<
"####### End of the POD for " << snapshots[0].name() <<
" #######" <<
1438 Eigen::MatrixXd modesEig = (SnapMatrix * eigenVectoreig);
1440 Eigen::MatrixXd normFact(nmodes, 1);
1442 for (label
i = 0;
i < nmodes;
i++)
1444 if (PODnorm ==
"L2")
1446 normFact(
i, 0) = std::sqrt((modesEig.col(
i).transpose() * V.asDiagonal() *
1447 modesEig.col(
i))(0, 0));
1449 if (Pstream::parRun())
1451 normFact(
i, 0) = (modesEig.col(
i).transpose() * V.asDiagonal() *
1452 modesEig.col(
i))(0, 0);
1455 else if (PODnorm ==
"Frobenius")
1457 normFact(
i, 0) = std::sqrt((modesEig.col(
i).transpose() * modesEig.col(
i))(0,
1460 if (Pstream::parRun())
1462 normFact(
i, 0) = (modesEig.col(
i).transpose() * modesEig.col(
i))(0, 0);
1467 if (Pstream::parRun())
1469 reduce(normFact, sumOp<Eigen::MatrixXd>());
1472 if (Pstream::parRun())
1474 normFact = normFact.cwiseSqrt();
1477 List<Eigen::MatrixXd> modesEigBC;
1478 modesEigBC.resize(NBC);
1480 for (label
i = 0;
i < NBC;
i++)
1482 modesEigBC[
i] = (SnapMatrixBC[
i] * eigenVectoreig);
1485 std::cout << normFact << std::endl;
1487 for (label
i = 0;
i < nmodes;
i++)
1489 modesEig.col(
i) = modesEig.col(
i).array() / normFact(
i, 0);
1491 for (label j = 0; j < NBC; j++)
1493 modesEigBC[j].col(
i) = modesEigBC[j].col(
i).array() / normFact(
i, 0);
1497 for (label
i = 0;
i < modes.size();
i++)
1499 GeometricField<Type, PatchField, GeoMesh> tmp2(snapshots[0].name(),
1501 Eigen::VectorXd vec = modesEig.col(
i);
1504 for (label k = 0; k < NBC; k++)
1509 modes.set(
i, tmp2.clone());
1512 eigenValueseig = eigenValueseig / eigenValueseig.sum();
1513 Eigen::VectorXd cumEigenValues(eigenValueseig);
1515 for (label j = 1; j < cumEigenValues.size(); ++j)
1517 cumEigenValues(j) += cumEigenValues(j - 1);
1520 Info <<
"####### Saving the POD bases for " << snapshots[0].name() <<
1524 for (label
i = 0;
i < modes.size();
i++)
1531 "./ITHACAoutput/DEIM/eigenValues_" + fieldName, para->
precision,
1534 "./ITHACAoutput/DEIM/cumEigenValues_" + fieldName, para->
precision,
1539 Info <<
"Reading the existing modes" << endl;
1546template<
class Field_type,
class Field_type_2>
1548 PtrList<Field_type>& snapshots, PtrList<Field_type>& modes,
1549 PtrList<Field_type_2>& fields2, word fieldName,
bool podex,
bool supex,
1550 bool sup, label nmodes,
bool correctBC)
1554 if ((podex == 0 && sup == 0) || (supex == 0 && sup == 1))
1560 nmodes = snapshots.size() - 2;
1563 M_Assert(nmodes <= snapshots.size() - 2,
1564 "The number of requested modes cannot be bigger than the number of Snapshots - 2");
1570 nmodes = snapshots.size();
1573 M_Assert(nmodes <= snapshots.size(),
1574 "The number of requested modes cannot be bigger than the number of Snapshots");
1579 label NBC = fields2[0].boundaryField().size();
1581 Eigen::MatrixXd _corMatrix = SnapMatrix.transpose() * VM.asDiagonal() *
1584 if (Pstream::parRun())
1586 reduce(_corMatrix, sumOp<Eigen::MatrixXd>());
1589 Eigen::VectorXd eigenValueseig;
1590 Eigen::MatrixXd eigenVectoreig;
1591 modes.resize(nmodes);
1592 Info <<
"####### Performing the POD using EigenDecomposition " <<
1593 fields2[0].name() <<
" #######" << endl;
1594 label ncv = snapshots.size();
1595 Spectra::DenseSymMatProd<double> op(_corMatrix);
1596 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> esEg;
1600 Spectra::SymEigsSolver<double, Spectra::LARGEST_ALGE, Spectra::DenseSymMatProd<double>>
1601 es(&op, nmodes, ncv);
1602 std::cout <<
"Using Spectra EigenSolver " << std::endl;
1604 es.compute(1000, 1e-10, Spectra::LARGEST_ALGE);
1605 M_Assert(es.info() == Spectra::SUCCESSFUL,
1606 "The Eigenvalue Decomposition did not succeed");
1607 eigenVectoreig = es.eigenvectors().real();
1608 eigenValueseig = es.eigenvalues().real();
1612 std::cout <<
"Using Eigen EigenSolver " << std::endl;
1613 esEg.compute(_corMatrix);
1614 M_Assert(esEg.info() == Eigen::Success,
1615 "The Eigenvalue Decomposition did not succeed");
1616 eigenVectoreig = esEg.eigenvectors().real().rowwise().reverse().leftCols(
1618 eigenValueseig = esEg.eigenvalues().real().reverse().head(nmodes);
1621 Info <<
"####### End of the POD for " << snapshots[0].name() <<
" #######" <<
1623 Eigen::VectorXd eigenValueseigLam =
1624 eigenValueseig.real().array().cwiseInverse().abs().sqrt() ;
1625 Eigen::MatrixXd modesEig = (SnapMatrix * eigenVectoreig) *
1626 eigenValueseigLam.asDiagonal();
1627 List<Eigen::MatrixXd> modesEigBC;
1628 modesEigBC.resize(NBC);
1630 for (label
i = 0;
i < modes.size();
i++)
1632 Field_type tmp2(snapshots[0].name(), snapshots[0] * 0);
1634 for (label j = 0; j < snapshots.size(); j++)
1636 tmp2 += eigenValueseigLam(
i) * eigenVectoreig.col(
i)(j) * snapshots[j];
1639 modes.set(
i, tmp2.clone());
1642 eigenValueseig = eigenValueseig / eigenValueseig.sum();
1643 Eigen::VectorXd cumEigenValues(eigenValueseig);
1645 for (label j = 1; j < cumEigenValues.size(); ++j)
1647 cumEigenValues(j) += cumEigenValues(j - 1);
1650 Info <<
"####### Saving the POD bases for " << snapshots[0].name() <<
1656 snapshots[0].name());
1664 "./ITHACAoutput/POD/Eigenvalues_" + snapshots[0].name(), para->
precision,
1667 "./ITHACAoutput/POD/CumEigenvalues_" + snapshots[0].name(), para->
precision,
1672 Info <<
"Reading the existing modes" << endl;
1686 PtrList<surfaceScalarField>& snapshots, PtrList<surfaceScalarField>& modes,
1687 PtrList<volVectorField>& fields2, word fieldName,
bool podex,
bool supex,
1688 bool sup, label nmodes,
bool correctBC);
1691 PtrList<volScalarField>& snapshots, PtrList<volScalarField>& modes,
1692 PtrList<volVectorField>& fields2, word fieldName,
bool podex,
bool supex,
1693 bool sup, label nmodes,
bool correctBC);
1695template PtrList<volScalarField>
1697 PtrList<volScalarField>& SnapShotsMatrix,
1699 word FunctionName, word FieldName);
1701template PtrList<volVectorField>
1703 PtrList<volVectorField>& SnapShotsMatrix,
1705 word FunctionName, word FieldName);
Header file of the EigenFunctions class.
Header file of the ITHACAPOD class.
#define M_Assert(Expr, Msg)
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 std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > LFvMatrix2LSM(PtrList< fvMatrix< Type > > &MatrixList)
Convert a PtrList of OpenFOAM fvMatrix into a tuple of lists of Eigen Sparse Matrices and source vect...
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)
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
word eigensolver
type of eigensolver used in the eigenvalue decomposition can be either be eigen or spectra
label precision
precision of the output Market Matrix objects (i.e. reduced matrices, eigenvalues,...
IOdictionary * ITHACAdict
Dictionary for input objects from file.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
std::_Ios_Fmtflags outytpe
type of output format can be fixed or scientific
bool saveMarketVector(const VectorType &vec, const std::string &filename, label prec, std::_Ios_Fmtflags outytpe=std::ios_base::scientific)
namespace for the computation of the POD, it exploits bot the method of snapshots and SVD.
void exportcumEigenvalues(scalarField cumEigenvalues, fileName name, bool sup)
Export the cumulative eigenvalues as a txt file.
void exportEigenvalues(scalarField Eigenvalues, fileName name, bool sup)
Exports the eigenvalues as a txt file.
void getModes(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Computes the bases or reads them for a field.
std::tuple< List< Eigen::SparseMatrix< double > >, List< Eigen::VectorXd > > DEIMmodes(List< Eigen::SparseMatrix< double > > &A, List< Eigen::VectorXd > &b, label nmodesA, label nmodesB, word MatrixName)
Get the DEIM modes for a generic non a parametrized matrix coming from a differential operator functi...
void getNestedSnapshotMatrix(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &ModesGlobal, word fieldName, label Npar, label NnestedOut)
Nested-POD approach.
void exportBases(PtrList< GeometricField< Type, PatchField, GeoMesh > > &s, PtrList< GeometricField< Type, PatchField, GeoMesh > > &bases, word fieldName, bool sup)
Export the Bases.
void getModesSVD(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Gets the bases for a scalar field using SVD instead of the method of snapshots.
Eigen::MatrixXd corMatrix(PtrList< volScalarField > &snapshots)
Construct the Correlation Matrix for Scalar Field.
void GrammSchmidt(Eigen::MatrixXd &Matrix)
Performs GrammSchmidt orthonormalization on an Eigen Matrix.
void getWeightedModes(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Computes the weighted bases (using the nested-pod approach) or read them for a vector field.
void SaveSparseMatrix(Eigen::SparseMatrix< T, Nrows, IND > &m, word folder, word MatrixName)
Export an Eigen sparse matrix into bynary format file.
void ReadDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Read a dense matrix from a binary format file.
void exportFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &field, word folder, word fieldname)
Function to export a scalar of vector field.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void exportList(T &list, word folder, word filename)
Export a list to file.
void SaveDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Save a dense matrix to a binary format file.
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
void ReadSparseMatrix(Eigen::SparseMatrix< T, Nrows, IND > &m, word folder, word MatrixName)
Read an Eigen sparse matrix from a bynary format file.
Eigen::VectorXd getMassMatrixFV(GeometricField< Type, PatchField, GeoMesh > &snapshot)
Gets a vector containing the volumes of each cell of the mesh.
Eigen::MatrixXd getMassMatrix(PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, label Nmodes, bool consider_volumes)
Gets the mass matrix using the eigen routine.
void assignBC(GeometricField< scalar, fvPatchField, volMesh > &s, label BC_ind, double value)
Assign uniform Boundary Condition to a volScalarField.
bool check_folder(word folder)
Checks if a folder exists.
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
Info<< endl;Info<< "*********************************************************"<< endl;Info<< "Performing test for the parameterized BC inverse solver"<< endl;Info<< endl;word outputFolder="./ITHACAoutput/parameterizedBCtest_RBFparameter/";volScalarField gTrueField=example_paramBC.list2Field(example_paramBC.gTrue);ITHACAstream::exportSolution(gTrueField, "1", outputFolder, "gTrue");Eigen::VectorXd rbfWidth=EigenFunctions::ExpSpaced(0.001, 100, rbfWidthTest_size);ITHACAstream::exportMatrix(rbfWidth, "rbfShapeParameters", "eigen", outputFolder);Eigen::VectorXd residualNorms;residualNorms.resize(rbfWidthTest_size);Eigen::VectorXd heatFluxL2norm(rbfWidthTest_size);Eigen::VectorXd heatFluxLinfNorm=heatFluxL2norm;Eigen::VectorXd condNumber=heatFluxL2norm;Eigen::MatrixXd singVal;for(int i=0;i< rbfWidthTest_size;i++){ Info<< "*********************************************************"<< endl;Info<< "RBF parameter "<< rbfWidth(i)<< endl;Info<< "Test "<< i+1<< endl;Info<< endl;example_paramBC.set_gParametrized("rbf", rbfWidth(i));example_paramBC.parameterizedBCoffline(1);example_paramBC.parameterizedBC("fullPivLU");volScalarField gParametrizedField=example_paramBC.list2Field(example_paramBC.g);ITHACAstream::exportSolution(gParametrizedField, std::to_string(i+1), outputFolder, "gParametrized");volScalarField &T(example_paramBC._T());ITHACAstream::exportSolution(T, std::to_string(i+1), outputFolder, "T");residualNorms(i)=Foam::sqrt(example_paramBC.residual.squaredNorm());Eigen::MatrixXd A=example_paramBC.Theta.transpose() *example_paramBC.Theta;Eigen::JacobiSVD< Eigen::MatrixXd > svd(A, Eigen::ComputeThinU|Eigen::ComputeThinV)