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 | Copyright (C) 2013-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "DESModelRegions.H"
27 #include "volFields.H"
28 #include "DESModelBase.H"
29 #include "turbulenceModel.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 defineTypeNameAndDebug(DESModelRegions, 0);
36 }
37 
38 
39 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
40 
42 {
43  writeHeader(os, "DES model region coverage (% volume)");
44 
45  writeCommented(os, "Time");
46  writeTabbed(os, "LES");
47  writeTabbed(os, "RAS");
48  os << endl;
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const word& name,
57  const objectRegistry& obr,
58  const dictionary& dict,
59  const bool loadFromFiles
60 )
61 :
62  functionObjectFile(obr, name, typeName, dict),
63  name_(name),
64  obr_(obr),
65  active_(true),
66  resultName_(name),
67  log_(true)
68 {
69  // Check if the available mesh is an fvMesh, otherwise deactivate
70  if (!isA<fvMesh>(obr_))
71  {
72  active_ = false;
74  << "No fvMesh available, deactivating " << name_ << nl
75  << endl;
76  }
77 
78  read(dict);
79 
80  if (active_)
81  {
82  const fvMesh& mesh = refCast<const fvMesh>(obr_);
83 
84  volScalarField* DESModelRegionsPtr
85  (
86  new volScalarField
87  (
88  IOobject
89  (
90  resultName_,
91  mesh.time().timeName(),
92  mesh,
95  ),
96  mesh,
97  dimensionedScalar("0", dimless, 0.0)
98  )
99  );
100 
101  mesh.objectRegistry::store(DESModelRegionsPtr);
102 
103  writeFileHeader(file());
104  }
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
109 
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  if (active_)
119  {
121 
122  log_.readIfPresent("log", dict);
123  dict.readIfPresent("resultName", resultName_);
124  }
125 }
126 
127 
129 {
130  if (active_)
131  {
132  const fvMesh& mesh = refCast<const fvMesh>(obr_);
133 
134  if (log_) Info<< type() << " " << name_ << " output:" << nl;
135 
137  const_cast<volScalarField&>
138  (
139  mesh.lookupObject<volScalarField>(resultName_)
140  );
141 
142 
144  {
145  const DESModelBase& model =
147  (
149  );
150 
151  DESModelRegions == model.LESRegion();
152 
153  scalar prc =
154  gSum(DESModelRegions.internalField()*mesh.V())
155  /gSum(mesh.V())*100.0;
156 
157  file() << obr_.time().value()
158  << token::TAB << prc
159  << token::TAB << 100.0 - prc
160  << endl;
161 
162  if (log_) Info
163  << " LES = " << prc << " % (volume)" << nl
164  << " RAS = " << 100.0 - prc << " % (volume)" << nl
165  << endl;
166  }
167  else
168  {
169  if (log_) Info
170  << " No DES turbulence model found in database" << nl
171  << endl;
172  }
173  }
174 }
175 
176 
178 {
179  // Do nothing
180 }
181 
182 
184 {
185  // Do nothing
186 }
187 
188 
190 {
191  if (active_)
192  {
194  obr_.lookupObject<volScalarField>(resultName_);
195 
196  if (log_) Info
197  << type() << " " << name_ << " output:" << nl
198  << " writing field " << DESModelRegions.name() << nl
199  << endl;
200 
202  }
203 }
204 
205 
206 // ************************************************************************* //
Foam::token::TAB
@ TAB
Definition: token.H:96
Foam::DESModelRegions::timeSet
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: DESModelRegions.C:183
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::DESModelRegions::execute
virtual void execute()
Execute, currently does nothing.
Definition: DESModelRegions.C:128
Foam::DESModelRegions::writeFileHeader
virtual void writeFileHeader(Ostream &os) const
File header information.
Definition: DESModelRegions.C:41
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::DESModelRegions::read
virtual void read(const dictionary &)
Read the DESModelRegions data.
Definition: DESModelRegions.C:116
Foam::functionObjectFile::writeTabbed
void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: functionObjectFile.C:230
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::functionObjectFile::read
void read(const dictionary &dict)
Read.
Definition: functionObjectFile.C:178
Foam::DESModelRegions::~DESModelRegions
virtual ~DESModelRegions()
Destructor.
Definition: DESModelRegions.C:110
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:199
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::functionObjectFile::writeHeader
void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: functionObjectFile.C:240
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
DESModelBase.H
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::DESModelRegions::DESModelRegions
DESModelRegions(const DESModelRegions &)
Disallow default bitwise copy construct.
Foam::DESModelRegions::write
virtual void write()
Calculate the DESModelRegions and write.
Definition: DESModelRegions.C:189
DESModelRegions.H
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
Foam::DESModelRegions
This function object writes out an indicator field for DES turbulence calculations,...
Definition: DESModelRegions.H:106
Foam::DESModelRegions::end
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: DESModelRegions.C:177
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::DESModelBase::LESRegion
virtual tmp< volScalarField > LESRegion() const =0
Return the LES field indicator.
Foam::functionObjectFile
Base class for output file data handling.
Definition: functionObjectFile.H:57
Foam::DESModelBase
Base class for DES models providing an interfaces to the LESRegion function.
Definition: DESModelBase.H:51
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::DESModelRegions::name
virtual const word & name() const
Return name of the set of DESModelRegions.
Definition: DESModelRegions.H:167
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::functionObjectFile::writeCommented
void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: functionObjectFile.C:219
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
turbulenceModel.H
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47