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-------------------------------------------------------------------------------
11License
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
37namespace Foam
38{
39namespace functionObjects
40{
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 (
78 (
81 mesh_,
84 ),
85 mesh_,
87 );
88
89 store(resultName_, tmodelRegions);
90
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// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Base class for DES models providing an interfaces to the LESRegion function.
Definition: DESModelBase.H:55
virtual tmp< volScalarField > LESRegion() const =0
Return the LES field indicator.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual bool read()
Re-read model coefficients if they have changed.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Definition: TAB.H:69
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
Computes an indicator field for detached eddy simulation (DES) turbulence calculations,...
word resultName_
Name of DES indicator field.
virtual void writeFileHeader(Ostream &os) const
File header information.
virtual bool write()
Calculate the DESModelRegions and write.
virtual bool read(const dictionary &)
Read the DESModelRegions data.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
const Time & time_
Reference to the time database.
Base class for writing single files from the function objects.
Definition: writeFile.H:120
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:285
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:295
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:233
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:269
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
engineTime & runTime
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimless
Dimensionless.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict