momentumError.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) 2020-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "momentumError.H"
29#include "fvcDiv.H"
30#include "fvcGrad.H"
31#include "fvcLaplacian.H"
32#include "turbulenceModel.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace functionObjects
42{
45}
46}
47
48
49// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50
53{
54 typedef compressible::turbulenceModel cmpTurbModel;
55 typedef incompressible::turbulenceModel icoTurbModel;
56
57 const auto& U = lookupObject<volVectorField>(UName_);
59
60 {
61 auto* turb = findObject<cmpTurbModel>
62 (
64 );
65
66 if (turb)
67 {
69 tmp<volScalarField> tnuEff = turb->nuEff();
70
72 {
73 const fvMeshSubset& subsetter = zoneSubSetPtr_->subsetter();
74
75 tU = subsetter.interpolate(U, false);
76 trho = subsetter.interpolate(turb->rho(), false);
77 tnuEff = subsetter.interpolate(turb->nuEff()(), false);
78 }
79
81 (
82 "divDevRhoReff",
83 - fvc::div
84 (
85 (trho()*tnuEff())
86 *dev2(T(fvc::grad(tU()))),
87 "div(((rho*nuEff)*dev2(T(grad(U)))))"
88 )
90 (
91 trho()*tnuEff(),
92 tU(),
93 "laplacian(nuEff,U)"
94 )
95 );
96 }
97 }
98
99 {
100 const auto* turb = findObject<icoTurbModel>
101 (
103 );
104
105 if (turb)
106 {
107 tmp<volScalarField> tnuEff = turb->nuEff();
108
109 if (zoneSubSetPtr_)
110 {
111 const fvMeshSubset& subsetter = zoneSubSetPtr_->subsetter();
112
113 tU = subsetter.interpolate(U, false);
114 tnuEff = subsetter.interpolate(turb->nuEff()(), false);
115 }
116
118 (
119 "divDevRhoReff",
120 - fvc::div
121 (
122 tnuEff()*dev2(T(fvc::grad(tU()))),
123 "div((nuEff*dev2(T(grad(U)))))"
124 )
125 - fvc::laplacian(tnuEff(), tU(), "laplacian(nuEff,U)")
126 );
127 }
128 }
129
130 return volVectorField::null();
131}
132
133
134// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
135
137(
138 const word& name,
139 const Time& runTime,
140 const dictionary& dict
141)
142:
144 pName_("p"),
145 UName_("U"),
146 phiName_("phi")
147{
148 read(dict);
149
150 const auto& phi =lookupObject<surfaceScalarField>(phiName_);
151
152 const dimensionSet momDims
153 (
155 );
156
157
158 volVectorField* momentPtr = nullptr;
159
160 word momName(scopedName("momentError"));
161
162 if (zoneSubSetPtr_)
163 {
164 const fvMesh& subMesh = zoneSubSetPtr_->subsetter().subMesh();
165
166 // momentErrorMap
167
168 momentPtr = new volVectorField
169 (
171 (
172 scopedName("momentErrorMap"),
173 subMesh.time().timeName(),
174 subMesh,
177 ),
178 subMesh,
179 dimensionedVector(momDims)
180 );
181
182 subMesh.objectRegistry::store(momentPtr);
183
184 momName = scopedName("momentErrorZone");
185 }
186
187 momentPtr = new volVectorField
188 (
190 (
191 momName,
192 time_.timeName(),
193 mesh_,
196 ),
197 mesh_,
198 dimensionedVector(momDims)
199 );
200
201 mesh_.objectRegistry::store(momentPtr);
202}
203
204
205// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
206
208{
210 {
211 Info<< type() << ' ' << name() << ':' << nl;
212
213 // Optional field name entries
214 if (dict.readIfPresent<word>("p", pName_))
215 {
216 Info<< " p: " << pName_ << endl;
217 }
218 if (dict.readIfPresent<word>("U", UName_))
219 {
220 Info<< " U: " << UName_ << endl;
221 }
222
223 if (dict.readIfPresent<word>("phi", phiName_))
224 {
225 Info<< " phi: " << phiName_ << endl;
226 }
227
228 if (dict.found("cellZones"))
229 {
230 zoneSubSetPtr_.reset(new Detail::zoneSubSet(mesh_, dict));
231 }
232 else
233 {
234 zoneSubSetPtr_.reset(nullptr);
235 }
236
237 return true;
238 }
239
240 return false;
241}
242
243
245{
246 const auto& p = lookupObject<volScalarField>(pName_);
247 const auto& U = lookupObject<volVectorField>(UName_);
248 const auto& phi = lookupObject<surfaceScalarField>(phiName_);
249
250 if (zoneSubSetPtr_)
251 {
252 const fvMeshSubset& subsetter = zoneSubSetPtr_->subsetter();
253
254 fvMesh& subMesh = zoneSubSetPtr_->subsetter().subMesh();
255
256 subMesh.fvSchemes::readOpt() = mesh_.fvSchemes::readOpt();
257 subMesh.fvSchemes::read();
258
259 auto& momentErrMap =
261 (
262 scopedName("momentErrorMap")
263 );
264
265 tmp<volScalarField> tp = subsetter.interpolate(p, false);
266 tmp<volVectorField> tU = subsetter.interpolate(U, false);
267 tmp<surfaceScalarField> tphi = subsetter.interpolate(phi, false);
268
269 momentErrMap =
270 (
271 divDevRhoReff()
272 + fvc::div(tphi, tU, "div(phi,U)")
273 + fvc::grad(tp, "grad(p)")
274 );
275 }
276 else
277 {
278 auto& momentErr =
279 lookupObjectRef<volVectorField>(scopedName("momentError"));
280
281 momentErr = fvc::div(phi, U) + fvc::grad(p) + divDevRhoReff();
282 }
283}
284
285
287{
288 calcMomentError();
289
290 return true;
291}
292
293
295{
296 Log << " functionObjects::" << type() << " " << name();
297
298 if (!zoneSubSetPtr_)
299 {
300 Log << " writing field: " << scopedName("momentError") << endl;
301
302 const auto& momentErr =
303 lookupObjectRef<volVectorField>(scopedName("momentError"));
304
305 momentErr.write();
306 }
307 else
308 {
309 Log << " writing field: " << scopedName("momentErrorMap") << endl;
310
311 const fvMeshSubset& subsetter = zoneSubSetPtr_->subsetter();
312 const fvMesh& subMesh = subsetter.subMesh();
313
314 const auto& momentErrMap =
316 (
317 scopedName("momentErrorMap")
318 );
319
320 tmp<volVectorField> mapMomErr =
321 zoneSubSetPtr_->mapToZone<vector>(momentErrMap);
322
323 mapMomErr().write();
324 }
325
326 return true;
327}
328
329
330// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
compressible::turbulenceModel & turb
const dimensionSet & dimensions() const
Return dimensions.
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated abstract base class for single-phase incompressible turbulence models.
virtual bool read()
Re-read model coefficients if they have changed.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Abstract base-class for Time/database function objects.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Computes balance terms for the steady momentum equation.
word UName_
Name of velocity field.
word phiName_
Name of flux field.
tmp< volVectorField > divDevRhoReff()
Return the effective viscous stress (laminar + turbulent).
Definition: momentumError.C:52
void calcMomentError()
Calculate the momentum error.
autoPtr< Detail::zoneSubSet > zoneSubSetPtr_
Sub-set mesh.
virtual bool read(const dictionary &)
Read the forces data.
const Time & time_
Reference to the time database.
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:80
static tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
Map volume internal (dimensioned) field.
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubsetI.H:48
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Type & lookupObjectRef(const word &name, const bool recursive=false) const
const Type & lookupObject(const word &name, const bool recursive=false) const
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & T
engineTime & runTime
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the laplacian of the given field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:47
Namespace for OpenFOAM.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
tmp< volScalarField > trho