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);
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  (
80  time_.timeName(),
81  mesh_,
85  ),
86  mesh_,
88  );
89 
90  store(resultName_, tmodelRegions);
91 
93 }
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
103  dict.readIfPresent("result", resultName_);
104 
105  return true;
106 }
107 
108 
110 {
111  Log << type() << " " << name() << " execute:" << nl;
112 
114  lookupObjectRef<volScalarField>(resultName_);
115 
116 
117  if (foundObject<DESModelBase>(turbulenceModel::propertiesName))
118  {
119  const DESModelBase& model =
120  lookupObject<DESModelBase>
121  (
123  );
124 
125  DESModelRegions == model.LESRegion();
126 
127  const scalar prc =
128  gSum(DESModelRegions.primitiveField()*mesh_.V())
129  /gSum(mesh_.V())*100.0;
130 
131  file() << time_.value()
132  << token::TAB << prc
133  << token::TAB << 100.0 - prc
134  << endl;
135 
136  Log << " LES = " << prc << " % (volume)" << nl
137  << " RAS = " << 100.0 - prc << " % (volume)" << nl
138  << endl;
139  }
140  else
141  {
142  Log << " No DES turbulence model found in database" << nl
143  << endl;
144  }
145 
146  return true;
147 }
148 
149 
151 {
152  const volScalarField& DESModelRegions =
153  lookupObject<volScalarField>(resultName_);
154 
155  Log << type() << " " << name() << " output:" << nl
156  << " writing field " << DESModelRegions.name() << nl
157  << endl;
158 
159  DESModelRegions.write();
160 
161  return true;
162 }
163 
164 
165 // ************************************************************************* //
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:339
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual bool read(const dictionary &)
Read the DESModelRegions data.
Tab [isspace].
Definition: token.H:129
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual bool write()
Calculate the DESModelRegions and write.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Abstract base-class for Time/database function objects.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
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
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
Base class for DES models providing an interfaces to DES fields.
Definition: DESModelBase.H:49
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Computes an indicator field for detached eddy simulation (DES) turbulence calculations, where the values of the indicator mean:
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
word resultName_
Name of DES indicator field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
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
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
virtual tmp< volScalarField > LESRegion() const =0
Return the LES field indicator.
virtual void writeFileHeader(Ostream &os) const
File header information.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
#define Log
Definition: PDRblock.C:28
virtual bool read(const dictionary &dict)
Read optional controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
DESModelRegions(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Base class for writing single files from the function objects.
Definition: writeFile.H:112
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329
const Time & time_
Reference to the time database.