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-2024 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 "momentumError.H"
29 #include "fvcDiv.H"
30 #include "fvcGrad.H"
31 #include "fvcLaplacian.H"
32 #include "turbulenceModel.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
43  defineTypeNameAndDebug(momentumError, 0);
44  addToRunTimeSelectionTable(functionObject, momentumError, dictionary);
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  {
68  tmp<volScalarField> trho = turb->rho();
69  tmp<volScalarField> tnuEff = turb->nuEff();
70 
71  if (zoneSubSetPtr_)
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 
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  (
154  phi.dimensions()*dimVelocity/dimVolume
155  );
156 
157  if (zoneSubSetPtr_)
158  {
159  const fvMesh& subMesh = zoneSubSetPtr_->subsetter().subMesh();
160 
161  // momentErrorMap
162 
163  volVectorField* fldPtr = new volVectorField
164  (
165  IOobject
166  (
167  scopedName("momentErrorMap"),
168  subMesh.time().timeName(),
169  subMesh.thisDb(),
173  ),
174  subMesh,
175  dimensionedVector(momDims, Zero)
176  );
177 
178  regIOobject::store(fldPtr);
179 
180  }
181 
182  const word momName =
183  (
185  ? scopedName("momentErrorZone")
186  : scopedName("momentError")
187  );
188 
189  volVectorField* fldPtr = new volVectorField
190  (
191  IOobject
192  (
193  momName,
194  time_.timeName(),
195  mesh_.thisDb(),
199  ),
200  mesh_,
201  dimensionedVector(momDims, Zero)
202  );
204  regIOobject::store(fldPtr);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
209 
211 {
213  {
214  Info<< type() << ' ' << name() << ':' << nl;
215 
216  // Optional field name entries
217  if (dict.readIfPresent<word>("p", pName_))
218  {
219  Info<< " p: " << pName_ << endl;
220  }
221  if (dict.readIfPresent<word>("U", UName_))
222  {
223  Info<< " U: " << UName_ << endl;
224  }
225 
226  if (dict.readIfPresent<word>("phi", phiName_))
227  {
228  Info<< " phi: " << phiName_ << endl;
229  }
230 
231  if (dict.found("cellZones"))
232  {
233  zoneSubSetPtr_.reset(new Detail::zoneSubSet(mesh_, dict));
234  }
235  else
236  {
237  zoneSubSetPtr_.reset(nullptr);
238  }
239 
240  return true;
241  }
242 
243  return false;
244 }
245 
246 
248 {
249  const auto& p = lookupObject<volScalarField>(pName_);
250  const auto& U = lookupObject<volVectorField>(UName_);
251  const auto& phi = lookupObject<surfaceScalarField>(phiName_);
252 
253  if (zoneSubSetPtr_)
254  {
255  const fvMeshSubset& subsetter = zoneSubSetPtr_->subsetter();
256 
257  fvMesh& subMesh = zoneSubSetPtr_->subsetter().subMesh();
258 
259  subMesh.fvSchemes::readOpt(mesh_.fvSchemes::readOpt());
260  subMesh.fvSchemes::read();
261 
262  auto& momentErrMap =
264  (
265  scopedName("momentErrorMap")
266  );
267 
268  tmp<volScalarField> tp = subsetter.interpolate(p, false);
269  tmp<volVectorField> tU = subsetter.interpolate(U, false);
270  tmp<surfaceScalarField> tphi = subsetter.interpolate(phi, false);
271 
272  momentErrMap =
273  (
274  divDevRhoReff()
275  + fvc::div(tphi, tU, "div(phi,U)")
276  + fvc::grad(tp, "grad(p)")
277  );
278  }
279  else
280  {
281  auto& momentErr =
282  lookupObjectRef<volVectorField>(scopedName("momentError"));
283 
284  momentErr = fvc::div(phi, U) + fvc::grad(p) + divDevRhoReff();
285  }
286 }
287 
288 
290 {
291  calcMomentError();
292 
293  return true;
294 }
295 
296 
298 {
299  Log << " functionObjects::" << type() << " " << name();
300 
301  if (!zoneSubSetPtr_)
302  {
303  Log << " writing field: " << scopedName("momentError") << endl;
304 
305  const auto& momentErr =
306  lookupObjectRef<volVectorField>(scopedName("momentError"));
307 
308  momentErr.write();
309  }
310  else
311  {
312  Log << " writing field: " << scopedName("momentErrorMap") << endl;
313 
314  const fvMeshSubset& subsetter = zoneSubSetPtr_->subsetter();
315  const fvMesh& subMesh = subsetter.subMesh();
316 
317  const auto& momentErrMap =
319  (
320  scopedName("momentErrorMap")
321  );
322 
323  tmp<volVectorField> mapMomErr =
324  zoneSubSetPtr_->mapToZone<vector>(momentErrMap);
325 
326  mapMomErr().write();
327  }
328 
329  return true;
330 }
331 
332 
333 // ************************************************************************* //
static tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
Map volume internal (dimensioned) field.
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
defineTypeNameAndDebug(ObukhovLength, 0)
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...
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
word UName_
Name of velocity field.
compressible::turbulenceModel & turb
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< volScalarField > trho
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
Ignore writing from objectRegistry::writeObject()
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:40
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
momentumError(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Macros for easy insertion into run-time selection tables.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubsetI.H:41
Templated abstract base class for single-phase incompressible turbulence models.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Calculate the gradient of the given field.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Calculate the laplacian of the given field.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
static const GeometricField< vector, fvPatchField, volMesh > & null() noexcept
Return a null GeometricField (reference to a nullObject).
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Calculate the divergence of the given field.
tmp< volVectorField > divDevRhoReff()
Return the effective viscous stress (laminar + turbulent).
Definition: momentumError.C:45
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word phiName_
Name of flux field.
void calcMomentError()
Calculate the momentum error.
virtual bool execute()
Execute.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
#define Log
Definition: PDRblock.C:28
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
volScalarField & p
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:180
Request registration (bool: true)
word scopedName(const word &name) const
Return a scoped (prefixed) name.
const fvMesh & mesh_
Reference to the fvMesh.
virtual bool read(const dictionary &)
Read the forces data.
Namespace for OpenFOAM.
autoPtr< Detail::zoneSubSet > zoneSubSetPtr_
Sub-set mesh.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity
const Time & time_
Reference to the time database.