All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
interQHDFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of QGDsolver library, based on OpenFOAM+.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Application
29  QHDFoam
30 
31 Description
32  Solver for unsteady 3D turbulent flow of 2 incompressible immisicible
33  fluids governed by quasi-hydrodynamic dynamic (QHD) equations.
34 
35  QHD system of equations has been developed by scientific group from
36  Keldysh Institute of Applied Mathematics,
37  see http://elizarova.imamod.ru/selection-of-papers.html
38 
39  A comprehensive description of QGD equations and their applications
40  can be found here:
41  \verbatim
42  Elizarova, T.G.
43  "Quasi-Gas Dynamic equations"
44  Springer, 2009
45  \endverbatim
46 
47  A brief of theory on QGD and QHD system of equations:
48  \verbatim
49  Elizarova, T.G. and Sheretov, Y.V.
50  "Theoretical and numerical analysis of quasi-gasdynamic and quasi-hydrodynamic
51  equations"
52  J. Computational Mathematics and Mathematical Physics, vol. 41, no. 2, pp 219-234,
53  2001
54  \endverbatim
55 
56  Developed by UniCFD group (www.unicfd.ru) of ISP RAS (www.ispras.ru).
57 
58 
59 \*---------------------------------------------------------------------------*/
60 
61 #include "fvCFD.H"
62 #include "wallFvPatch.H"
63 #include "fvsc.H"
64 #include "twoPhaseIcoQGDThermo.H"
65 #include "QGDInterpolate.H"
66 #include "MULES.H"
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 int main(int argc, char *argv[])
71 {
72  #define NO_CONTROL
73  #include "postProcess.H"
74 
75  #include "setRootCase.H"
76  #include "createTime.H"
77  #include "createMesh.H"
78  #include "createFields.H"
79  #include "createFaceFields.H"
80  #include "createFaceFluxes.H"
81  #include "createTimeControls.H"
82 
83  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85  // Courant numbers used to adjust the time-step
86  scalar CoNum = 0.0;
87  scalar meanCoNum = 0.0;
88 
89  Info<< "\nStarting time loop\n" << endl;
90 
91  while (runTime.run())
92  {
93  /*
94  *
95  * Update physical properties
96  *
97  */
98  thermo.correct();
99  //turbulence.correct();
100 
101  /*
102  *
103  * Update fields
104  *
105  */
106  #include "updateFields.H"
107 
108  /*
109  *
110  * Update fluxes
111  *
112  */
113  #include "updateFluxes.H"
114  /*
115  *
116  * Update time step
117  *
118  */
119  #include "readTimeControls.H"
120  #include "QHDCourantNo.H"
121  #include "setDeltaT-QGDQHD.H"
122 
123  runTime++;
124  Info<< "Time = " << runTime.timeName() << nl << endl;
125 
126  // --- Store old time values
127  U.oldTime();
128  alpha1.oldTime();
129  rho.oldTime();
130 
131  //Continuity equation
133  surfaceScalarField tphi = phiu*da1dtf*(Tau1-Tau2);
134  tphi.setOriented(true);
135  phiwm += tphi;
136  coeffp = alpha1f*Tau1/rho1 + alpha2f*Tau2/rho2;
137  p.correctBoundaryConditions();
138 
139  //solve for pressure
140  {
141  fvScalarMatrix pEqn1
142  (
143  -fvm::laplacian(Tau1/rho1,p)
144  );
145  fvScalarMatrix pEqn2
146  (
147  -fvm::laplacian(Tau2/rho2,p)
148  );
149  fvScalarMatrix pEqn
150  (
151  fvc::div(phiu)
152  +fvc::div(phiwm)
153  -fvm::laplacian(coeffp,p)
154  );
155  pEqn.setReference(pRefCell, getRefCellValue(p, pRefCell));
156  pEqn.solve();
157 
158  phiw1 = phiwo1 - pEqn1.flux();
159  phiw2 = phiwo2 - pEqn2.flux();
160  phi1 = phiu - phiw1;
161  phi2 = phiu - phiw2;
162 
163  phi = phiu+phiwm+pEqn.flux();
164  }
165 
166  gradpf = fvsc::grad(p);
167  W1 = ((Uf & gradUf) + (1./rho1)*gradpf - g - cFrcf/rho1)*Tau1;
168  W2 = ((Uf & gradUf) + (1./rho2)*gradpf - g - cFrcf/rho2)*Tau2;
169 
170  phiWr =
171  qgdFlux
172  (
173  -phiw2+phiw1,
174  alpha2,
175  alpha2f,
176  "div(phi,alpha1)"
177  );
178  phiAlpha1f =
179  qgdFlux
180  (
181  phi,
182  alpha1,
183  alpha1f,
184  "div(phi,alpha1)"
185  )
186  +
187  //add qgd terms from Wr
188  qgdFlux
189  (
190  -phiWr,
191  alpha1,
192  alpha1f,
193  "div(phi,alpha1)"
194  );
195 
196  //add terms from da1dt
197  {
198  surfaceScalarField DeltaTauFlux =
199  phiu*da1dtf*(Tau1 - alpha1f*(Tau1-Tau2));
200  //phiu*da1dtf*Tau1;
201  DeltaTauFlux.setOriented(true);
202  phiAlpha1f += DeltaTauFlux;
203  }
204 
205  //apply compressive fluxes and MULES limiter
206  if (thermo.cAlpha() > SMALL)
207  {
208  surfaceScalarField phic(thermo.cAlpha()*mag(phi/mesh.magSf()));
209 
210  // Do not compress interface at non-coupled boundary faces
211  // (inlets, outlets etc.)
212  {
213  surfaceScalarField::Boundary& phicBf =
214  phic.boundaryFieldRef();
215 
216  forAll(phic.boundaryField(), patchi)
217  {
218  fvsPatchScalarField& phicp = phicBf[patchi];
219 
220  if (!phicp.coupled())
221  {
222  phicp == 0;
223  }
224  }
225  }
226  surfaceScalarField phir(phic*thermo.nHatf());
227 
228  phiAlpha1f +=
229  fvc::flux
230  (
231  -fvc::flux(-phir, alpha2, "div(phir,alphar)"),
232  alpha1,
233  "div(phir,alphar)"
234  );
235 
236  //limit flux with MULES limiter
237  MULES::limit
238  (
239  1.0 / runTime.deltaTValue(),
240  geometricOneField(),
241  alpha1,
242  phi,
243  phiAlpha1f,
244  zeroField(),
245  zeroField(),
246  oneField(),
247  zeroField(),
248  false //return total flux
249  );
250  }
251 
252  // --- Solve for volume fraction
253  {
254 
255  solve
256  (
257  fvm::ddt(alpha1)
258  + fvc::div(phiAlpha1f)
259  );
260  //
261  Info << "max/min alpha1: " << max(alpha1).value() << "/" << min(alpha1).value() << endl;
262  alpha1 = max(min(alpha1,1.0),0.0);
263  alpha2 = 1.0 - alpha1;
264  }
265 
266  // --- Update density field
267  rho = thermo.rho();
268 
269  // --- Update mass fluxes
270  {
272  rhoPhi = phiAlpha1f*rho1 + phiAlpha2f*rho2;
273 
274  phiRhofWf = phiu*(alpha1f*rho1*W1 + alpha2f*rho2*W2);
275  phiRhofWf.setOriented(true);
277  (
278  rhoPhi,
279  U,
280  Uf
281  )
282  -
283  phiRhofWf;
284  }
285 
286  // --- Solve U
287  if (implicitDiffusion)
288  {
289  solve
290  (
291  fvm::ddt(rho,U)
292  +
294  +
295  fvc::grad(p)
296  *(1.0 + da1dt*(Tau1-Tau2))
297  -
298  fvm::laplacian(muf,U)
299  -
300  fvc::div(muf*mesh.Sf() & qgdInterpolate(Foam::T(fvc::grad(U))))
301  -
302  BdFrc
303  -
304  cFrc*(1.0 + da1dt*(Tau1-Tau2))
305  );
306  }
307  else
308  {
309  solve
310  (
311  fvm::ddt(rho,U)
312  +
314  +
315  fvc::grad(p)
316  *(1.0 + da1dt*(Tau1-Tau2))
317  -
318  fvc::laplacian(muf,U)
319  -
320  fvc::div(muf*mesh.Sf() & qgdInterpolate(Foam::T(fvc::grad(U))))
321  -
322  BdFrc
323  -
324  cFrc*(1.0 + da1dt*(Tau1-Tau2))
325  );
326  }
327 
328  runTime.write();
329 
330  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
331  << " ClockTime = " << runTime.elapsedClockTime() << " s"
332  << nl << endl;
333 
334  }
335 
336  Info<< "End\n" << endl;
337 
338  return 0;
339 }
340 
341 
342 // ************************************************************************* //
surfaceScalarField phiAlpha1f("phiAlpha1f", phi *alpha1f)
surfaceScalarField phi1("phi1", mesh.Sf()&linearInterpolate(U))
tmp< GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > > qgdFlux(const GeometricField< scalar, Foam::fvsPatchField, Foam::surfaceMesh > &flux, const GeometricField< T, Foam::fvPatchField, Foam::volMesh > &psi, const GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > &psif, const word fluxName)
Calculates the mean and maximum wave speed based Courant Numbers.
surfaceScalarField phiw2("phiw2", phi *alpha1f)
Creates interpolation instances templated for QGD solver.
surfaceScalarField phiw1("phiw1", phi *alpha1f)
volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho())
Creates the face fields: linear interpolation of fields from volumes to face centers.
muf
Definition: updateFields.H:38
int main(int argc, char *argv[])
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
surfaceVectorField phiRhofWf("phiUf", phi *Uf *rho1)
EEqn solve()
gradUf
Definition: updateFields.H:34
phiwo2
Definition: updateFluxes.H:43
surfaceScalarField phi2("phi2", mesh.Sf()&linearInterpolate(U))
phiwo1
Definition: updateFluxes.H:42
surfaceVectorField phiUfRhof("phiUf", phi *Uf *rho1)
tmp< GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > > qgdInterpolate(const GeometricField< T, Foam::fvPatchField, Foam::volMesh > &psi)
Reset the timestep to maintain a constant maximum courant Number. Reduction of time-step is immediate...
surfaceScalarField phiwm("phiwm", phiwo1 *0.0)
surfaceVectorField gradpf("gradpf", fvsc::grad(p))
forAll(Y, i)
Definition: QGDYEqn.H:36
BdFrc
Definition: updateFields.H:35
volScalarField & p
Definition: createFields.H:14
Switch implicitDiffusion(thermo.implicitDiffusion())
Creates the face-flux fields.
tmp< surfaceScalarField > div(const volVectorField &vF)
surfaceVectorField cFrcf("cFrcf",)
surfaceScalarField rhoPhi("rhoPhi", rho1 *phi)
Updates fields.
surfaceScalarField alpha1f("alpha1f",)
Uf
Definition: updateFields.H:30
surfaceVectorField W1("Wf", linearInterpolate(U)*0.0)
phiu
——–Start———
Definition: updateFluxes.H:33
surfaceScalarField alpha2f("alpha2f", 1.0-alpha1f)
surfaceScalarField phiWr("phiWr", phiwo1 *0.0)
psiQGDThermo & thermo
Definition: createFields.H:7
tmp< surfaceVectorField > grad(const volScalarField &vF)
fvScalarMatrix pEqn(fvc::div(phiu)-fvc::div(phiwo)-fvm::laplacian(taubyrhof, p))
surfaceVectorField W2("Wf", linearInterpolate(U)*0.0)
const volScalarField & T
Definition: createFields.H:15
surfaceScalarField phiAlpha2f("phiAlpha2f", phi *alpha1f)
Updates fluxes for continuity equation.
surfaceScalarField coeffp("coeffp",)