Loading...
Searching...
No Matches
17YJunction.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2021 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Description
25 Example of an unsteady NS Reduction Problem with time-dependent boundary
26 conditions.
27SourceFiles
28 17YJunction.C
29\*---------------------------------------------------------------------------*/
30#include "unsteadyNS.H"
31#include "ITHACAPOD.H"
32#include "ReducedUnsteadyNS.H"
33#include "ITHACAstream.H"
34#include <chrono>
35#include <math.h>
36#include <iomanip>
37
39{
40 public:
41 explicit tutorial17(int argc, char* argv[])
42 :
43 unsteadyNS(argc, argv),
44 U(_U()),
45 p(_p())
46 {}
47
48 // Fields To Perform
49 volVectorField& U;
50 volScalarField& p;
51
53 {
54 List<scalar> mu_now(1);
55 volVectorField U0 = U;
56 volScalarField P0 = p;
57
58 if (offline)
59 {
60 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
61 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
62 }
63 else
64 {
65 U = U0;
66 p = P0;
67 mu_now[0] = mu(0, 0);
68 truthSolve(mu_now);
69 }
70 }
71
72 // Method to compute the lifting functions for this tutorial
73 void liftSolve(Eigen::MatrixXd BCs)
74 {
75 for (label k = 0; k < inletPatch.rows(); k++)
76 {
77 Time& runTime = _runTime();
78 surfaceScalarField& phi = _phi();
79 fvMesh& mesh = _mesh();
80 volScalarField p = _p();
81 volVectorField U = _U();
82 IOMRFZoneList& MRF = _MRF();
83 label BCind = inletPatch(k, 0);
84 volVectorField Ulift("Ulift" + name(k), U);
85 instantList Times = runTime.times();
86 runTime.setTime(Times[1], 1);
87 pisoControl potentialFlow(mesh, "potentialFlow");
88 Info << "Solving a lifting Problem" << endl;
89 Vector<double> v1(0, 0, 0);
90 v1[0] = BCs(0, k);
91 v1[1] = BCs(1, k);
92 Vector<double> v0(0, 0, 0);
93
94 for (label j = 0; j < U.boundaryField().size(); j++)
95 {
96 if (j == BCind)
97 {
98 assignBC(Ulift, j, v1);
99 }
100 else if (U.boundaryField()[BCind].type() == "fixedValue")
101 {
102 assignBC(Ulift, j, v0);
103 }
104 else
105 {
106 }
107
108 assignIF(Ulift, v0);
109 phi = linearInterpolate(Ulift) & mesh.Sf();
110 }
111
112 Info << "Constructing velocity potential field Phi\n" << endl;
113 volScalarField Phi
114 (
115 IOobject
116 (
117 "Phi",
118 runTime.timeName(),
119 mesh,
120 IOobject::READ_IF_PRESENT,
121 IOobject::NO_WRITE
122 ),
123 mesh,
124 dimensionedScalar("Phi", dimLength * dimVelocity, 0),
125 p.boundaryField().types()
126 );
127 label PhiRefCell = 0;
128 scalar PhiRefValue = 0;
130 (
131 Phi,
132 potentialFlow.dict(),
133 PhiRefCell,
134 PhiRefValue
135 );
136 mesh.setFluxRequired(Phi.name());
137 runTime.functionObjects().start();
138 MRF.makeRelative(phi);
139 adjustPhi(phi, Ulift, p);
140
141 while (potentialFlow.correctNonOrthogonal())
142 {
143 fvScalarMatrix PhiEqn
144 (
145 fvm::laplacian(dimensionedScalar("1", dimless, 1), Phi)
146 ==
147 fvc::div(phi)
148 );
149 PhiEqn.setReference(PhiRefCell, PhiRefValue);
150 PhiEqn.solve();
151
152 if (potentialFlow.finalNonOrthogonalIter())
153 {
154 phi -= PhiEqn.flux();
155 }
156 }
157
158 MRF.makeAbsolute(phi);
159 Info << "Continuity error = "
160 << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
161 << endl;
162 Ulift = fvc::reconstruct(phi);
163 Ulift.correctBoundaryConditions();
164 Info << "Interpolated velocity error = "
165 << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
166 / sum(mesh.magSf())).value()
167 << endl;
168 Ulift.write();
169 volVectorField Uzero
170 (
171 IOobject
172 (
173 "Uzero",
174 U.time().timeName(),
175 U.mesh(),
176 IOobject::NO_READ,
177 IOobject::AUTO_WRITE
178 ),
179 U.mesh(),
180 dimensionedVector("zero", U.dimensions(), vector::zero)
181 );
182 volVectorField Uliftx("Uliftx" + name(k), Uzero);
183 Uliftx.replace(0, Ulift.component(0));
184 Uliftx.write();
185 liftfield.append((Uliftx).clone());
186 volVectorField Ulifty("Ulifty" + name(k), Uzero);
187 Ulifty.replace(1, Ulift.component(1));
188 Ulifty.write();
189 liftfield.append((Ulifty).clone());
190 }
191 }
192
193};
194
195/*---------------------------------------------------------------------------*\
196 Starting the MAIN
197\*---------------------------------------------------------------------------*/
198int main(int argc, char* argv[])
199{
200 // Construct the tutorial17 object
201 tutorial17 example(argc, argv);
202 // the offline samples for the boundary conditions
203 word par_offline_BC("./timeBCoff");
204 // Convert boundary values x-direction and y-direction to x,y,z
205 Eigen::MatrixXd timeBCoff2D = ITHACAstream::readMatrix(par_offline_BC);
206 Eigen::MatrixXd timeBCoff3D = Eigen::MatrixXd::Zero(6, timeBCoff2D.cols());
207 timeBCoff3D.row(0) = timeBCoff2D.row(0);
208 timeBCoff3D.row(1) = timeBCoff2D.row(1);
209 timeBCoff3D.row(3) = timeBCoff2D.row(2);
210 timeBCoff3D.row(4) = timeBCoff2D.row(3);
211 example.timeBCoff = timeBCoff3D;
212 Eigen::MatrixXd par_on_BC;
213 word par_online_BC("./timeBCon");
214 par_on_BC = ITHACAstream::readMatrix(par_online_BC);
215 // Read parameters from ITHACAdict file
217 example._runTime());
218 int NmodesUout = para->ITHACAdict->lookupOrDefault<int>("NmodesUout", 15);
219 int NmodesPout = para->ITHACAdict->lookupOrDefault<int>("NmodesPout", 15);
220 int NmodesSUPout = para->ITHACAdict->lookupOrDefault<int>("NmodesSUPout", 15);
221 int NmodesUproj = para->ITHACAdict->lookupOrDefault<int>("NmodesUproj", 10);
222 int NmodesPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesPproj", 10);
223 int NmodesSUPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesSUPproj", 10);
224 int NmodesOut = para->ITHACAdict->lookupOrDefault<int>("NmodesOut", 20);
226 example.Pnumber = 1;
228 example.Tnumber = 1;
230 example.setParameters();
231 // Set the parameter ranges
232 example.mu_range(0, 0) = 0.01;
233 example.mu_range(0, 1) = 0.01;
234 // Generate equispaced samples inside the parameter range
235 example.genEquiPar();
236 // Set the inlet boundaries where we have non homogeneous boundary conditions
237 example.inletIndex.resize(4, 2); // rows: total number of patches
238 example.inletIndex(0, 0) = 1; // Patch inlet 1
239 example.inletIndex(0, 1) = 0; // Patch inlet 1: x-direction
240 example.inletIndex(1, 0) = 1; // Patch inlet 1: y-direction
241 example.inletIndex(1, 1) = 1; // Patch inlet 2
242 example.inletIndex(2, 0) = 2; // Patch inlet 2: x-direction
243 example.inletIndex(2, 1) = 0; // Patch inlet 2: y-direction
244 example.inletIndex(3, 0) = 2; // Patch inlet 2: x-direction
245 example.inletIndex(3, 1) = 1; // Patch inlet 2: y-direction
246 example.inletPatch.resize(2, 1);
247 example.inletPatch(0, 0) = example.inletIndex(0, 0); // Patch inlet 1
248 example.inletPatch(1, 0) = example.inletIndex(2, 0); // Patch inlet 2
249 // Time parameters
250 example.startTime = 0;
251 example.finalTime = 12;
252 example.timeStep = 0.0005;
253 example.writeEvery = 0.03;
254 // Perform The Offline Solve;
255 example.offlineSolve();
256
257 // Perform POD
258 if (example.bcMethod == "lift")
259 {
260 // Search the lift function
261 Eigen::MatrixXd BCs;
262 BCs.resize(2, 2);
263 BCs(0, 0) = 1; // Patch inlet 1
264 BCs(1, 0) = -1; // Patch inlet 2
265 BCs(0, 1) = 1; // Patch inlet 1
266 BCs(1, 1) = 1; // Patch inlet 2
267 example.liftSolve(BCs);
268 // Normalize the lifting function
270 // Create homogeneous basis functions for velocity
271 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
272 // Perform a POD decomposition for velocity and pressure
273 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
274 example.podex, 0, 0,
275 NmodesUout);
276 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
277 example.podex, 0, 0,
278 NmodesPout);
279 }
280 else if (example.bcMethod == "penalty")
281 {
282 // Perform a POD decomposition for velocity and pressure
283 ITHACAPOD::getModes(example.Ufield, example.Umodes, example._U().name(),
284 example.podex, 0, 0,
285 NmodesUout);
286 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
287 example.podex, 0, 0,
288 NmodesPout);
289 }
290
291 // Reduced Matrices
292 example.projectPPE("./Matrices", NmodesUproj, NmodesPproj, NmodesSUPproj);
293 reducedUnsteadyNS reduced(example);
294 // Set values of the online phase
295 reduced.nu = 0.01;
296 reduced.tstart = 0;
297 reduced.finalTime = 18;
298 reduced.dt = 0.0005;
299 reduced.storeEvery = 0.0005;
300 reduced.exportEvery = 0.03;
301 // Set values velocity boundary conditions of the online phase
302 Eigen::MatrixXd vel_now = par_on_BC;
303
304 if (example.bcMethod == "penalty")
305 {
306 // Set values of the iterative penalty method
307 reduced.maxIterPenalty = 100;
308 reduced.tolerancePenalty = 1e-5;
309 reduced.timeStepPenalty = 5;
310 // Set initial quess for penalty factors
311 reduced.tauIter = Eigen::MatrixXd::Zero(4, 1);
312 reduced.tauIter << 1e-6, 1e-6, 1e-6, 1e-6;
313 // Solve for the penalty factors with the iterative solver
314 reduced.tauU = reduced.penalty_PPE(vel_now, reduced.tauIter);
315 }
316
317 // Set the online temperature BC and solve reduced model
318 reduced.solveOnline_PPE(vel_now);
319 reduced.reconstruct(false, "./ITHACAoutput/Reconstruction/");
320 exit(0);
321}
322
510
511
512
513
int main(int argc, char *argv[])
Header file of the ITHACAPOD class.
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the reducedUnsteadyNS class.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
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.
scalar startTime
Start Time (initial time to start storing the snapshots)
scalar writeEvery
Time step of the writing procedure.
scalar timeStep
Time step of the simulation.
scalar finalTime
Final time (final time of the simulation and consequently of the acquisition of the snapshots)
scalar nu
Reduced viscosity in case of parametrized viscosity.
Eigen::MatrixXd tauU
Penalty Factor.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
double exportEvery
A variable for exporting the fields.
scalar tolerancePenalty
Tolerance for the residual of the boundary values, there is the same tolerance for velocity and tempe...
scalar finalTime
Scalar to store the final time if the online simulation.
int timeStepPenalty
Number of timesteps calculated for the iterative penalty method.
void solveOnline_PPE(Eigen::MatrixXd vel_now, int startSnap=0)
Method to perform an online solve using a PPE stabilisation method.
scalar tstart
Scalar to store the initial time if the online simulation.
void reconstruct(bool exportFields=false, fileName folder="./online_rec")
Method to reconstruct the solutions from an online solve with a supremizer stabilisation technique.
Eigen::MatrixXd tauIter
Penalty Factor determined with iterative solver.
scalar maxIterPenalty
Maximum number of iterations to be done for the computation of the penalty factor.
double dt
Scalar to store the time increment.
Eigen::MatrixXd penalty_PPE(Eigen::MatrixXd &vel_now, Eigen::MatrixXd &tauIter, int startSnap=0)
Method to determine the penalty factors iteratively.
double storeEvery
A variable for storing the reduced coefficients.
Eigen::MatrixXi inletPatch
Matrix that contains informations about the inlet boundaries without specifing the direction Rows = N...
label Pnumber
Number of parameters.
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.
label Tnumber
Dimension of the training set (used only when gerating parameters without input)
void computeLift(T &Lfield, T &liftfield, T &omfield)
Homogenize the snapshot matrix, it works with PtrList of volVectorField and volScalarField.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
Eigen::MatrixXd mu
Row matrix of parameters.
Eigen::MatrixXd mu_range
Range of the parameter spaces.
void setParameters()
Set Parameters Problems.
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 genEquiPar()
Generate Equidistributed Numbers.
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Poisson Equation for pressure.
Definition steadyNS.C:353
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:290
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition steadyNS.H:101
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition steadyNS.H:89
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:281
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
PtrList< volVectorField > Uomfield
List of pointers used to form the homogeneous velocity snapshots.
Definition steadyNS.H:113
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:305
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
word bcMethod
Boundary Method.
Definition steadyNS.H:317
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
tutorial17(int argc, char *argv[])
Definition 17YJunction.C:41
void liftSolve(Eigen::MatrixXd BCs)
Definition 17YJunction.C:73
void offlineSolve()
Definition 17YJunction.C:52
volScalarField & p
Definition 17YJunction.C:50
volVectorField & U
Definition 17YJunction.C:49
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
Definition unsteadyNS.H:62
Eigen::MatrixXd timeBCoff
Definition unsteadyNS.H:96
volScalarField v0(v)
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.
Definition ITHACAPOD.C:93
List< Eigen::MatrixXd > readMatrix(word folder, word mat_name)
Read a three dimensional matrix from a txt file in Eigen format.
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 normalizeFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &fields)
Normalize list of Geometric fields.
surfaceScalarField & phi
volVectorField U0(U)
setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue)
adjustPhi(phiHbyA, U, p)
Header file of the unsteadyNS class.