viscousDissipation.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) 2017-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "viscousDissipation.H"
29 #include "fvMatrices.H"
32 #include "basicThermo.H"
34 
35 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace fv
40 {
41  defineTypeNameAndDebug(viscousDissipation, 0);
42  addToRunTimeSelectionTable(option, viscousDissipation, dictionary);
43 }
44 }
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 Foam::tmp<Foam::volScalarField> Foam::fv::viscousDissipation::rho() const
49 {
51  (
52  IOobject
53  (
54  "trho",
55  mesh_.time().timeName(),
56  mesh_,
59  ),
60  mesh_,
61  rho_
62  );
63 
64  if (rho_.value() > 0)
65  {
66  return trho;
67  }
68  else if (rhoName_ != "none")
69  {
70  trho.ref() = mesh_.lookupObject<volScalarField>(rhoName_);
71  return trho;
72  }
73 
75  << "Neither rhoName nor rho are specified."
76  << exit(FatalError);
77 
78  return nullptr;
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
85 (
86  const word& sourceName,
87  const word& modelType,
88  const dictionary& dict,
89  const fvMesh& mesh
90 )
91 :
92  fv::option(sourceName, modelType, dict, mesh),
93  UName_(coeffs_.getOrDefault<word>("U", "U")),
94  rhoName_(coeffs_.getOrDefault<word>("rho", "none")),
95  rho_
96  (
97  coeffs_.getOrDefault
98  (
99  "rhoInf",
101  )
102  )
103 {
104  const auto* thermoPtr =
106 
107  if (thermoPtr)
108  {
109  fieldNames_.resize(1, thermoPtr->he().name());
110  }
111 
112  if (fieldNames_.empty())
113  {
114  coeffs_.readEntry("fields", fieldNames_);
115  }
116 
117  if (fieldNames_.size() != 1)
118  {
120  << "settings are:" << fieldNames_ << exit(FatalError);
121  }
122 
124 }
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
130 Foam::fv::viscousDissipation::devRhoReff() const
131 {
132  // Incompressible
133  {
134  const auto* turbPtr =
135  mesh_.findObject<incompressible::turbulenceModel>
136  (
138  );
139 
140  if (turbPtr)
141  {
142  return tmp<volSymmTensorField>(rho()*turbPtr->devRhoReff());
143  }
144  }
145 
146  // Compressible
147  {
148  const auto* turbPtr =
149  mesh_.findObject<compressible::turbulenceModel>
150  (
152  );
153 
154  if (turbPtr)
155  {
156  return tmp<volSymmTensorField>(turbPtr->devRhoReff());
157  }
158  }
159 
161  << " The turbulence model is not found in the database."
163 
164  return nullptr;
165 }
166 
167 
169 (
170  const volScalarField& rho,
171  fvMatrix<scalar>& eqn,
172  const label fieldi
173 )
174 {
175  typedef typename outerProduct<vector, vector>::type GradType;
177 
178  const word gradUName("grad(" + UName_ + ')');
179 
180  auto tgradU = tmp<GradFieldType>::New
181  (
182  IOobject
183  (
184  "gradU",
185  mesh_.time().timeName(),
186  mesh_.time(),
189  ),
190  mesh_,
192  );
193 
194  // Cached?
195  const auto* gradUPtr = mesh_.findObject<GradFieldType>(gradUName);
196 
197  if (gradUPtr)
198  {
199  tgradU.ref() = *gradUPtr;
200  }
201  else
202  {
203  const auto& U = mesh_.lookupObject<volVectorField>(UName_);
204  tgradU.ref() = fvc::grad(U);
205  }
206 
207  const volScalarField D("D", devRhoReff() && tgradU.ref());
208 
209  eqn -= D;
210 }
211 
212 
213 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
const Type & value() const noexcept
Return const reference to value.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
tmp< volScalarField > trho
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible energy equation.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:118
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
dynamicFvMesh & mesh
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const dimensionSet dimDensity
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const dimensionedScalar & D
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:152
Namespace for OpenFOAM.
viscousDissipation(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:123
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127