36#include "viscosityModel.H"
44 _args = autoPtr<argList>
46 new argList(argc, argv)
49 if (!
_args->checkRootCase())
51 Foam::FatalError.exit();
54 argList& args =
_args();
55#include "createTime.H"
56#include "createMesh.H"
57 _simple = autoPtr<simpleControl>
64 simpleControl& simple =
_simple();
65#include "createFields.H"
66#include "createFvOptions.H"
83 "The BC method must be set to lift or penalty or none in ITHACAdict");
86 "The flux method must be set to inconsistent or consistent in ITHACAdict");
100 fvMesh& mesh =
_mesh();
101 volScalarField& p =
_p();
102 volVectorField& U =
_U();
103 surfaceScalarField& phi =
_phi();
105 simpleControl& simple =
_simple();
106 IOMRFZoneList& MRF =
_MRF();
108#include "NLsolvesteadyNS.H"
118 for (label i = 0; i < mu_now.size(); i++)
132 "./ITHACAoutput/Offline");
139 M_Assert(type ==
"modes"
140 || type ==
"snapshots",
141 "You must specify the variable type with either snapshots or modes");
142 PtrList<volScalarField> P_sup;
144 if (type ==
"snapshots")
155 volVectorField U =
_U();
167 dimensionedVector(
"zero", U.dimensions(), vector::zero)
170 if (type ==
"snapshots")
181 volVectorField U =
_U();
193 dimensionedVector(
"zero", U.dimensions(), vector::zero)
195 dimensionedScalar nu_fake
198 dimensionSet(0, 2, -1, 0, 0, 0, 0),
201 Vector<double> v(0, 0, 0);
203 for (label i = 0; i < Usup.boundaryField().size(); i++)
205 if (Usup.boundaryField()[i].type() !=
"processor")
213 if (type ==
"snapshots")
215 for (label i = 0; i < P_sup.size(); i++)
217 fvVectorMatrix u_sup_eqn
219 - fvm::laplacian(nu_fake, Usup)
223 u_sup_eqn == fvc::grad(P_sup[i])
233 for (label i = 0; i < P_sup.size(); i++)
235 fvVectorMatrix u_sup_eqn
237 - fvm::laplacian(nu_fake, Usup)
241 u_sup_eqn == fvc::grad(P_sup[i])
258 surfaceScalarField& phi =
_phi();
259 fvMesh& mesh =
_mesh();
260 volScalarField p =
_p();
261 volVectorField U =
_U();
262 IOMRFZoneList& MRF =
_MRF();
264 volVectorField Ulift(
"Ulift" + name(k), U);
265 instantList Times = runTime.times();
266 runTime.setTime(Times[1], 1);
267 pisoControl potentialFlow(mesh,
"potentialFlow");
268 Info <<
"Solving a lifting Problem" << endl;
269 Vector<double> v1(0, 0, 0);
271 Vector<double> v0(0, 0, 0);
273 for (label j = 0; j < U.boundaryField().size(); j++)
279 else if (U.boundaryField()[BCind].type() ==
"fixedValue")
288 phi = linearInterpolate(Ulift) & mesh.Sf();
291 Info <<
"Constructing velocity potential field Phi\n" << endl;
299 IOobject::READ_IF_PRESENT,
303 dimensionedScalar(
"Phi", dimLength * dimVelocity, 0),
304 p.boundaryField().types()
306 label PhiRefCell = 0;
307 scalar PhiRefValue = 0;
311 potentialFlow.dict(),
315 mesh.setFluxRequired(Phi.name());
316 runTime.functionObjects().start();
317 MRF.makeRelative(phi);
318 adjustPhi(phi, Ulift, p);
320 while (potentialFlow.correctNonOrthogonal())
322 fvScalarMatrix PhiEqn
324 fvm::laplacian(dimensionedScalar(
"1", dimless, 1), Phi)
328 PhiEqn.setReference(PhiRefCell, PhiRefValue);
331 if (potentialFlow.finalNonOrthogonalIter())
333 phi -= PhiEqn.flux();
337 MRF.makeAbsolute(phi);
338 Info <<
"Continuity error = "
339 << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
341 Ulift = fvc::reconstruct(phi);
342 Ulift.correctBoundaryConditions();
343 Info <<
"Interpolated velocity error = "
344 << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
345 / sum(mesh.magSf())).value()
363 for (label k = 0; k <
liftfield.size(); k++)
371 for (label k = 0; k <
NUmodes; k++)
379 word B_str =
"B_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
391 word K_str =
"K_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
403 word M_str =
"M_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
415 word D_str =
"D_" + name(
NPmodes);
426 word bc3_str =
"BC3_" + name(
liftfield.size()) +
"_" + name(
438 word bc4_str =
"BC4_" + name(
liftfield.size()) +
"_" + name(
450 word C_str =
"C_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
462 word G_str =
"G_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
499 if (para->exportPython)
506 "./ITHACAoutput/Matrices/");
508 "./ITHACAoutput/Matrices/");
513 if (para->exportMatlab)
520 "./ITHACAoutput/Matrices/");
522 "./ITHACAoutput/Matrices/");
534 "./ITHACAoutput/Matrices/");
536 "./ITHACAoutput/Matrices/");
538 "./ITHACAoutput/Matrices/C");
540 "./ITHACAoutput/Matrices/G");
545 cnpy::save(
B_matrix,
"./ITHACAoutput/Matrices/B.npy");
546 cnpy::save(
K_matrix,
"./ITHACAoutput/Matrices/K.npy");
547 cnpy::save(
D_matrix,
"./ITHACAoutput/Matrices/D.npy");
548 cnpy::save(
M_matrix,
"./ITHACAoutput/Matrices/M.npy");
549 cnpy::save(
BC3_matrix,
"./ITHACAoutput/Matrices/BC3.npy");
550 cnpy::save(
BC4_matrix,
"./ITHACAoutput/Matrices/BC4.npy");
551 cnpy::save(
C_tensor,
"./ITHACAoutput/Matrices/C.npy");
552 cnpy::save(
gTensor,
"./ITHACAoutput/Matrices/G.npy");
565 for (label k = 0; k <
liftfield.size(); k++)
573 for (label k = 0; k <
NUmodes; k++)
589 word B_str =
"B_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
601 word K_str =
"K_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
613 word P_str =
"P_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
625 word M_str =
"M_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
637 word C_str =
"C_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
671 if (para->exportPython)
680 if (para->exportMatlab)
696 "./ITHACAoutput/Matrices/C");
701 cnpy::save(
B_matrix,
"./ITHACAoutput/Matrices/B.npy");
702 cnpy::save(
K_matrix,
"./ITHACAoutput/Matrices/K.npy");
703 cnpy::save(
P_matrix,
"./ITHACAoutput/Matrices/P.npy");
704 cnpy::save(
M_matrix,
"./ITHACAoutput/Matrices/M.npy");
705 cnpy::save(
C_tensor,
"./ITHACAoutput/Matrices/C.npy");
716 Vector<double> inl(0, 0, 0);
720 for (label k = 0; k <
NUmodes; k++)
735 for (label k = 0; k <
NUmodes; k++)
738 forAll(
L_PHImodes[k].boundaryFieldRef()[0], l)
747 Uinl = autoPtr<volVectorField> (
new volVectorField(
_U()));
750 dt_dummy = autoPtr<dimensionedScalar> (
new dimensionedScalar
753 dimensionSet(0, 0, 1, 0, 0, 0, 0),
756 nu_dummy = autoPtr<dimensionedScalar> (
new dimensionedScalar
759 dimensionSet(0, 2, -1, 0, 0, 0, 0),
765 word B_str =
"B_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
777 word K_str =
"K_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
789 word P_str =
"P_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
801 word C_str =
"C_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
813 word Cf_str =
"Cf_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
825 word BP_str =
"BP_" + name(
liftfield.size()) +
"_" + name(
NUmodes) +
"_" + name(
891 for (label i = 0; i < Bsize; i++)
893 for (label j = 0; j < Bsize; j++)
896 dimensionedScalar(
"1", dimless, 1),
L_U_SUPmodes[j])).value();
900 if (Pstream::master())
917 for (label i = 0; i < Bsize; i++)
919 for (label j = 0; j < Bsize; j++)
926 if (Pstream::master())
940 Eigen::MatrixXd
K_matrix(K1size, K2size);
943 for (label i = 0; i < K1size; i++)
945 for (label j = 0; j < K2size; j++)
952 if (Pstream::master())
969 for (label j = 0; j < Csize; j++)
974 for (label i = 0; i < Csize; i++)
976 for (label j = 0; j < Csize; j++)
978 for (label k = 0; k < Csize; k++)
987 if (Pstream::master())
1003 C_tensor.resize(Csize, Csize, Csize);
1005 for (label i = 0; i < Csize; i++)
1007 for (label j = 0; j < Csize; j++)
1009 for (label k = 0; k < Csize; k++)
1027 if (Pstream::master())
1043 Eigen::Tensor<double, 3>
C_tensor(Csize, Csize, Csize);
1045 for (label j = 0; j < Csize; ++j)
1047 surfaceScalarField SfUj = linearInterpolate(
L_U_SUPmodes[j]) &
1050 List<tmp<volVectorField>> divRow(Csize);
1052 for (label k = 0; k < Csize; ++k)
1064 for (label i = 0; i < Csize; ++i)
1066 for (label k = 0; k < Csize; ++k)
1068 const volVectorField& divField = divRow[k]();
1074 if (Pstream::master())
1089 Eigen::MatrixXd
M_matrix(Msize, Msize);
1092 for (label i = 0; i < Msize; i++)
1094 for (label j = 0; j < Msize; j++)
1101 if (Pstream::master())
1117 Eigen::MatrixXd
P_matrix(P1size, P2size);
1120 for (label i = 0; i < P1size; i++)
1122 for (label j = 0; j < P2size; j++)
1129 if (Pstream::master())
1147 for (label j = 0; j < G1size; j++)
1149 G_matrix[j].resize(G2size, G2size);
1152 for (label i = 0; i < G1size; i++)
1154 for (label j = 0; j < G2size; j++)
1156 for (label k = 0; k < G2size; k++)
1158 G_matrix[i](j, k) = fvc::domainIntegrate(fvc::grad(
Pmodes[i]) & (fvc::div(
1165 if (Pstream::master())
1180 Eigen::Tensor<double, 3>
gTensor;
1181 gTensor.resize(g1Size, g2Size, g2Size);
1183 for (label i = 0; i < g1Size; i++)
1185 for (label j = 0; j < g2Size; j++)
1187 for (label k = 0; k < g2Size; k++)
1189 gTensor(i, j, k) = fvc::domainIntegrate(fvc::grad(
Pmodes[i]) & (fvc::div(
1196 if (Pstream::master())
1212 Eigen::Tensor<double, 3>
gTensor (g1Size, g2Size, g2Size);
1213 List<tmp<volVectorField>> PmodesGrad(g1Size);
1215 for (label i = 0; i < g1Size; i++)
1217 PmodesGrad[i] = fvc::grad(
Pmodes[i]);
1220 for (label j = 0; j < g2Size; j++)
1222 surfaceScalarField interpUj = fvc::interpolate(
L_U_SUPmodes[j]) &
1225 List<tmp<volVectorField>> divRow(g2Size);
1227 for (label k = 0; k < g2Size; ++k)
1232 for (label i = 0; i < g1Size; i++)
1234 const volVectorField& gradPi = PmodesGrad[i]();
1236 for (label k = 0; k < g2Size; k++)
1238 const volVectorField& divField = divRow[k]();
1239 gTensor(i, j, k) = fvc::domainIntegrate(gradPi & divField).value();
1244 if (Pstream::master())
1260 Eigen::MatrixXd
L_matrix(Lsize, Lsize);
1262 for (label i = 0; i < Lsize; i++)
1264 for (label j = 0; j < Lsize; j++)
1267 fvc::interpolate(vls) & vls.mesh().Sf(),
1276 label
NUmodes, volVectorField vls)
1280 Eigen::MatrixXd
L_D_matrix(LDsize1, LDsize2);
1282 for (label i = 0; i < LDsize1; i++)
1284 for (label j = 0; j < LDsize2; j++)
1287 fvc::interpolate(vls) & vls.mesh().Sf(),
1298 Eigen::MatrixXd
D_matrix(Dsize, Dsize);
1301 for (label i = 0; i < Dsize; i++)
1303 for (label j = 0; j < Dsize; j++)
1305 D_matrix(i, j) = fvc::domainIntegrate(fvc::grad(
Pmodes[i]) & fvc::grad(
1310 if (Pstream::master())
1323 Eigen::MatrixXd
BC1_matrix(P_BC1size, P_BC2size);
1326 for (label i = 0; i < P_BC1size; i++)
1328 for (label j = 0; j < P_BC2size; j++)
1330 surfaceScalarField lpl((fvc::interpolate(fvc::laplacian(
1334 for (label k = 0; k < lpl.boundaryField().size(); k++)
1336 if (lpl.boundaryField()[k].type() !=
"processor")
1338 s += gSum(lpl.boundaryField()[k]);
1346 if (Pstream::master())
1364 for (label j = 0; j < P2_BC1size; j++)
1366 BC2_matrix[j].resize(P2_BC2size, P2_BC2size);
1369 for (label i = 0; i < P2_BC1size; i++)
1371 for (label j = 0; j < P2_BC2size; j++)
1373 for (label k = 0; k < P2_BC2size; k++)
1375 surfaceScalarField div_m(fvc::interpolate(fvc::div(fvc::interpolate(
1380 for (label k = 0; k < div_m.boundaryField().size(); k++)
1382 if (div_m.boundaryField()[k].type() !=
"processor")
1384 s += gSum(div_m.boundaryField()[k]);
1398 label pressureBC1Size =
NPmodes;
1402 bc2Tensor.resize(pressureBC1Size, pressureBC2Size, pressureBC2Size);
1404 for (label i = 0; i < pressureBC1Size; i++)
1406 for (label j = 0; j < pressureBC2Size; j++)
1408 for (label k = 0; k < pressureBC2Size; k++)
1410 surfaceScalarField div_m(fvc::interpolate(fvc::div(fvc::interpolate(
1415 for (label k = 0; k < div_m.boundaryField().size(); k++)
1417 if (div_m.boundaryField()[k].type() !=
"processor")
1419 s += gSum(div_m.boundaryField()[k]);
1428 if (Pstream::master())
1443 Eigen::MatrixXd
BC3_matrix(P3_BC1size, P3_BC2size);
1445 surfaceVectorField n(mesh.Sf() / mesh.magSf());
1447 for (label i = 0; i < P3_BC1size; i++)
1449 for (label j = 0; j < P3_BC2size; j++)
1451 surfaceVectorField BC3 = fvc::interpolate(fvc::curl(
L_U_SUPmodes[j])).ref();
1452 surfaceVectorField BC4 = (n ^ fvc::interpolate(fvc::grad(
Pmodes[i]))).ref();
1453 surfaceScalarField BC5 = ((BC3 & BC4) * mesh.magSf()).ref();
1456 for (label k = 0; k < BC5.boundaryField().size(); k++)
1458 if (BC5.boundaryField()[k].type() !=
"processor")
1460 s += gSum(BC5.boundaryField()[k]);
1468 if (Pstream::master())
1481 Eigen::MatrixXd
BC4_matrix(P4_BC1size, P4_BC2size);
1483 surfaceVectorField n(mesh.Sf() / mesh.magSf());
1485 for (label i = 0; i < P4_BC1size; i++)
1487 for (label j = 0; j < P4_BC2size; j++)
1489 surfaceScalarField BC3 = fvc::interpolate(
Pmodes[i]).ref();
1490 surfaceScalarField BC4 = (n & fvc::interpolate(
L_U_SUPmodes[j])).ref();
1491 surfaceScalarField BC5 = ((BC3 * BC4) * mesh.magSf()).ref();
1494 for (label k = 0; k < BC5.boundaryField().size(); k++)
1496 if (BC5.boundaryField()[k].type() !=
"processor")
1498 s += gSum(BC5.boundaryField()[k]);
1506 if (Pstream::master())
1521 for (label j = 0; j <
inletIndex.rows(); j++)
1526 for (label k = 0; k <
inletIndex.rows(); k++)
1531 for (label i = 0; i < BCsize; i++)
1538 if (Pstream::master())
1541 "./ITHACAoutput/Matrices/bcVelVec");
1552 List <Eigen::MatrixXd>
bcVelMat(BCUsize);
1554 for (label j = 0; j <
inletIndex.rows(); j++)
1556 bcVelMat[j].resize(BCsize, BCsize);
1559 for (label k = 0; k <
inletIndex.rows(); k++)
1564 for (label i = 0; i < BCsize; i++)
1566 for (label j = 0; j < BCsize; j++)
1570 L_U_SUPmodes[j].boundaryField()[BCind].component(BCcomp));
1575 if (Pstream::master())
1578 "./ITHACAoutput/Matrices/bcVelMat");
1589 Eigen::MatrixXd
BP_matrix(BPsize1, BPsize2);
1592 for (label i = 0; i < BPsize1; i++)
1594 for (label j = 0; j < BPsize2; j++)
1596 L_U_SUPmodesaux =
dt_dummy * fvc::laplacian(
1599 fvc::div(L_U_SUPmodesaux)).value();
1603 if (Pstream::master())
1615 Eigen::VectorXd ModeVector;
1618 List <Eigen::MatrixXd>
RD_matrix(BCsize);
1619 Eigen::SparseMatrix<double> A;
1622 for (label i = 0; i < BCsize; i++)
1625 Vector<double> v(0, 0, 0);
1627 volVectorField Upara(
Uinl());
1637 for (label j = 0; j < RDsize; j++)
1640 RD_matrix[i](j, 0) = ModeVector.dot(b.col(0));
1643 if (Pstream::parRun())
1645 reduce(
RD_matrix[i], sumOp<Eigen::MatrixXd>());
1648 if (Pstream::master())
1661 Eigen::VectorXd ModeVector;
1664 List <Eigen::MatrixXd>
RC_matrix(BCsize);
1665 Eigen::SparseMatrix<double> A;
1668 for (label i = 0; i < BCsize; i++)
1671 Vector<double> v(0, 0, 0);
1673 volVectorField Upara(
Uinl());
1678 fvm::div(fvc::flux(Upara), Upara)
1683 for (label j = 0; j < RCsize; j++)
1686 RC_matrix[i](j, 0) = ModeVector.dot(b.col(0));
1689 if (Pstream::parRun())
1691 reduce(
RC_matrix[i], sumOp<Eigen::MatrixXd>());
1694 if (Pstream::master())
1711 Cf_tensor.resize(Csize2, Csize1, Csize1);
1713 for (label i = 0; i < Csize2; i++)
1715 for (label j = 0; j < Csize1; j++)
1717 for (label k = 0; k < Csize1; k++)
1721 L_U_SUPmodesaux =
dt_dummy * (fvc::div(
1727 L_U_SUPmodesaux =
dt_dummy * (fvc::div(
1733 * fvc::div(L_U_SUPmodesaux)).value();
1738 if (Pstream::master())
1751 List<Eigen::MatrixXd> LinSysDivDummy;
1753 volScalarField p(
_p());
1755 for (label i = 0; i < BCsize; i++)
1758 Vector<double> v(0, 0, 0);
1760 volVectorField Upara(
Uinl());
1764 fvm::laplacian(p) == (1.0 /
dt_dummy) * fvc::div(Upara)
1766 pEqn.setReference(0, 0);
1784 List<Eigen::MatrixXd> LinSysConvDummy;
1786 volScalarField p(
_p());
1788 for (label i = 0; i < BCsize; i++)
1791 Vector<double> v(0, 0, 0);
1793 volVectorField Upara(
Uinl());
1796 Caux =
dt_dummy * (-fvc::div(fvc::flux(Upara), Upara));
1799 fvm::laplacian(p) == (1.0 /
dt_dummy) * fvc::div(Caux)
1801 pEqn.setReference(0, 0);
1819 List<Eigen::MatrixXd> LinSysDiffDummy;
1821 volScalarField p(
_p());
1823 for (label i = 0; i < BCsize; i++)
1826 Vector<double> v(0, 0, 0);
1828 volVectorField Upara(
Uinl());
1834 fvm::laplacian(p) == (1.0 /
dt_dummy) * fvc::div(Daux)
1836 pEqn.setReference(0, 0);
1858 for (label i = 0; i <
NUmodes; i++)
1860 volVectorField CoeffA = (fvc::reconstruct(
L_PHImodes[i])).ref();
1862 for (label j = 0; j < Isize; j++)
1864 surfaceScalarField B = fvc::flux(
L_U_SUPmodes[j]).ref();
1865 volVectorField CoeffB = fvc::reconstruct(B).ref();
1866 I_matrix(i, j) = fvc::domainIntegrate(CoeffA & CoeffB).value();
1870 if (Pstream::master())
1886 surfaceScalarField phi_tmp(
"Phi_tmp",
_phi());
1888 for (label i = 0; i <
NUmodes; i++)
1890 for (label j = 0; j < DFsize; j++)
1893 dimensionedScalar(
"1", dimless, 1),
1895 volVectorField CoeffB = fvc::reconstruct(phi_tmp).ref();
1896 volVectorField CoeffA = fvc::reconstruct(
L_PHImodes[i]).ref();
1897 DF_matrix(i, j) = fvc::domainIntegrate(CoeffA & CoeffB).value();
1901 if (Pstream::master())
1904 "DF_" + name(
NUmodes) +
"_" + name(
1917 Eigen::MatrixXd
KF_matrix(KF1size, KF2size);
1919 for (label i = 0; i < KF1size; i++)
1921 for (label j = 0; j < KF2size; j++)
1923 volVectorField CoeffA = (fvc::reconstruct(
dt_dummy * fvc::snGrad(
1925 volVectorField CoeffB = fvc::reconstruct(
L_PHImodes[i]).ref();
1926 KF_matrix(i, j) = fvc::domainIntegrate(CoeffA & CoeffB).value();
1930 if (Pstream::master())
1933 "KF_" + name(
NUmodes) +
"_" + name(
1948 surfaceScalarField phi_tmp(
"Phi_tmp",
_phi());
1950 for (label i = 0; i <
NUmodes; i++)
1952 for (label j = 0; j < Csize1; j++)
1954 for (label k = 0; k < Csize1; k++)
1957 volVectorField CoeffA = fvc::reconstruct(phi_tmp).ref();
1958 volVectorField CoeffB = fvc::reconstruct(
L_PHImodes[j]).ref();
1959 Ci_tensor(i, j, k) = fvc::domainIntegrate(CoeffB & CoeffA).value();
1964 if (Pstream::master())
1980 List <Eigen::MatrixXd>
SD_matrix(BCsize);
1981 surfaceScalarField phi_tmp(
"Phi_tmp",
_phi());
1983 for (label i = 0; i < BCsize; i++)
1986 Vector<double> v(0, 0, 0);
1988 volVectorField Upara(
Uinl());
1993 for (label j = 0; j < SDsize; j++)
1995 volVectorField CoeffB = fvc::reconstruct(
L_PHImodes[j]).ref();
1997 volVectorField CoeffA = fvc::reconstruct(phi_tmp).ref();
1998 SD_matrix[i](j, 0) = fvc::domainIntegrate(CoeffA & CoeffB).value();
2014 List <Eigen::MatrixXd>
SC_matrix(BCsize);
2015 surfaceScalarField phi_tmp(
"Phi_tmp",
_phi());
2017 for (label i = 0; i < BCsize; i++)
2020 Vector<double> v(0, 0, 0);
2022 volVectorField Upara(
Uinl());
2027 for (label j = 0; j < SCsize; j++)
2029 volVectorField CoeffB = fvc::reconstruct(
L_PHImodes[j]).ref();
2030 phi_tmp =
dt_dummy * fvc::flux(fvc::div(fvc::flux(Upara), Upara));
2031 volVectorField CoeffA = fvc::reconstruct(phi_tmp).ref();
2032 SC_matrix[i](j, 0) = fvc::domainIntegrate(CoeffA & CoeffB).value();
2049 for (label i = 0; i <
NUmodes; i++)
2051 volVectorField A = fvc::reconstruct(
L_PHImodes[i]).ref();
2053 for (label j = 0; j <
NUmodes; j++)
2055 volVectorField B = fvc::reconstruct(
L_PHImodes[j]).ref();
2056 W_matrix(i, j) = fvc::domainIntegrate(A & B).value();
2069 volScalarField& ciao =
const_cast<volScalarField&
>(nu);
2072 for (label i = 0; i < ciao.boundaryFieldRef().size(); i++)
2086 instantList Times = runTime.times();
2087 fvMesh& mesh =
_mesh();
2088 volScalarField& p =
_p();
2089 volVectorField& U =
_U();
2091 IOdictionary FORCESdict
2098 IOobject::MUST_READ,
2102 IOdictionary transportProperties
2106 "transportProperties",
2109 IOobject::MUST_READ,
2113 word pName(FORCESdict.lookup(
"pName"));
2114 word UName(FORCESdict.lookup(
"UName"));
2115 functionObjects::ITHACAforces f(
"Forces", mesh, FORCESdict);
2123 f.calcForcesMoment();
2125 for (label j = 0; j < 3; j++)
2131 for (label i = 0; i <
NPmodes; i++)
2137 f.calcForcesMoment();
2139 for (label j = 0; j < 3; j++)
2141 nMatrix(i, j) = f.forcePressure()[j];
2145 if (para->exportPython)
2148 "./ITHACAoutput/Matrices/");
2152 if (para->exportMatlab)
2155 "./ITHACAoutput/Matrices/");
2159 if (para->exportTxt)
2162 "./ITHACAoutput/Matrices/");
2174 instantList Times = runTime.times();
2175 fvMesh& mesh =
_mesh();
2176 volScalarField& p =
_p();
2177 volVectorField& U =
_U();
2179 IOdictionary FORCESdict
2186 IOobject::MUST_READ,
2190 IOdictionary transportProperties
2194 "transportProperties",
2197 IOobject::MUST_READ,
2201 word pName(FORCESdict.lookup(
"pName"));
2202 word UName(FORCESdict.lookup(
"UName"));
2203 functionObjects::ITHACAforces f(
"Forces", mesh, FORCESdict);
2205 for (label i = 0; i < nModes; i++)
2211 f.calcForcesMoment();
2213 for (label j = 0; j < 3; j++)
2219 for (label i = 0; i < nModes; i++)
2225 f.calcForcesMoment();
2227 for (label j = 0; j < 3; j++)
2229 nMatrix(i, j) = f.forcePressure()[j];
2233 if (para->exportPython)
2236 "./ITHACAoutput/Matrices/");
2240 if (para->exportMatlab)
2243 "./ITHACAoutput/Matrices/");
2247 if (para->exportTxt)
2250 "./ITHACAoutput/Matrices/");
2256 const Eigen::MatrixXd& pressureCoeffs, fileName folder)
2258 M_Assert(velCoeffs.cols() ==
tauMatrix.rows(),
2259 "The number of velocity modes in the coefficients matrix is not equal to the number of modes in the viscous forces matrix.");
2260 M_Assert(pressureCoeffs.cols() ==
nMatrix.rows(),
2261 "The number of pressure modes in the coefficients matrix is not equal to the number of modes in the pressure forces matrix.");
2263 system(
"ln -s ../../constant " + folder +
"/constant");
2264 system(
"ln -s ../../0 " + folder +
"/0");
2265 system(
"ln -s ../../system " + folder +
"/system");
2267 IOdictionary FORCESdict
2274 IOobject::MUST_READ,
2278 Eigen::MatrixXd fTau;
2280 fTau.setZero(velCoeffs.rows(), 3);
2281 fN.setZero(pressureCoeffs.rows(), 3);
2283 fN = pressureCoeffs *
nMatrix;
2286 if (para->exportPython)
2292 if (para->exportMatlab)
2298 if (para->exportTxt)
2317 argList& args =
_args();
2319 runTime.setTime(0, 1);
2320 Foam::fvMesh& mesh =
_mesh();
2328 simpleControl& simple =
_simple();
2329 Info <<
"ReReading field p\n" << endl;
2330 _p = autoPtr<volScalarField>
2339 IOobject::MUST_READ,
2340 IOobject::AUTO_WRITE
2345 volScalarField& p =
_p();
2346 Info <<
"ReReading field U\n" << endl;
2347 _U = autoPtr<volVectorField>
2356 IOobject::MUST_READ,
2357 IOobject::AUTO_WRITE
2362 volVectorField& U =
_U();
2363 Info <<
"ReReading/calculating face flux field phi\n" << endl;
2364 _phi = autoPtr<surfaceScalarField>
2366 new surfaceScalarField
2373 IOobject::READ_IF_PRESENT,
2374 IOobject::AUTO_WRITE
2376 linearInterpolate(U) & mesh.Sf()
2379 surfaceScalarField& phi =
_phi();
2385 new singlePhaseTransportModel( U, phi )
2388 turbulence = autoPtr<incompressible::turbulenceModel>
2390 incompressible::turbulenceModel::New(U, phi, laminarTransport)
2392 _MRF = autoPtr<IOMRFZoneList>
2394 new IOMRFZoneList(mesh)
2396 _fvOptions = autoPtr<fv::options>(
new fv::options(mesh));
static Eigen::MatrixXd field2Eigen(GeometricField< Type, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
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.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
void writeMu(List< scalar > mu_now)
Write out a list of scalar corresponding to the parameters used in the truthSolve.
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
label counter
Counter used for the output of the full order solutions.
void assignBC(volVectorField &s, label BC_ind, Vector< double > &value)
Assign Boundary Condition to a volVectorField.
void assignIF(T &s, G &value)
Assign internal field condition.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
IOdictionary * ITHACAdict
dictionary to store input output infos
Eigen::MatrixXd mu
Row matrix of parameters.
autoPtr< argList > _args
argList
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
void change_viscosity(double mu)
Function to change the viscosity.
List< Eigen::MatrixXd > bcVelVec
Boundary term for penalty method - vector.
scalar maxIter
Number of maximum iterations to be done for the computation of the truth solution.
void restart()
set U and P back to the values into the 0 folder
List< Eigen::MatrixXd > pressure_gradient_term_linsys_div(label NPmodes)
Laplacian of pressure Linear System - Divergence term.
label NPmodes
Number of pressure modes used for the projection.
List< Eigen::MatrixXd > pressure_gradient_term_linsys_conv(label NPmodes)
Laplacian of pressure Linear System - Convection term.
bool supex
Boolean variable to check the existence of the supremizer modes.
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Poisson Equation for pressure.
Eigen::MatrixXd nMatrix
Pressure forces.
Eigen::MatrixXd BC1_matrix
PPE BC1.
Eigen::MatrixXd diffusive_term(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Term.
autoPtr< surfaceScalarField > _phi
Flux.
Eigen::MatrixXd BC4_matrix
PPE BC4.
Eigen::MatrixXd tauMatrix
Viscous forces.
Eigen::MatrixXd BC3_matrix
PPE BC3.
Eigen::Tensor< double, 3 > C_tensor
Diffusion term.
Eigen::Tensor< double, 3 > convective_term_tens_cache(label NUmodes, label NPmodes, label NSUPmodes)
Export convective term as a tensor using the cached procedure.
void forcesMatrices(label NUmodes, label NPmodes, label NSUPmodes)
Compute lift and drag matrices.
autoPtr< simpleControl > _simple
simpleControl
surfaceScalarModes L_PHImodes
List of pointers containing the total number of flux modes.
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
volVectorModes supmodes
List of pointers used to form the supremizer modes.
List< Eigen::MatrixXd > bcVelocityVec(label NUmodes, label NSUPmodes)
Boundary integral modes on boundary used by the penaly method.
List< Eigen::MatrixXd > boundary_vector_convection_consistent(label NUmodes, label NSUPmodes)
Boundary vector convection term - Consistent Flux Method.
Eigen::MatrixXd W_matrix
Mass Matrix New Time Step - Consistent Flux Method.
Eigen::MatrixXd mass_matrix_oldtime_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Mass Matrix old time step (consistent flux method).
List< Eigen::MatrixXd > pressure_gradient_term_linsys_diff(label NPmodes)
Laplacian of pressure Linear System - Diffusion term.
autoPtr< fv::options > _fvOptions
fvOptions
Eigen::MatrixXd pressure_BC4(label NPmodes, label NUmodes)
Term N° 4 given by the additional boundary condition using a PPE approach for time-dependent BCs.
void solvesupremizer(word type="snapshots")
solve the supremizer either with the use of the pressure snaphots or the pressure modes
autoPtr< Time > _runTime
Time.
Eigen::MatrixXd I_matrix
Mass Matrix Old Time Step - Consistent Flux Method.
volVectorModes Umodes
List of pointers used to form the velocity modes.
Eigen::MatrixXd diffusive_term_sym(label NUmodes, label NPmodes, label NSUPmodes)
Symetric diffusive Term.
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
List< Eigen::MatrixXd > boundary_vector_diffusion_consistent(label NUmodes, label NSUPmodes)
Boundary vector diffusion term (consistent flux method).
List< Eigen::MatrixXd > LinSysDiff
Projection Peqn onto Pressure modes - Diffusion term.
autoPtr< dimensionedScalar > dt_dummy
Dummy time step including unit.
Eigen::MatrixXd diffusive_term_flux_method(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Flux Method.
Eigen::MatrixXd diffusive_term_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Diffusion Term (consistent flux method).
Eigen::MatrixXd pressure_BC1(label NPmodes, label NUmodes)
Term N° 1 given by the additional boundary condition using a PPE approach.
steadyNS()
Null constructor.
Eigen::MatrixXd divergence_term(label NUmodes, label NPmodes, label NSUPmodes)
Divergence Term (supremizer approach).
List< Eigen::MatrixXd > G_matrix
Divergence of momentum PPE.
scalar tolerance
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and p...
Eigen::MatrixXd mass_term(label NUmodes, label NPmodes, label NSUPmodes)
Mass Term.
Eigen::Tensor< double, 3 > divMomentum_cache(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach) using the cached procedure.
autoPtr< fvMesh > _mesh
Mesh.
autoPtr< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model).
label NUmodes
Number of velocity modes used for the projection.
autoPtr< dimensionedScalar > nu_dummy
Dummy viscocity including unit.
List< Eigen::MatrixXd > boundary_vector_diffusion(label NUmodes, label NPmodes, label NSUPmodes)
Boundary vector diffusion term.
List< Eigen::MatrixXd > LinSysDiv
Projection Peqn onto Pressure modes - Divergence term.
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Eigen::Tensor< double, 3 > convective_term_flux_tens(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term.
Eigen::Tensor< double, 3 > convective_term_consistent_tens(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term (consistent flux method).
scalar pRefValue
Reference pressure value.
label pRefCell
Reference pressure cell.
word fluxMethod
Flux Method.
volVectorModes L_U_SUPmodes
List of pointers containing the total number of lift, supremizer and velocity modes.
Eigen::MatrixXd pressure_BC3(label NPmodes, label NUmodes)
Term N° 3 given by the additional boundary condition using a PPE approach.
Eigen::MatrixXd L_matrix
Convective background / Large scale advection term.
List< Eigen::MatrixXd > boundary_vector_convection(label NUmodes, label NPmodes, label NSUPmodes)
Boundary vector convection term.
Eigen::MatrixXd B_matrix
Diffusion term.
List< Eigen::MatrixXd > RD_matrix
Boundary term for diffusion term.
Eigen::MatrixXd D_matrix
Laplacian term PPE.
autoPtr< IOMRFZoneList > _MRF
MRF variable.
List< Eigen::MatrixXd > SD_matrix
Boundary term for diffusion term - Consistent Flux Method.
label NSUPmodes
Number of supremizer modes used for the projection.
Eigen::Tensor< double, 3 > gTensor
Divergence of momentum PPE.
Eigen::MatrixXd KF_matrix
Pressure Gradient Term - Consistent Flux Method.
Eigen::MatrixXd K_matrix
Gradient of pressure matrix.
List< Eigen::MatrixXd > C_matrix
Non linear term.
Eigen::Tensor< double, 3 > convective_term_tens(label NUmodes, label NPmodes, label NSUPmodes)
Export convective term as a tensor.
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
List< Eigen::MatrixXd > convective_term(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term.
Eigen::MatrixXd pressure_gradient_term_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Pressure Gradient Term (consistent flux method).
List< Eigen::MatrixXd > pressure_BC2(label NPmodes, label NUmodes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Eigen::MatrixXd divergent_convective_background(label NPmodes, label NUmodes, volVectorField vls)
Divergent of Large Scale / Background Advection.
List< Eigen::MatrixXd > LinSysConv
Projection Peqn onto Pressure modes - Convection term.
PtrList< volVectorField > supfield
List of pointers used to form the supremizer snapshots matrix.
Eigen::Tensor< double, 3 > divMomentum(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach).
List< Eigen::MatrixXd > SC_matrix
Boundary term for convection term - Consistent Flux Method.
Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPmodes, label NSUPmodes)
Gradient of pressure.
Eigen::MatrixXd DF_matrix
Diffusion Term - Consistent Flux Method.
List< Eigen::MatrixXd > bcVelMat
Boundary term for penalty method - matrix.
Eigen::MatrixXd P_matrix
Div of velocity.
List< Eigen::MatrixXd > RC_matrix
Boundary vector for convection term.
Eigen::MatrixXd BP_matrix
Diffusion term for flux method PPE.
Eigen::Tensor< double, 3 > Cf_tensor
Convection term for flux method.
Eigen::MatrixXd L_D_matrix
Divergent convective background / Large scale advection term.
Eigen::MatrixXd M_matrix
Mass Matrix.
List< Eigen::MatrixXd > BC2_matrix
PPE BC2.
autoPtr< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
Eigen::MatrixXd convective_background(label NUmodes, volVectorField vls)
Large Scale / Background Advection.
void discretizeThenProject(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Discretize-then-project approach.
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes)
Project using a supremizer approach.
Eigen::MatrixXd laplacian_pressure(label NPmodes)
Laplacian of pressure term (PPE approach).
surfaceScalarModes Phimodes
List of pointers used to form the flux modes.
Eigen::MatrixXd mass_matrix_newtime_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Mass Matrix new time step (consistent flux method).
List< Eigen::MatrixXd > bcVelocityMat(label NUmodes, label NSUPmodes)
Boundary integral modes on boundary used by the penaly method.
Eigen::Tensor< double, 3 > bc2Tensor
PPE BC2.
PtrList< volScalarField > liftfieldP
List of pointer used to form the list of lifting functions for the pressure.
autoPtr< volVectorField > _U
Velocity field.
Eigen::Tensor< double, 3 > pressureBC2(label NPmodes, label NUmodes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Eigen::Tensor< double, 3 > Ci_tensor
Convection term - Consistent Flux Method.
void liftSolve()
Perform a lift solve.
List< Eigen::MatrixXd > div_momentum(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach).
volScalarModes Pmodes
List of pointers used to form the pressure modes.
word bcMethod
Boundary Method.
void reconstructLiftAndDrag(const Eigen::MatrixXd &velCoeffs, const Eigen::MatrixXd &pressureCoeffs, fileName folder)
Method to reconstruct the forces using velocity and pressure coefficients.
autoPtr< volScalarField > _p
Pressure field.
void ReadDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Read a dense matrix from a binary format file.
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 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 ...
void SaveDenseMatrix(MatrixType &Matrix, word folder, word MatrixName)
Save a dense matrix to a binary format file.
void ReadDenseTensor(TensorType &Tensor, word folder, word MatrixName)
Read a dense tensor from 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 exportTensor(Eigen::Tensor< T, 3 > tensor, word Name, word type, word folder)
Export the reduced tensor in numpy (tipo=python), matlab (tipo=matlab) and txt (tipo=eigen) format.
void SaveDenseTensor(TensorType &Tensor, word folder, word MatrixName)
Save a dense tensor to file.
bool check_pod()
Check if the POD data folder "./ITHACAoutput/POD" exists.
bool check_off()
Check if the offline data folder "./ITHACAoutput/Offline" exists.
void assignZeroDirichlet(GeometricField< Type, fvPatchField, volMesh > &field)
Assign zero internal field.
bool check_folder(word folder)
Checks if a folder exists.
bool check_file(std::string fileName)
Function that returns true if a file exists.
bool check_sup()
Check if the supremizer folder exists.
void createSymLink(word folder)
Creates symbolic links to 0, system and constant.
void changeBCtype(GeometricField< Type, fvPatchField, volMesh > &field, word BCtype, label BC_ind)
Change the boundary condition type for a GeometricField.
Header file of the steadyNS class.