All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
twoPhaseConstTau.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 \*---------------------------------------------------------------------------*/
29 
30 #include "twoPhaseConstTau.H"
31 #include "QGDThermo.H"
32 #include "twoPhaseIcoQGDThermo.H"
33 #include "addToRunTimeSelectionTable.H"
34 #include "linear.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace qgd
41 {
42  defineTypeNameAndDebug(twoPhaseConstTau,0);
44  (
45  QGDCoeffs,
46  twoPhaseConstTau,
47  dictionary
48  );
49 }
50 }
51 
54 (
55  const IOobject& io,
56  const fvMesh& mesh,
57  const dictionary& dict
58 )
59 :
60  QGDCoeffs(io, mesh, dict)
61 {
62  scalar ScQGD = 0.0, PrQGD = 1.0;
63 //
64  ScQGD_.primitiveFieldRef() = ScQGD;
65  PrQGD_.primitiveFieldRef() = PrQGD;
66  muQGD_.primitiveFieldRef() = 0.0;
67  alphauQGD_.primitiveFieldRef() = 0.0;
68 
69  ScQGD_.boundaryFieldRef() = ScQGD;
70  PrQGD_.boundaryFieldRef() = PrQGD;
71  muQGD_.boundaryFieldRef() = 0.0;
72  alphauQGD_.boundaryFieldRef() = 0.0;
73 }
74 
77 {
78 }
79 
80 void Foam::qgd::
81 twoPhaseConstTau::correct(const QGDThermo& qgdThermo)
82 {
83  const volScalarField nu = qgdThermo.mu()/qgdThermo.rho();
84  if (isA<twoPhaseIcoQGDThermo>(qgdThermo))
85  {
86  const twoPhaseIcoQGDThermo& refThermo =
87  refCast<const twoPhaseIcoQGDThermo>(qgdThermo);
88 
89  const dimensionedScalar& Tau1 = refThermo.Tau1();
90  const dimensionedScalar& Tau2 = refThermo.Tau2();
91  const volScalarField& alpha1 = refThermo.alpha1();
92  this->tauQGD_ = alpha1*Tau1 + (1.0 - alpha1)*Tau2;
93  Info << "max/min tauQGD: " << max(tauQGD_).value() << "/" << min(tauQGD_).value() << endl;
94  }
95  else
96  {
97  FatalErrorInFunction
98  << "twoPhaseConstTau::correct(const QGDThermo& qgdThermo)"
99  << "MUST be only with two-phase solver interQHDFoam"
100  << exit(FatalError);
101  }
102 
103  this->tauQGDf_= linearInterpolate(this->tauQGD_);
104 }
105 
106 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
107 
108 
109 // ************************************************************************* //
volScalarField PrQGD_
Definition: QGDCoeffs.H:115
void correct(const QGDThermo &)
Thermodynamics and mechanics class for incompressible two-phase mixture of immiscible components...
const dimensionedScalar & Tau1() const
volScalarField alphauQGD_
Definition: QGDCoeffs.H:80
virtual tmp< volScalarField > rho() const =0
const dimensionedScalar & Tau2() const
virtual tmp< volScalarField > mu() const =0
Abstract base class for classes implementing thermophysical properties of gases and fluids governed b...
Definition: QGDThermo.H:53
volScalarField muQGD_
Definition: QGDCoeffs.H:75
twoPhaseConstTau(const IOobject &io, const fvMesh &mesh, const dictionary &dict)
addToRunTimeSelectionTable(QGDCoeffs, constScPrModel1, dictionary)
defineTypeNameAndDebug(constScPrModel1, 0)
volScalarField ScQGD_
Definition: QGDCoeffs.H:120