All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
twoPhaseIcoQGDThermo.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 
13 License
14  This file is part of QGDsolver, based on OpenFOAM library.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU 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  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 Group
30  grpTwoPhaseIcoQGDThermo
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "twoPhaseIcoQGDThermo.H"
35 #include "QGDCoeffs.H"
36 #include "twoPhaseConstTau.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 Foam::twoPhaseIcoQGDThermo::twoPhaseIcoQGDThermo(const fvMesh& mesh, const volVectorField& U)
49 :
50  IOdictionary
51  (
52  IOobject
53  (
54  "thermophysicalProperties",
55  mesh.time().constant(),
56  mesh.time().db(),
57  IOobject::MUST_READ_IF_MODIFIED,
58  IOobject::NO_WRITE
59  )
60  ),
61  twoPhaseMixture(mesh, *this),
62  constTwoPhaseProperties(*this, phase1Name(), phase2Name()),
63  qInterfaceProperties(alpha1(),U,*this),
64  QGDThermo(mesh, *this),
65  p_
66  (
67  IOobject
68  (
69  "p",
70  mesh.time().timeName(),
71  mesh,
72  IOobject::MUST_READ,
73  IOobject::AUTO_WRITE
74  ),
75  mesh
76  )
77 {
78  validateQGDCoeffs();
79  this->read();
80 }
81 
82 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
83 
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  if (!isA<qgd::twoPhaseConstTau>(this->qgdCoeffs()))
93  {
94  FatalErrorInFunction
95  << "twoPhaseIcoQGDThermo works only with twoPhaseConstTau QGD Coeffs Model"
96  << exit(FatalError);
97  }
98 }
99 
101 {
102 
103  if (!QGDThermo::read())
104  {
105  return false;
106  }
107 
108  return true;
109 }
110 
111 const Foam::volScalarField& Foam::twoPhaseIcoQGDThermo::p() const
112 {
113  return p_;
114 }
115 
116 Foam::volScalarField& Foam::twoPhaseIcoQGDThermo::p()
117 {
118  return p_;
119 }
121 const Foam::volScalarField& Foam::twoPhaseIcoQGDThermo::c() const
122 {
123  notImplemented ("const volScalarField& Foam::twoPhaseIcoQGDThermo::c() const");
124  return volScalarField::null();
125 }
126 
127 Foam::tmp<Foam::volScalarField> Foam::twoPhaseIcoQGDThermo::rho() const
128 {
129  return (rho1()-rho2())*this->alpha1() + rho2();
130 }
131 
132 Foam::tmp<Foam::volScalarField> Foam::twoPhaseIcoQGDThermo::mu() const
133 {
134  return (rho1()*nu1()-rho2()*nu2())*this->alpha1() + rho2()*nu2();
135 }
136 
138 {
139  qgdCoeffs().correct(*this);
141 }
142 
143 
144 
145 // ************************************************************************* //
virtual tmp< volScalarField > mu() const
Thermodynamics and mechanics class for incompressible two-phase mixture of immiscible components...
virtual bool read()
Definition: QGDThermo.C:65
virtual tmp< volScalarField > rho() const
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
virtual const volScalarField & c() const
Deprecated for incompressible models.
virtual const volScalarField & p() const
virtual ~twoPhaseIcoQGDThermo()
Destructor.
Abstract base class for classes implementing thermophysical properties of gases and fluids governed b...
Definition: QGDThermo.H:53
Contains the interface properties.
defineTypeNameAndDebug(psiQGDReactionThermo, 0)
twoPhaseIcoQGDThermo(const twoPhaseIcoQGDThermo &)
Construct as copy (not implemented)