Loading...
Searching...
No Matches
CompressibleUnSteadyRhoPimple.H
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-------------------------------------------------------------------------------
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/>.
24Class
25 CompressibleUnSteadyRhoPimple
26Description
27 full-order class for a Compressible UnSteady problem with Pimple algorithm.
28SourceFiles
29 CompressibleUnSteadyRhoPimple.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef CompressibleUnSteadyRhoPimple_H
38#define CompressibleUnSteadyRhoPimple_H
39#include "fvCFD.H"
40#include "dynamicFvMesh.H"
41#include "fluidThermo.H"
42#include "turbulentFluidThermoModel.H"
44#include "CorrectPhi.H"
45#include "UnsteadyProblem.H"
46#include "pimpleControl.H"
47#include "pressureControl.H"
48#include "fvOptions.H"
49#include "reductionProblem.H"
50#include "ITHACAstream.H"
51#include "ITHACAparameters.H"
52#include "ITHACAforces.H"
53#include <iostream>
54
55
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58/*---------------------------------------------------------------------------*\
59 Class CompressibleUnSteadyRhoPimple Declaration
60\*---------------------------------------------------------------------------*/
61
63
67 public UnsteadyProblem
68{
69 public:
70 // Constructors
73
75 CompressibleUnSteadyRhoPimple(int argc, char* argv[]);
77
78 // Functions
79
80 //--------------------------------------------------------------------------
86 //void truthSolve(List<scalar> mu_now);
87 void truthSolve(word folder);
88
89 // Momuntum equation terms
90 fvVectorMatrix getUmatrix(volVectorField& U);//, Vector<double>& uresidual_v);
91 fvVectorMatrix getNLTerm(volVectorField& U);
92 fvVectorMatrix getViscTerm(volVectorField& U);
93 volVectorField getGradP(volScalarField& p);
94
95 // Energy equation terms
96 fvScalarMatrix getEmatrix(volVectorField& U,
97 volScalarField& p);//, scalar& eresidual);
98 fvScalarMatrix getFluxTerm();
99 volScalarField getKinEnTerm(volVectorField& U, volScalarField& p);
100 fvScalarMatrix getDiffTerm();
101
102 // Pressure equation terms
103 fvScalarMatrix getPmatrix(fvVectorMatrix& Ueqn, volVectorField& U,
104 volScalarField& p);
105 surfaceScalarField getPhiHbyA(fvVectorMatrix& Ueqn, volVectorField& U,
106 volScalarField& p);
107 volScalarField getDivPhiHbyA(fvVectorMatrix& Ueqn, volVectorField& U,
108 volScalarField& p);
109 surfaceScalarField getRhorAUf(fvVectorMatrix& Ueqn);
110 fvScalarMatrix getPoissonTerm(fvVectorMatrix& Ueqn, volScalarField& p);
111
112 //--------------------------------------------------------------------------
117 void changeViscosity(double mu_new);
118
119 //--------------------------------------------------------------------------
122 void restart();
123
124 autoPtr<fluidThermo> pThermo;
125 // /// Dynamic mesh field
126 autoPtr<Foam::dynamicFvMesh> meshPtr;
128 autoPtr<pimpleControl> _pimple;
129
130 autoPtr<volScalarField> _rho;
131
132 autoPtr<pressureControl> _pressureControl;
133
134 autoPtr<compressible::turbulenceModel> turbulence;
135
136 autoPtr<dimensionedScalar> _initialMass;
137
138 autoPtr<volScalarField> _psi;
139
140 autoPtr<volScalarField> _E;
141 //autoPtr<surfaceVectorField> _rhoUf;
142
144 autoPtr<volScalarField> _p0;
145
147 autoPtr<volScalarField> _E0;
148
150 autoPtr<volVectorField> _U0;
151
153 autoPtr<volScalarField> _rho0;
154
156 autoPtr<surfaceScalarField> _phi0, _phi;
158 PtrList<volScalarField> Efield;
159
161 volScalarModes Emodes;
162
163 autoPtr<fvVectorMatrix> Ueqn_global;
164
165 autoPtr<fvScalarMatrix> Peqn_global;
166
167 autoPtr<fvScalarMatrix> Eeqn_global;
168
169 autoPtr<surfaceScalarField> phiHbyA;
170
171 autoPtr<surfaceScalarField> rhorAUf;
172
173 autoPtr<volVectorField> HbyA;
174
177 //label middleStepInit;
178
181
182 bool correctPhi;
183 bool checkMeshCourantNo;
184 bool moveMeshOuterCorrectors;
185 //const dimensionedScalar rhoMax;
186 //const dimensionedScalar rhoMin;
187
188
189
190};
191
192#endif
193
194
195
196
197
198
199
200
201
202
Header file of the steadyNS class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
CompressibleSteadyNS()
Null constructor.
autoPtr< volVectorField > _U0
Initial Velocity field (for restart purposes).
autoPtr< volScalarField > _p0
Initial Pressure field (for restart purposes).
PtrList< volScalarField > Efield
List of pointers used to store the energy solutions.
volScalarModes Emodes
List of pointers used to form the energy modes.
autoPtr< volScalarField > _E0
Initial Energy field (for restart purposes).
autoPtr< surfaceScalarField > _phi0
Initial Flux field (for restart purposes).
void changeViscosity(double mu_new)
Function to change the viscosity.
void restart()
set all variables back to the values into the 0 folder
bool middleExport
Export also intermediate fields.
label middleStep
Distancing between intermediate steps (for turbulent case only).
autoPtr< pimpleControl > _pimple
pimpleControl
autoPtr< volScalarField > _rho0
Initial Density field (for restart purposes).
void truthSolve()
Perform a TruthSolve.
Header file of the reductionProblem class.