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-2016 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 
44  (
45  functionObject,
46  DESModelRegions,
47  dictionary
48  );
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
54 
56 {
57  writeHeader(os, "DES model region coverage (% volume)");
58 
59  writeCommented(os, "Time");
60  writeTabbed(os, "LES");
61  writeTabbed(os, "RAS");
62  os << endl;
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
69 (
70  const word& name,
71  const Time& runTime,
72  const dictionary& dict
73 )
74 :
76  writeFile(obr_, name, typeName, dict),
77  resultName_(name)
78 {
79  read(dict);
80 
81  auto tmodelRegions = tmp<volScalarField>::New
82  (
83  IOobject
84  (
85  resultName_,
86  time_.timeName(),
87  mesh_,
90  ),
91  mesh_,
93  );
94 
95  store(resultName_, tmodelRegions);
96 
97  writeFileHeader(file());
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
102 
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
113 
114  dict.readIfPresent("result", resultName_);
115 
116  return true;
117 }
118 
119 
121 {
122  Log << type() << " " << name() << " execute:" << nl;
123 
125  lookupObjectRef<volScalarField>(resultName_);
126 
127 
128  if (foundObject<DESModelBase>(turbulenceModel::propertiesName))
129  {
130  const DESModelBase& model =
131  lookupObject<DESModelBase>
132  (
134  );
135 
136  DESModelRegions == model.LESRegion();
137 
138  scalar prc =
139  gSum(DESModelRegions.primitiveField()*mesh_.V())
140  /gSum(mesh_.V())*100.0;
141 
142  file() << time_.value()
143  << token::TAB << prc
144  << token::TAB << 100.0 - prc
145  << endl;
146 
147  Log << " LES = " << prc << " % (volume)" << nl
148  << " RAS = " << 100.0 - prc << " % (volume)" << nl
149  << endl;
150  }
151  else
152  {
153  Log << " No DES turbulence model found in database" << nl
154  << endl;
155  }
156 
157  return true;
158 }
159 
160 
162 {
164  lookupObject<volScalarField>(resultName_);
165 
166  Log << type() << " " << name() << " output:" << nl
167  << " writing field " << DESModelRegions.name() << nl
168  << endl;
169 
171 
172  return true;
173 }
174 
175 
176 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
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:104
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
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:62
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::functionObjects::DESModelRegions::DESModelRegions
DESModelRegions(const DESModelRegions &)=delete
No copy construct.
Foam::functionObjects::DESModelRegions
This function object writes out an indicator field for DES turbulence calculations,...
Definition: DESModelRegions.H:104
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:337
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::functionObjects::DESModelRegions::~DESModelRegions
virtual ~DESModelRegions()
Destructor.
Definition: DESModelRegions.C:103
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, add, dictionary)
Foam::functionObjects::writeFile::writeHeader
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:284
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:212
Foam::functionObjects::DESModelRegions::execute
virtual bool execute()
Execute.
Definition: DESModelRegions.C:120
Foam::functionObjects::DESModelRegions::read
virtual bool read(const dictionary &)
Read the DESModelRegions data.
Definition: DESModelRegions.C:109
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::functionObjects::DESModelRegions::write
virtual bool write()
Calculate the DESModelRegions and write.
Definition: DESModelRegions.C:161
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
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:121
Foam::functionObjects::DESModelRegions::writeFileHeader
virtual void writeFileHeader(Ostream &os) const
File header information.
Definition: DESModelRegions.C:55
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:166
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::defineTypeNameAndDebug
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
Foam::functionObjects::writeFile::writeCommented
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:258
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::functionObject::name
const word & name() const
Return the name of this functionObject.
Definition: functionObject.C:131
Foam::token::TAB
Tab [isspace].
Definition: token.H:113
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::writeFile
functionObject base class for writing single files
Definition: writeFile.H:59
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:123
Foam::functionObjects::writeFile::writeTabbed
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:274
turbulenceModel.H
Log
#define Log
Report write to Foam::Info if the local log switch is true.
Definition: messageStream.H:332