Loading...
Searching...
No Matches
SteadyNSTurbIntrusive.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) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
31#include "steadyNS.H"
33#include "viscosityModel.H"
34
35// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
36
37// Constructor
40{
41 _args = autoPtr<argList>
42 (
43 new argList(argc, argv)
44 );
45
46 if (!_args->checkRootCase())
47 {
48 Foam::FatalError.exit();
49 }
50
51 argList& args = _args();
52#include "createTime.H"
53#include "createMesh.H"
54 _simple = autoPtr<simpleControl>
55 (
56 new simpleControl
57 (
58 mesh
59 )
60 );
61 simpleControl& simple = _simple();
62#include "createFields.H"
63#include "createFvOptions.H"
64 //
65#pragma GCC diagnostic push
66#pragma GCC diagnostic ignored "-Wunused-variable"
67#include "initContinuityErrs.H"
68#pragma GCC diagnostic pop
69 //
70 turbulence->validate();
71 ITHACAdict = new IOdictionary
72 (
73 IOobject
74 (
75 "ITHACAdict",
76 runTime.system(),
77 mesh,
78 IOobject::MUST_READ,
79 IOobject::NO_WRITE
80 )
81 );
82 tolerance = ITHACAdict->lookupOrDefault<scalar>("tolerance", 1e-5);
83 maxIter = ITHACAdict->lookupOrDefault<scalar>("maxIter", 1000);
84 bcMethod = ITHACAdict->lookupOrDefault<word>("bcMethod", "lift");
85 M_Assert(bcMethod == "lift" || bcMethod == "penalty",
86 "The BC method must be set to lift or penalty in ITHACAdict");
91}
92
93
94// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
95
96// Method to performa a truthSolve
97void SteadyNSTurbIntrusive::truthSolve(List<scalar> mu_now)
98{
99 Time& runTime = _runTime();
100 fvMesh& mesh = _mesh();
101 volScalarField& p = _p();
102 volVectorField& U = _U();
103 surfaceScalarField& phi = _phi();
104 fv::options& fvOptions = _fvOptions();
105 simpleControl& simple = _simple();
106 IOMRFZoneList& MRF = _MRF();
107 singlePhaseTransportModel& laminarTransport = _laminarTransport();
109 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
110 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
111 volScalarField _nut(turbulence->nut());
112 ITHACAstream::exportSolution(_nut, name(counter), "./ITHACAoutput/Offline/");
113 Ufield.append(U.clone());
114 Pfield.append(p.clone());
115 nutFields.append(_nut.clone());
116 counter++;
117 writeMu(mu_now);
118 // --- Fill in the mu_samples with parameters (mu) to be used for the PODI sample points
119 mu_samples.conservativeResize(mu_samples.rows() + 1, mu_now.size());
120
121 for (label i = 0; i < mu_now.size(); i++)
122 {
123 mu_samples(mu_samples.rows() - 1, i) = mu_now[i];
124 }
125
126 // Resize to Unitary if not initialized by user (i.e. non-parametric problem)
127 if (mu.cols() == 0)
128 {
129 mu.resize(1, 1);
130 }
131
132 if (mu_samples.rows() == mu.cols())
133 {
134 ITHACAstream::exportMatrix(mu_samples, "mu_samples", "eigen",
135 "./ITHACAoutput/Offline");
136 }
137}
138
139Eigen::Tensor<double, 3> SteadyNSTurbIntrusive::turbulenceTensor1(label nModes)
140{
141 Eigen::Tensor<double, 3> ct1Tensor;
142 ct1Tensor.resize(nModes, nModes, nModes);
143
144 for (label i = 0; i < nModes; i++)
145 {
146 for (label j = 0; j < nModes; j++)
147 {
148 for (label k = 0; k < nModes; k++)
149 {
150 ct1Tensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & fvc::laplacian(
151 nutModes[j], Umodes[k])).value();
152 }
153 }
154 }
155
156 if (Pstream::parRun())
157 {
158 reduce(ct1Tensor, sumOp<Eigen::Tensor<double, 3>>());
159 }
160
161 // Export the tensor
162 ITHACAstream::SaveDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/",
163 "ct1_" + name(nModes) + "_t");
164 return ct1Tensor;
165}
166
167Eigen::Tensor<double, 3> SteadyNSTurbIntrusive::turbulenceTensor2(label nModes)
168{
169 Eigen::Tensor<double, 3> ct2Tensor;
170 ct2Tensor.resize(nModes, nModes, nModes);
171
172 for (label i = 0; i < nModes; i++)
173 {
174 for (label j = 0; j < nModes; j++)
175 {
176 for (label k = 0; k < nModes; k++)
177 {
178 ct2Tensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & (fvc::div(
179 nutModes[j] * dev((fvc::grad(Umodes[k]))().T())))).value();
180 }
181 }
182 }
183
184 if (Pstream::parRun())
185 {
186 reduce(ct2Tensor, sumOp<Eigen::Tensor<double, 3>>());
187 }
188
189 // Export the tensor
190 ITHACAstream::SaveDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/",
191 "ct2_" + name(nModes) + "_t");
192 return ct2Tensor;
193}
194
195Eigen::MatrixXd SteadyNSTurbIntrusive::btTurbulence(label nModes)
196{
197 Eigen::MatrixXd btMatrix(nModes, nModes);
198 btMatrix = btMatrix * 0;
199
200 // Project everything
201 for (label i = 0; i < nModes; i++)
202 {
203 for (label j = 0; j < nModes; j++)
204 {
205 btMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & (fvc::div(dev((T(
206 fvc::grad(
207 Umodes[j]))))))).value();
208 }
209 }
210
211 // Export the matrix
212 ITHACAstream::SaveDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/",
213 "bt_" + name(nModes) + "_t");
214 return btMatrix;
215}
216
217void SteadyNSTurbIntrusive::project(fileName folder, label nModes)
218{
219 nModesOnline = nModes;
220
221 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
222 {
223 word bStr = "b_" + name(nModes);
224
225 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bStr))
226 {
227 ITHACAstream::ReadDenseMatrix(bMatrix, "./ITHACAoutput/Matrices/", bStr);
228 }
229 else
230 {
231 bMatrix = diffusiveTerm(nModes);
232 }
233
234 word btStr = "bt_" + name(nModes);
235
236 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + btStr))
237 {
238 ITHACAstream::ReadDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/", btStr);
239 }
240 else
241 {
242 btMatrix = btTurbulence(nModes);
243 }
244
245 word kStr = "k_" + name(nModes);
246
247 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + kStr))
248 {
249 ITHACAstream::ReadDenseMatrix(kMatrix, "./ITHACAoutput/Matrices/", kStr);
250 }
251 else
252 {
254 }
255
256 word cStr = "c_" + name(nModes) + "_t";
257
258 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + cStr))
259 {
260 ITHACAstream::ReadDenseTensor(convTensor, "./ITHACAoutput/Matrices/", cStr);
261 }
262 else
263 {
264 convTensor = convectiveTerm(nModes);
265 }
266
267 word ct1Str = "ct1_" + name(nModes) + "_t";
268
269 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct1Str))
270 {
271 ITHACAstream::ReadDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/", ct1Str);
272 }
273 else
274 {
276 }
277
278 word ct2Str = "ct2_" + name(nModes) + "_t";
279
280 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct2Str))
281 {
282 ITHACAstream::ReadDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/", ct2Str);
283 }
284 else
285 {
287 }
288
289 if (bcMethod == "penalty")
290 {
291 bcVelVec = bcVelocityVec(nModes);
292 bcVelMat = bcVelocityMat(nModes);
293 }
294 }
295 else
296 {
297 bMatrix = diffusiveTerm(nModes);
298 convTensor = convectiveTerm(nModes);
300 btMatrix = btTurbulence(nModes);
303
304 if (bcMethod == "penalty")
305 {
306 bcVelVec = bcVelocityVec(nModes);
307 bcVelMat = bcVelocityMat(nModes);
308 }
309 }
310
311 // Export the matrices
312 if (para->exportPython)
313 {
314 ITHACAstream::exportMatrix(bMatrix, "b", "python", "./ITHACAoutput/Matrices/");
315 ITHACAstream::exportMatrix(kMatrix, "k", "python", "./ITHACAoutput/Matrices/");
316 ITHACAstream::exportMatrix(btMatrix, "bt", "python",
317 "./ITHACAoutput/Matrices/");
319 "./ITHACAoutput/Matrices/");
320 ITHACAstream::exportTensor(ct1Tensor, "ct1", "python",
321 "./ITHACAoutput/Matrices/");
322 ITHACAstream::exportTensor(ct2Tensor, "ct2", "python",
323 "./ITHACAoutput/Matrices/");
324 }
325
326 if (para->exportMatlab)
327 {
328 ITHACAstream::exportMatrix(bMatrix, "b", "matlab", "./ITHACAoutput/Matrices/");
329 ITHACAstream::exportMatrix(kMatrix, "k", "matlab", "./ITHACAoutput/Matrices/");
330 ITHACAstream::exportMatrix(btMatrix, "bt", "matlab",
331 "./ITHACAoutput/Matrices/");
333 "./ITHACAoutput/Matrices/");
334 ITHACAstream::exportTensor(ct1Tensor, "ct1", "matlab",
335 "./ITHACAoutput/Matrices/");
336 ITHACAstream::exportTensor(ct2Tensor, "ct2", "matlab",
337 "./ITHACAoutput/Matrices/");
338 }
339
340 if (para->exportTxt)
341 {
342 ITHACAstream::exportMatrix(bMatrix, "b", "eigen", "./ITHACAoutput/Matrices/");
343 ITHACAstream::exportMatrix(kMatrix, "k", "eigen", "./ITHACAoutput/Matrices/");
344 ITHACAstream::exportMatrix(btMatrix, "bt", "eigen", "./ITHACAoutput/Matrices/");
346 "./ITHACAoutput/Matrices/c");
347 ITHACAstream::exportTensor(ct1Tensor, "ct1_", "eigen",
348 "./ITHACAoutput/Matrices/ct1");
349 ITHACAstream::exportTensor(ct2Tensor, "ct2_", "eigen",
350 "./ITHACAoutput/Matrices/ct2");
351 }
352
354 cTotalTensor.resize(nModes, nModes, nModes);
356}
357
358Eigen::MatrixXd SteadyNSTurbIntrusive::diffusiveTerm(label nModes)
359{
360 Eigen::MatrixXd bMatrix;
361 bMatrix.resize(nModes, nModes);
362
363 // Project everything
364 for (label i = 0; i < nModes; i++)
365 {
366 for (label j = 0; j < nModes; j++)
367 {
368 bMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & fvc::laplacian(
369 dimensionedScalar("1", dimless, 1), Umodes[j])).value();
370 }
371 }
372
373 if (Pstream::parRun())
374 {
375 reduce(bMatrix, sumOp<Eigen::MatrixXd>());
376 }
377
378 ITHACAstream::SaveDenseMatrix(bMatrix, "./ITHACAoutput/Matrices/",
379 "b_" + name(nModes));
380 return bMatrix;
381}
382
383Eigen::Tensor<double, 3> SteadyNSTurbIntrusive::convectiveTerm(label nModes)
384{
385 Eigen::Tensor<double, 3> convTensor;
386 convTensor.resize(nModes, nModes, nModes);
387
388 for (label i = 0; i < nModes; i++)
389 {
390 for (label j = 0; j < nModes; j++)
391 {
392 for (label k = 0; k < nModes; k++)
393 {
394 convTensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & fvc::div(
395 linearInterpolate(Umodes[j]) & Umodes[j].mesh().Sf(),
396 Umodes[k])).value();
397 }
398 }
399 }
400
401 if (Pstream::parRun())
402 {
403 reduce(convTensor, sumOp<Eigen::Tensor<double, 3>>());
404 }
405
406 // Export the tensor
407 ITHACAstream::SaveDenseTensor(convTensor, "./ITHACAoutput/Matrices/",
408 "c_" + name(nModes) + "_t");
409 return convTensor;
410}
411
413{
414 Eigen::MatrixXd kMatrix(nModes, nModes);
415
416 // Project everything
417 for (label i = 0; i < nModes; i++)
418 {
419 for (label j = 0; j < nModes; j++)
420 {
421 kMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & fvc::grad(
422 Pmodes[j])).value();
423 }
424 }
425
426 if (Pstream::parRun())
427 {
428 reduce(kMatrix, sumOp<Eigen::MatrixXd>());
429 }
430
431 ITHACAstream::SaveDenseMatrix(kMatrix, "./ITHACAoutput/Matrices/",
432 "k_" + name(nModes));
433 return kMatrix;
434}
435
436List< Eigen::MatrixXd > SteadyNSTurbIntrusive::bcVelocityVec(label nModes)
437{
438 List < Eigen::MatrixXd > bcVelVec(inletIndex.rows());
439
440 for (label j = 0; j < inletIndex.rows(); j++)
441 {
442 bcVelVec[j].resize(nModes, 1);
443 }
444
445 for (label k = 0; k < inletIndex.rows(); k++)
446 {
447 label BCind = inletIndex(k, 0);
448 label BCcomp = inletIndex(k, 1);
449
450 for (label i = 0; i < nModes; i++)
451 {
452 bcVelVec[k](i, 0) = gSum(Umodes[i].boundaryField()[BCind].component(
453 BCcomp));
454 }
455 }
456
457 ITHACAstream::exportMatrix(bcVelVec, "bcVelVec", "eigen",
458 "./ITHACAoutput/Matrices/bcVelVec");
459 return bcVelVec;
460}
461
462List< Eigen::MatrixXd > SteadyNSTurbIntrusive::bcVelocityMat(label nModes)
463{
464 label BCUsize = inletIndex.rows();
465 List < Eigen::MatrixXd > bcVelMat(BCUsize);
466
467 for (label j = 0; j < inletIndex.rows(); j++)
468 {
469 bcVelMat[j].resize(nModes, nModes);
470 }
471
472 for (label k = 0; k < inletIndex.rows(); k++)
473 {
474 label BCind = inletIndex(k, 0);
475 label BCcomp = inletIndex(k, 1);
476
477 for (label i = 0; i < nModes; i++)
478 {
479 for (label j = 0; j < nModes; j++)
480 {
481 bcVelMat[k](i, j) = gSum(Umodes[i].boundaryField()[BCind].component(
482 BCcomp) *
483 Umodes[j].boundaryField()[BCind].component(BCcomp));
484 }
485 }
486 }
487
488 ITHACAstream::exportMatrix(bcVelMat, "bcVelMat", "eigen",
489 "./ITHACAoutput/Matrices/bcVelMat");
490 return bcVelMat;
491}
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
#define M_Assert(Expr, Msg)
Header file of the SteadyNSTurbIntrusive class.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
Eigen::MatrixXd diffusiveTerm(label nModes)
Diffusive Term.
Eigen::MatrixXd bTotalMatrix
Total B Matrix.
PtrList< volScalarField > nutModes
List of POD modes for eddy viscosity.
Eigen::MatrixXd btTurbulence(label nModes)
bt added matrix for the turbulence treatement
Eigen::Tensor< double, 3 > turbulenceTensor2(label nModes)
Method to compute one of the turbulence eddy viscosity tensors.
Eigen::Tensor< double, 3 > convTensor
Convective tensor.
PtrList< volScalarField > nutFields
List of snapshots for the solution for eddy viscosity.
List< Eigen::MatrixXd > bcVelocityVec(label nModes)
Boundary integral modes on boundary used by the penaly method.
Eigen::Tensor< double, 3 > cTotalTensor
Total Turbulent tensor.
Eigen::Tensor< double, 3 > convectiveTerm(label nModes)
The method for computing the convective term tensor.
autoPtr< volScalarField > _nut
Eddy viscosity field.
Eigen::MatrixXd kMatrix
Pressure Gradient matrix.
Eigen::MatrixXd btMatrix
Turbulent viscosity matrix.
List< Eigen::MatrixXd > bcVelocityMat(label nModes)
Boundary integral modes on boundary used by the penaly method.
Eigen::Tensor< double, 3 > turbulenceTensor1(label nModes)
Method to compute one of the turbulence eddy viscosity tensors.
Eigen::MatrixXd pressureGradientTerm(label nModes)
The method for computing the pressure gradient term with number of modes of pressure being equal to t...
label nModesOnline
Number of modes used in the online stage.
Eigen::MatrixXd bMatrix
Diffusive matrix.
Eigen::Tensor< double, 3 > ct2Tensor
Turbulent viscosity tensor.
Eigen::Tensor< double, 3 > ct1Tensor
Turbulent viscosity tensor.
void project()
General projection operation.
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.
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.
List< Eigen::MatrixXd > bcVelVec
Boundary term for penalty method - vector.
Definition steadyNS.H:218
scalar maxIter
Number of maximum iterations to be done for the computation of the truth solution.
Definition steadyNS.H:125
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:256
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:284
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
autoPtr< fv::options > _fvOptions
fvOptions
Definition steadyNS.H:287
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
scalar tolerance
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and p...
Definition steadyNS.H:122
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:281
autoPtr< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model)
Definition steadyNS.H:302
ITHACAparameters * para
Definition steadyNS.H:82
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:305
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition steadyNS.H:299
List< Eigen::MatrixXd > bcVelMat
Boundary term for penalty method - matrix.
Definition steadyNS.H:221
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
volScalarField & T
Definition createT.H:46
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 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.
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.
simpleControl simple(mesh)
surfaceScalarField & phi
volVectorField & U
volScalarField & p
singlePhaseTransportModel & laminarTransport
label i
Definition pEqn.H:46
Header file of the steadyNS class.