Loading...
Searching...
No Matches
CompressibleUnSteadyRhoPimple.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
13 License
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
34
36//#include "viscosityModel.H"
37
38// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
39// Constructor
42 char* argv[])
43 : UnsteadyProblem()
44{
45 _args = autoPtr<argList>
46 (
47 new argList(argc, argv)
48 );
49
50 if (!_args->checkRootCase())
51 {
52 Foam::FatalError.exit();
53 }
54
55 argList& args = _args();
56#include "createTime.H"
57 Info << "Create a dynamic mesh for time = "
58 << runTime.timeName() << nl << endl;
59 meshPtr = autoPtr<dynamicFvMesh> (dynamicFvMesh::New(args, runTime));
60 dynamicFvMesh& mesh = meshPtr();
61 _pimple = autoPtr<pimpleControl>
62 (
63 new pimpleControl
64 (
65 mesh
66 )
67 );
68 pimpleControl& pimple = _pimple();
69#include "createFields.H"
70 //#include "createFvOptions.H"
71#include "initContinuityErrs.H"
72 turbulence->validate();
73 ITHACAdict = new IOdictionary
74 (
75 IOobject
76 (
77 "ITHACAdict",
78 runTime.system(),
79 mesh,
80 IOobject::MUST_READ,
81 IOobject::NO_WRITE
82 )
83 );
87 //middleExport = para->ITHACAdict->lookupOrDefault<bool>("middleExport", true);
88 //middleStep = para->ITHACAdict->lookupOrDefault<label>("middleStep", 20);
89}
90
91
92// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
93
96{
97 Time& runTime = _runTime();
98 volVectorField& U = _U();
99 volScalarField& p = pThermo().p();
100 volScalarField& E = pThermo().he();
101 pimpleControl& pimple = _pimple();
102 volScalarField _nut(turbulence->nut());
103 instantList Times = runTime.times();
104 runTime.setEndTime(finalTime);
105 runTime.setTime(Times[1], 1);
106 runTime.setDeltaT(timeStep);
107 nextWrite = startTime; // timeStep initialization
108 int counter = 1;
109#include "NLsolve.H"
110}
111
113{
114 const volScalarField& mu = pThermo().mu();
115 volScalarField& mu_field = const_cast<volScalarField&>(mu);
116 this->assignIF(mu_field, mu_new);
117
118 for (int i = 0; i < mu_field.boundaryFieldRef().size(); i++)
119 {
120 this->assignBC(mu_field, i, mu_new);
121 }
122}
123
124fvVectorMatrix CompressibleUnSteadyRhoPimple::getNLTerm(volVectorField& U)
125{
126 surfaceScalarField& phi = _phi();
127 fvVectorMatrix NLTerm = fvm::div(phi, U);
128 return NLTerm;
129}
130
131fvVectorMatrix CompressibleUnSteadyRhoPimple::getViscTerm(volVectorField& U)
132{
133 volScalarField& rho = _rho();
134 volScalarField nuEff = turbulence->nuEff();
135 //fvVectorMatrix viscTerm = turbulence->divDevRhoReff(U);
136 fvVectorMatrix viscTerm = - fvc::div((rho * nuEff) * dev2(T(fvc::grad(
137 U)))) - fvm::laplacian(rho * nuEff, U);
138 return viscTerm;
139}
140
141volVectorField CompressibleUnSteadyRhoPimple::getGradP(volScalarField& p)
142{
143 volVectorField gradP = fvc::grad(p);
144 return gradP;
145}
146
147
148fvVectorMatrix CompressibleUnSteadyRhoPimple::getUmatrix(volVectorField& U)
149{
150 volScalarField& rho = _rho();
151 fv::options& fvOptions = _fvOptions();
152 Ueqn_global.reset(
153 new fvVectorMatrix
154 ( fvm::ddt(rho, U)
155 +
156 getNLTerm(U)
157 +
158 getViscTerm(U)
159 ==
160 fvOptions(rho, U)
161 )
162 );
163 Ueqn_global().relax();
164 fvOptions.constrain(Ueqn_global());
165 return Ueqn_global();
166}
167
168fvScalarMatrix CompressibleUnSteadyRhoPimple::getFluxTerm()
169{
170 volScalarField& he = pThermo().he();
171 surfaceScalarField& phi = _phi();
172 fvScalarMatrix fluxTerm = fvm::div(phi, he);
173 return fluxTerm;
174}
175
176volScalarField CompressibleUnSteadyRhoPimple::getKinEnTerm(volVectorField& U,
177 volScalarField& p)
178{
179 surfaceScalarField& phi = _phi();
180 volScalarField& rho = _rho();
181 volScalarField kinEn = fvc::div(phi, volScalarField("Ekp",
182 0.5 * magSqr(U) + p / rho));
183 return kinEn;
184}
185
186fvScalarMatrix CompressibleUnSteadyRhoPimple::getDiffTerm()
187{
188 volScalarField& he = pThermo().he();
189 fvScalarMatrix diffTerm = fvm::laplacian(turbulence->alphaEff(), he);
190 return diffTerm;
191}
192
193fvScalarMatrix CompressibleUnSteadyRhoPimple::getEmatrix(volVectorField& U,
194 volScalarField& p)
195{
196 volScalarField& he = pThermo().he();
197 volScalarField& rho = _rho();
198 fv::options& fvOptions = _fvOptions();
199 Eeqn_global.reset(new fvScalarMatrix(
200 getFluxTerm() + getKinEnTerm(U, p) - getDiffTerm()
201 ==
202 fvOptions(rho, he)
203 ));
204 Eeqn_global().relax();
205 fvOptions.constrain(Eeqn_global());
206 return Eeqn_global();
207}
208
209surfaceScalarField CompressibleUnSteadyRhoPimple::getPhiHbyA(
210 fvVectorMatrix& Ueqn,
211 volVectorField& U, volScalarField& p)
212{
213 volScalarField& rho = _rho();
214 volScalarField rAU(1.0 /
215 Ueqn.A()); // Inverse of the diagonal part of the U equation matrix
216 HbyA.reset(new volVectorField(constrainHbyA(rAU * Ueqn.H(), U,
217 p))); // H is the extra diagonal part summed to the r.h.s. of the U equation
218 phiHbyA.reset(new surfaceScalarField("phiHbyA",
219 fvc::interpolate(rho)*fvc::flux(HbyA.ref())));
220 return phiHbyA;
221}
222
223volScalarField CompressibleUnSteadyRhoPimple::getDivPhiHbyA(
224 fvVectorMatrix& Ueqn,
225 volVectorField& U, volScalarField& p)
226{
227 volScalarField divPhiHbyA = fvc::div(getPhiHbyA(Ueqn, U, p));
228 return divPhiHbyA;
229}
230
231surfaceScalarField CompressibleUnSteadyRhoPimple::getRhorAUf(
232 fvVectorMatrix& Ueqn)
233{
234 volScalarField& rho = _rho();
235 volScalarField rAU(1.0 /
236 Ueqn.A()); // Inverse of the diagonal part of the U equation matrix
237 rhorAUf.reset(new surfaceScalarField("rhorAUf", fvc::interpolate(rho * rAU)));
238 return rhorAUf;
239}
240
241fvScalarMatrix CompressibleUnSteadyRhoPimple::getPoissonTerm(
242 fvVectorMatrix& Ueqn, volScalarField& p)
243{
244 fvScalarMatrix poissonTerm = fvm::laplacian(getRhorAUf(Ueqn), p);
245 return poissonTerm;
246}
247
248fvScalarMatrix CompressibleUnSteadyRhoPimple::getPmatrix(fvVectorMatrix& Ueqn,
249 volVectorField& U, volScalarField& p)
250{
251 volScalarField& rho = _rho();
252 fv::options& fvOptions = _fvOptions();
253 volScalarField& psi = _psi();
254 pressureControl& pressureControl = _pressureControl();
255 Peqn_global.reset(new fvScalarMatrix(
256 getDivPhiHbyA(Ueqn, U, p)
257 - getPoissonTerm(Ueqn, p)
258 ==
259 fvOptions(psi, p, rho.name())
260 ));
261 Peqn_global().setReference
262 (
263 pressureControl.refCell(),
264 pressureControl.refValue()
265 );
266 return Peqn_global();
267}
268
270{
271 _runTime().objectRegistry::clear();
272 _mesh().objectRegistry::clear();
273 // _mesh.clear();
274 // _runTime.clear();
275 _simple.clear();
276 pThermo.clear();
277 _p.clear();
278 _rho.clear();
279 _E.clear();
280 _U.clear();
281 _phi.clear();
282 turbulence.clear();
283 _initialMass.clear();
284 _pressureControl.clear();
285 _psi.clear();
286 _fvOptions.clear();
287 argList& args = _args();
288 Info << "ReCreate time\n" << Foam::endl;
289 Time& runTime = _runTime();
290 runTime.setTime(0, 1);
291 Foam::fvMesh& mesh = _mesh();
292 _simple = autoPtr<simpleControl>
293 (
294 new simpleControl
295 (
296 mesh
297 )
298 );
299 simpleControl& simple = _simple();
300 pThermo = autoPtr<fluidThermo>
301 (
302 fluidThermo::New(mesh)
303 );
304 fluidThermo& thermo = pThermo();
305 thermo.validate(_args().executable(), "h", "e");
306 _p = autoPtr<volScalarField>
307 (
308 new volScalarField(thermo.p())
309 );
310 volScalarField& p = thermo.p();
311 _rho = autoPtr<volScalarField>
312 (
313 new volScalarField
314 (
315 IOobject
316 (
317 "rho",
318 runTime.timeName(),
319 mesh,
320 IOobject::READ_IF_PRESENT,
321 IOobject::AUTO_WRITE
322 ),
323 thermo.rho()
324 )
325 );
326 volScalarField& rho = _rho();
327 _E = autoPtr<volScalarField>
328 (
329 new volScalarField(thermo.he())
330 );
331 Info << "ReReading field U\n" << endl;
332 _U = autoPtr<volVectorField>
333 (
334 new volVectorField
335 (
336 IOobject
337 (
338 "U",
339 runTime.timeName(),
340 mesh,
341 IOobject::MUST_READ,
342 IOobject::AUTO_WRITE
343 ),
344 mesh
345 )
346 );
347 volVectorField& U = _U();
348 Info << "ReReading/calculating face flux field phi\n" << endl;
349 _phi = autoPtr<surfaceScalarField>
350 (
351 new surfaceScalarField
352 (
353 IOobject
354 (
355 "phi",
356 runTime.timeName(),
357 mesh,
358 IOobject::READ_IF_PRESENT,
359 IOobject::AUTO_WRITE
360 ),
361 linearInterpolate(rho * U) & mesh.Sf()
362 )
363 );
364 surfaceScalarField& phi = _phi();
365 _pressureControl = autoPtr<pressureControl>
366 (
367 new pressureControl(p, rho, _simple().dict())
368 );
369 mesh.setFluxRequired(p.name());
370 Info << "ReCreating turbulence model\n" << endl;
371 turbulence = autoPtr<compressible::turbulenceModel>
372 (
373 compressible::turbulenceModel::New
374 (
375 rho,
376 U,
377 phi,
378 thermo
379 )
380 );
381 _initialMass = autoPtr<dimensionedScalar>
382 (
383 new dimensionedScalar(fvc::domainIntegrate(rho))
384 );
385 _psi = autoPtr<volScalarField>
386 (
387 new volScalarField(thermo.psi())
388 );
389 _fvOptions = autoPtr<fv::options>(new fv::options(mesh));
390#include "initContinuityErrs.H"
391 //turbulence->validate();
392}
Header file of the CompressibleUnSteadyPimple class.
void changeViscosity(double mu_new)
Function to change the viscosity.
void restart()
set all variables back to the values into the 0 folder
autoPtr< pimpleControl > _pimple
pimpleControl
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
autoPtr< volScalarField > _nut
Eddy viscosity field.
scalar startTime
Start Time (initial time to start storing the snapshots).
scalar timeStep
Time step of the simulation.
scalar nextWrite
Auxiliary variable to store the next writing instant.
scalar finalTime
Final time (final time of the simulation and consequently of the acquisition of the snapshots).
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
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:302
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:293
autoPtr< fv::options > _fvOptions
fvOptions
Definition steadyNS.H:296
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:299
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:290
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:272
bool check_pod()
Check if the POD data folder "./ITHACAoutput/POD" exists.
bool check_off()
Check if the offline data folder "./ITHACAoutput/Offline" exists.