DESModelRegions.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) 2013-2015 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "DESModelRegions.H"
30 #include "volFields.H"
31 #include "DESModelBase.H"
32 #include "turbulenceModel.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
41  defineTypeNameAndDebug(DESModelRegions, 0);
42  addToRunTimeSelectionTable(functionObject, DESModelRegions, dictionary);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
48 
50 {
51  writeHeader(os, "DES model region coverage (% volume)");
52 
53  writeCommented(os, "Time");
54  writeTabbed(os, "LES");
55  writeTabbed(os, "RAS");
56  os << endl;
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const word& name,
65  const Time& runTime,
66  const dictionary& dict
67 )
68 :
70  writeFile(obr_, name, typeName, dict),
71  resultName_(scopedName("regions"))
72 {
73  read(dict);
74 
75  auto tmodelRegions = tmp<volScalarField>::New
76  (
77  IOobject
78  (
79  resultName_,
80  time_.timeName(),
81  mesh_,
84  ),
85  mesh_,
87  );
88 
89  store(resultName_, tmodelRegions);
90 
91  writeFileHeader(file());
92 }
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
98 {
101 
102  dict.readIfPresent("result", resultName_);
103 
104  return true;
105 }
106 
107 
109 {
110  Log << type() << " " << name() << " execute:" << nl;
111 
113  lookupObjectRef<volScalarField>(resultName_);
114 
115 
116  if (foundObject<DESModelBase>(turbulenceModel::propertiesName))
117  {
118  const DESModelBase& model =
119  lookupObject<DESModelBase>
120  (
122  );
123 
124  DESModelRegions == model.LESRegion();
125 
126  const scalar prc =
127  gSum(DESModelRegions.primitiveField()*mesh_.V())
128  /gSum(mesh_.V())*100.0;
129 
130  file() << time_.value()
131  << token::TAB << prc
132  << token::TAB << 100.0 - prc
133  << endl;
134 
135  Log << " LES = " << prc << " % (volume)" << nl
136  << " RAS = " << 100.0 - prc << " % (volume)" << nl
137  << endl;
138  }
139  else
140  {
141  Log << " No DES turbulence model found in database" << nl
142  << endl;
143  }
144 
145  return true;
146 }
147 
148 
150 {
152  lookupObject<volScalarField>(resultName_);
153 
154  Log << type() << " " << name() << " output:" << nl
155  << " writing field " << DESModelRegions.name() << nl
156  << endl;
157 
159 
160  return true;
161 }
162 
163 
164 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
volFields.H
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Log
#define Log
Definition: PDRblock.C:35
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::functionObjects::DESModelRegions
Computes an indicator field for detached eddy simulation (DES) turbulence calculations,...
Definition: DESModelRegions.H:153
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::functionObjects::writeFile::writeHeader
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:298
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:213
Foam::functionObjects::DESModelRegions::execute
virtual bool execute()
Execute.
Definition: DESModelRegions.C:108
Foam::functionObjects::DESModelRegions::read
virtual bool read(const dictionary &)
Read the DESModelRegions data.
Definition: DESModelRegions.C:97
Foam::functionObjects::DESModelRegions::write
virtual bool write()
Calculate the DESModelRegions and write.
Definition: DESModelRegions.C:149
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::functionObjects::DESModelRegions::writeFileHeader
virtual void writeFileHeader(Ostream &os) const
File header information.
Definition: DESModelRegions.C:49
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
DESModelBase.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
DESModelRegions.H
Foam::functionObjects::DESModelRegions::DESModelRegions
DESModelRegions(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: DESModelRegions.C:63
Foam::functionObject::name
const word & name() const noexcept
Return the name of this functionObject.
Definition: functionObject.C:143
Foam::functionObjects::writeFile::writeCommented
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:272
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::token::TAB
Tab [isspace].
Definition: token.H:123
Foam::DESModelBase::LESRegion
virtual tmp< volScalarField > LESRegion() const =0
Return the LES field indicator.
Foam::DESModelBase
Base class for DES models providing an interfaces to the LESRegion function.
Definition: DESModelBase.H:54
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::functionObjects::writeFile::writeTabbed
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:288
turbulenceModel.H
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189