viscousHeating.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "viscousHeating.H"
27 #include "basicThermo.H"
29 #include "fvMatrices.H"
30 #include "fvcDiv.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
40 
42  (
43  fvModel,
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 void Foam::fv::viscousHeating::readCoeffs()
54 {
55  phaseName_ = coeffs().lookupOrDefault<word>("phase", word::null);
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
62 (
63  const word& name,
64  const word& modelType,
65  const fvMesh& mesh,
66  const dictionary& dict
67 )
68 :
69  fvModel(name, modelType, mesh, dict),
70  phaseName_(word::null)
71 {
72  readCoeffs();
73 }
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 {
80  const basicThermo& thermo =
81  mesh().lookupObject<basicThermo>
82  (
83  IOobject::groupName(physicalProperties::typeName, phaseName_)
84  );
85 
86  return wordList(1, thermo.he().name());
87 }
88 
89 
91 (
92  const volScalarField& rho,
93  const volScalarField& he,
94  fvMatrix<scalar>& eqn
95 ) const
96 {
97  const compressible::momentumTransportModel& momentumTransport =
98  mesh().lookupType<compressible::momentumTransportModel>();
99 
100  volVectorField& U = const_cast<volVectorField&>(momentumTransport.U());
101 
102  mesh().schemes().setFluxRequired(U.name());
103 
104  eqn -= fvc::div
105  (
106  fvc::dotInterpolate(momentumTransport.divDevTau(U)->flux(), U)
107  );
108 }
109 
110 
112 (
113  const volScalarField& alpha,
114  const volScalarField& rho,
115  const volScalarField& he,
116  fvMatrix<scalar>& eqn
117 ) const
118 {
119  const compressible::momentumTransportModel& momentumTransport =
120  mesh().lookupType<compressible::momentumTransportModel>(phaseName_);
121 
122  volVectorField& U = const_cast<volVectorField&>(momentumTransport.U());
123 
124  mesh().schemes().setFluxRequired(U.name());
125 
126  eqn -= fvc::div
127  (
128  fvc::dotInterpolate(momentumTransport.divDevTau(U)->flux(), U)
129  );
130 }
131 
132 
134 {
135  return true;
136 }
137 
138 
140 {}
141 
142 
144 {}
145 
146 
148 {}
149 
150 
152 {
153  if (fvModel::read(dict))
154  {
155  readCoeffs();
156  return true;
157  }
158  else
159  {
160  return false;
161  }
162 }
163 
164 
165 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static word groupName(Name name, const word &group)
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
Base class for single-phase compressible turbulence models.
virtual tmp< fvVectorMatrix > divDevTau(volVectorField &U) const =0
Return the stress matrix for the momentum equation.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
Applies the viscous heating source to the total energy equation.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual void addSup(const volScalarField &rho, const volScalarField &he, fvMatrix< scalar > &eqn) const
Add explicit contribution to compressible energy equation.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
viscousHeating(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
const volVectorField & U() const
Access function to velocity field.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the divergence of the given field.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
static tmp< SurfaceField< typename innerProduct< vector, Type >::type > > dotInterpolate(const surfaceVectorField &Sf, const VolField< Type > &tvf)
Interpolate field onto faces.
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
thermo he()
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31