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
39
41{
42 _args = autoPtr<argList>
43 (
44 new argList(argc, argv)
45 );
46
47 if (!_args->checkRootCase())
48 {
49 Foam::FatalError.exit();
50 }
51
52 argList& args = _args();
53#include "createTime.H"
54#include "createMesh.H"
55 _simple = autoPtr<simpleControl>
56 (
57 new simpleControl
58 (
59 mesh
60 )
61 );
62 simpleControl& simple = _simple();
63#include "createFields.H"
64#include "createFvOptions.H"
65 //
66#pragma GCC diagnostic push
67#pragma GCC diagnostic ignored "-Wunused-variable"
68#include "initContinuityErrs.H"
69#pragma GCC diagnostic pop
70 //
71 turbulence->validate();
72 ITHACAdict = new IOdictionary
73 (
74 IOobject
75 (
76 "ITHACAdict",
77 runTime.system(),
78 mesh,
79 IOobject::MUST_READ,
80 IOobject::NO_WRITE
81 )
82 );
83 tolerance = ITHACAdict->lookupOrDefault<scalar>("tolerance", 1e-5);
84 maxIter = ITHACAdict->lookupOrDefault<scalar>("maxIter", 1000);
85 bcMethod = ITHACAdict->lookupOrDefault<word>("bcMethod", "lift");
86 M_Assert(bcMethod == "lift" || bcMethod == "penalty",
87 "The BC method must be set to lift or penalty in ITHACAdict");
92}
93
94
95// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
96
97// Method to performa a truthSolve
98void SteadyNSTurbIntrusive::truthSolve(List<scalar> mu_now)
99{
100 Time& runTime = _runTime();
101 fvMesh& mesh = _mesh();
102 volScalarField& p = _p();
103 volVectorField& U = _U();
104 surfaceScalarField& phi = _phi();
105 fv::options& fvOptions = _fvOptions();
106 simpleControl& simple = _simple();
107 IOMRFZoneList& MRF = _MRF();
108 singlePhaseTransportModel& laminarTransport = _laminarTransport();
110 ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
111 ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
112 volScalarField _nut(turbulence->nut());
113 ITHACAstream::exportSolution(_nut, name(counter), "./ITHACAoutput/Offline/");
114 Ufield.append(U.clone());
115 Pfield.append(p.clone());
116 nutFields.append(_nut.clone());
117 counter++;
118 writeMu(mu_now);
119 // --- Fill in the mu_samples with parameters (mu) to be used for the PODI sample points
120 mu_samples.conservativeResize(mu_samples.rows() + 1, mu_now.size());
121
122 for (label i = 0; i < mu_now.size(); i++)
123 {
124 mu_samples(mu_samples.rows() - 1, i) = mu_now[i];
125 }
126
127 // Resize to Unitary if not initialized by user (i.e. non-parametric problem)
128 if (mu.cols() == 0)
129 {
130 mu.resize(1, 1);
131 }
132
133 if (mu_samples.rows() == mu.cols())
134 {
135 ITHACAstream::exportMatrix(mu_samples, "mu_samples", "eigen",
136 "./ITHACAoutput/Offline");
137 }
138}
139
140Eigen::Tensor<double, 3> SteadyNSTurbIntrusive::turbulenceTensor1(label nModes)
141{
142 Eigen::Tensor<double, 3> ct1Tensor;
143 ct1Tensor.resize(nModes, nModes, nModes);
144
145 for (label i = 0; i < nModes; i++)
146 {
147 for (label j = 0; j < nModes; j++)
148 {
149 for (label k = 0; k < nModes; k++)
150 {
151 ct1Tensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & fvc::laplacian(
152 nutModes[j], Umodes[k])).value();
153 }
154 }
155 }
156
157 if (Pstream::parRun())
158 {
159 reduce(ct1Tensor, sumOp<Eigen::Tensor<double, 3 >> ());
160 }
161
162 // Export the tensor
163 ITHACAstream::SaveDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/",
164 "ct1_" + name(nModes) + "_t");
165 return ct1Tensor;
166}
167
168Eigen::Tensor<double, 3> SteadyNSTurbIntrusive::turbulenceTensor2(label nModes)
169{
170 Eigen::Tensor<double, 3> ct2Tensor;
171 ct2Tensor.resize(nModes, nModes, nModes);
172
173 for (label i = 0; i < nModes; i++)
174 {
175 for (label j = 0; j < nModes; j++)
176 {
177 for (label k = 0; k < nModes; k++)
178 {
179 ct2Tensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & (fvc::div(
180 nutModes[j] * dev((fvc::grad(Umodes[k]))().T())))).value();
181 }
182 }
183 }
184
185 if (Pstream::parRun())
186 {
187 reduce(ct2Tensor, sumOp<Eigen::Tensor<double, 3 >> ());
188 }
189
190 // Export the tensor
191 ITHACAstream::SaveDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/",
192 "ct2_" + name(nModes) + "_t");
193 return ct2Tensor;
194}
195
196Eigen::MatrixXd SteadyNSTurbIntrusive::btTurbulence(label nModes)
197{
198 Eigen::MatrixXd btMatrix(nModes, nModes);
199 btMatrix = btMatrix * 0;
200
201 // Project everything
202 for (label i = 0; i < nModes; i++)
203 {
204 for (label j = 0; j < nModes; j++)
205 {
206 btMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & (fvc::div(dev((T(
207 fvc::grad(
208 Umodes[j]))))))).value();
209 }
210 }
211
212 // Export the matrix
213 ITHACAstream::SaveDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/",
214 "bt_" + name(nModes) + "_t");
215 return btMatrix;
216}
217
218void SteadyNSTurbIntrusive::project(fileName folder, label nModes)
219{
220 nModesOnline = nModes;
221 if (ITHACAutilities::check_folder("./ITHACAoutput/Matrices/"))
222 {
223 word bStr = "b_" + name(nModes);
224 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + bStr))
225 {
226 ITHACAstream::ReadDenseMatrix(bMatrix, "./ITHACAoutput/Matrices/", bStr);
227 }
228
229 else
230 {
231 bMatrix = diffusiveTerm(nModes);
232 }
233 word btStr = "bt_" + name(nModes);
234 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + btStr))
235 {
236 ITHACAstream::ReadDenseMatrix(btMatrix, "./ITHACAoutput/Matrices/", btStr);
237 }
238
239 else
240 {
241 btMatrix = btTurbulence(nModes);
242 }
243 word kStr = "k_" + name(nModes);
244 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + kStr))
245 {
246 ITHACAstream::ReadDenseMatrix(kMatrix, "./ITHACAoutput/Matrices/", kStr);
247 }
248
249 else
250 {
252 }
253 word cStr = "c_" + name(nModes) + "_t";
254 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + cStr))
255 {
256 ITHACAstream::ReadDenseTensor(convTensor, "./ITHACAoutput/Matrices/", cStr);
257 }
258
259 else
260 {
261 convTensor = convectiveTerm(nModes);
262 }
263 word ct1Str = "ct1_" + name(nModes) + "_t";
264 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct1Str))
265 {
266 ITHACAstream::ReadDenseTensor(ct1Tensor, "./ITHACAoutput/Matrices/", ct1Str);
267 }
268
269 else
270 {
272 }
273 word ct2Str = "ct2_" + name(nModes) + "_t";
274 if (ITHACAutilities::check_file("./ITHACAoutput/Matrices/" + ct2Str))
275 {
276 ITHACAstream::ReadDenseTensor(ct2Tensor, "./ITHACAoutput/Matrices/", ct2Str);
277 }
278
279 else
280 {
282 }
283 if (bcMethod == "penalty")
284 {
285 bcVelVec = bcVelocityVec(nModes);
286 bcVelMat = bcVelocityMat(nModes);
287 }
288 }
289
290 else
291 {
292 bMatrix = diffusiveTerm(nModes);
293 convTensor = convectiveTerm(nModes);
295 btMatrix = btTurbulence(nModes);
298
299 if (bcMethod == "penalty")
300 {
301 bcVelVec = bcVelocityVec(nModes);
302 bcVelMat = bcVelocityMat(nModes);
303 }
304 }
305
306 // Export the matrices
307 if (para->exportPython)
308 {
309 ITHACAstream::exportMatrix(bMatrix, "b", "python", "./ITHACAoutput/Matrices/");
310 ITHACAstream::exportMatrix(kMatrix, "k", "python", "./ITHACAoutput/Matrices/");
311 ITHACAstream::exportMatrix(btMatrix, "bt", "python",
312 "./ITHACAoutput/Matrices/");
314 "./ITHACAoutput/Matrices/");
315 ITHACAstream::exportTensor(ct1Tensor, "ct1", "python",
316 "./ITHACAoutput/Matrices/");
317 ITHACAstream::exportTensor(ct2Tensor, "ct2", "python",
318 "./ITHACAoutput/Matrices/");
319 }
320
321 if (para->exportMatlab)
322 {
323 ITHACAstream::exportMatrix(bMatrix, "b", "matlab", "./ITHACAoutput/Matrices/");
324 ITHACAstream::exportMatrix(kMatrix, "k", "matlab", "./ITHACAoutput/Matrices/");
325 ITHACAstream::exportMatrix(btMatrix, "bt", "matlab",
326 "./ITHACAoutput/Matrices/");
328 "./ITHACAoutput/Matrices/");
329 ITHACAstream::exportTensor(ct1Tensor, "ct1", "matlab",
330 "./ITHACAoutput/Matrices/");
331 ITHACAstream::exportTensor(ct2Tensor, "ct2", "matlab",
332 "./ITHACAoutput/Matrices/");
333 }
334
335 if (para->exportTxt)
336 {
337 ITHACAstream::exportMatrix(bMatrix, "b", "eigen", "./ITHACAoutput/Matrices/");
338 ITHACAstream::exportMatrix(kMatrix, "k", "eigen", "./ITHACAoutput/Matrices/");
339 ITHACAstream::exportMatrix(btMatrix, "bt", "eigen", "./ITHACAoutput/Matrices/");
341 "./ITHACAoutput/Matrices/c");
342 ITHACAstream::exportTensor(ct1Tensor, "ct1_", "eigen",
343 "./ITHACAoutput/Matrices/ct1");
344 ITHACAstream::exportTensor(ct2Tensor, "ct2_", "eigen",
345 "./ITHACAoutput/Matrices/ct2");
346 }
347
349 cTotalTensor.resize(nModes, nModes, nModes);
351}
352
353Eigen::MatrixXd SteadyNSTurbIntrusive::diffusiveTerm(label nModes)
354{
355 Eigen::MatrixXd bMatrix;
356 bMatrix.resize(nModes, nModes);
357
358 // Project everything
359 for (label i = 0; i < nModes; i++)
360 {
361 for (label j = 0; j < nModes; j++)
362 {
363 bMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & fvc::laplacian(
364 dimensionedScalar("1", dimless, 1), Umodes[j])).value();
365 }
366 }
367
368 if (Pstream::parRun())
369 {
370 reduce(bMatrix, sumOp<Eigen::MatrixXd>());
371 }
372
373 ITHACAstream::SaveDenseMatrix(bMatrix, "./ITHACAoutput/Matrices/",
374 "b_" + name(nModes));
375 return bMatrix;
376}
377
378Eigen::Tensor<double, 3> SteadyNSTurbIntrusive::convectiveTerm(label nModes)
379{
380 Eigen::Tensor<double, 3> convTensor;
381 convTensor.resize(nModes, nModes, nModes);
382
383 for (label i = 0; i < nModes; i++)
384 {
385 for (label j = 0; j < nModes; j++)
386 {
387 for (label k = 0; k < nModes; k++)
388 {
389 convTensor(i, j, k) = fvc::domainIntegrate(Umodes[i] & fvc::div(
390 linearInterpolate(Umodes[j]) & Umodes[j].mesh().Sf(),
391 Umodes[k])).value();
392 }
393 }
394 }
395
396 if (Pstream::parRun())
397 {
398 reduce(convTensor, sumOp<Eigen::Tensor<double, 3 >> ());
399 }
400
401 // Export the tensor
402 ITHACAstream::SaveDenseTensor(convTensor, "./ITHACAoutput/Matrices/",
403 "c_" + name(nModes) + "_t");
404 return convTensor;
405}
406
408{
409 Eigen::MatrixXd kMatrix(nModes, nModes);
410
411 // Project everything
412 for (label i = 0; i < nModes; i++)
413 {
414 for (label j = 0; j < nModes; j++)
415 {
416 kMatrix(i, j) = fvc::domainIntegrate(Umodes[i] & fvc::grad(
417 Pmodes[j])).value();
418 }
419 }
420
421 if (Pstream::parRun())
422 {
423 reduce(kMatrix, sumOp<Eigen::MatrixXd>());
424 }
425
426 ITHACAstream::SaveDenseMatrix(kMatrix, "./ITHACAoutput/Matrices/",
427 "k_" + name(nModes));
428 return kMatrix;
429}
430
431List< Eigen::MatrixXd > SteadyNSTurbIntrusive::bcVelocityVec(label nModes)
432{
433 List < Eigen::MatrixXd > bcVelVec(inletIndex.rows());
434
435 for (label j = 0; j < inletIndex.rows(); j++)
436 {
437 bcVelVec[j].resize(nModes, 1);
438 }
439
440 for (label k = 0; k < inletIndex.rows(); k++)
441 {
442 label BCind = inletIndex(k, 0);
443 label BCcomp = inletIndex(k, 1);
444
445 for (label i = 0; i < nModes; i++)
446 {
447 bcVelVec[k](i, 0) = gSum(Umodes[i].boundaryField()[BCind].component(
448 BCcomp));
449 }
450 }
451
452 ITHACAstream::exportMatrix(bcVelVec, "bcVelVec", "eigen",
453 "./ITHACAoutput/Matrices/bcVelVec");
454 return bcVelVec;
455}
456
457List< Eigen::MatrixXd > SteadyNSTurbIntrusive::bcVelocityMat(label nModes)
458{
459 label BCUsize = inletIndex.rows();
460 List < Eigen::MatrixXd > bcVelMat(BCUsize);
461
462 for (label j = 0; j < inletIndex.rows(); j++)
463 {
464 bcVelMat[j].resize(nModes, nModes);
465 }
466
467 for (label k = 0; k < inletIndex.rows(); k++)
468 {
469 label BCind = inletIndex(k, 0);
470 label BCcomp = inletIndex(k, 1);
471
472 for (label i = 0; i < nModes; i++)
473 {
474 for (label j = 0; j < nModes; j++)
475 {
476 bcVelMat[k](i, j) = gSum(Umodes[i].boundaryField()[BCind].component(
477 BCcomp) *
478 Umodes[j].boundaryField()[BCind].component(BCcomp));
479 }
480 }
481 }
482
483 ITHACAstream::exportMatrix(bcVelMat, "bcVelMat", "eigen",
484 "./ITHACAoutput/Matrices/bcVelMat");
485 return bcVelMat;
486}
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
#define M_Assert(Expr, Msg)
fv::options & fvOptions
Definition NLsolve.H:25
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:227
scalar maxIter
Number of maximum iterations to be done for the computation of the truth solution.
Definition steadyNS.H:128
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:265
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:302
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:293
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:296
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:299
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:125
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:290
autoPtr< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model)
Definition steadyNS.H:311
ITHACAparameters * para
Definition steadyNS.H:82
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:314
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition steadyNS.H:308
List< Eigen::MatrixXd > bcVelMat
Boundary term for penalty method - matrix.
Definition steadyNS.H:230
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:272
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
word bcMethod
Boundary Method.
Definition steadyNS.H:326
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:269
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.