blendingFactor.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-2016 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 "blendingFactor.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace functionObjects
38{
41}
42}
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
48{
49 writeHeader(os, "Blending factor");
50 writeCommented(os, "Time");
51 writeTabbed(os, "Scheme1");
52 writeTabbed(os, "Scheme2");
53 writeTabbed(os, "Blended");
54 os << endl;
55}
56
57
58bool Foam::functionObjects::blendingFactor::calc()
59{
60 bool processed = false;
61
62 processed = processed || calcScheme<scalar>();
63 processed = processed || calcScheme<vector>();
64
65 return processed;
66}
67
68
69// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70
72(
73 const word& name,
74 const Time& runTime,
75 const dictionary& dict
76)
77:
79 writeFile(obr_, name, typeName, dict),
80 phiName_("phi"),
81 tolerance_(0.001)
82{
83 read(dict);
85 setResultName(typeName, "");
86
87 auto indicatorPtr = tmp<volScalarField>::New
88 (
90 (
93 mesh_,
96 ),
97 mesh_,
99 zeroGradientFvPatchScalarField::typeName
100 );
101
102 store(resultName_, indicatorPtr);
103}
104
105
106// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107
109{
111 {
112 phiName_ = dict.getOrDefault<word>("phi", "phi");
113
114 tolerance_ =
115 dict.getCheckOrDefault
116 (
117 "tolerance",
118 0.001,
119 [](scalar tol){ return (tol > 0 && tol < 1); }
120 );
121
122 return true;
123 }
124
125 return false;
126}
127
128
130{
132 {
133 const volScalarField& indicator =
134 lookupObject<volScalarField>(resultName_);
135
136 // Generate scheme statistics
137 label nCellsScheme1 = 0;
138 label nCellsScheme2 = 0;
139 label nCellsBlended = 0;
140 for (const auto i : indicator)
141 {
142 if (i < tolerance_)
143 {
144 nCellsScheme1++;
145 }
146 else if (i > (1 - tolerance_))
147 {
148 nCellsScheme2++;
149 }
150 else
151 {
152 nCellsBlended++;
153 }
154 }
155
156 reduce(nCellsScheme1, sumOp<label>());
157 reduce(nCellsScheme2, sumOp<label>());
158 reduce(nCellsBlended, sumOp<label>());
159
160 Log << " scheme 1 cells : " << nCellsScheme1 << nl
161 << " scheme 2 cells : " << nCellsScheme2 << nl
162 << " blended cells : " << nCellsBlended << nl
163 << endl;
164
165 writeCurrentTime(file());
166
167 file()
168 << token::TAB << nCellsScheme1
169 << token::TAB << nCellsScheme2
170 << token::TAB << nCellsBlended
171 << endl;
172 }
173
174 return true;
175}
176
177
178// ************************************************************************* //
#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.
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.
Computes blending factor as an indicator about which of the schemes is active across the domain.
virtual void writeFileHeader(Ostream &os) const
Write the file header.
virtual bool write()
Write the blendingFactor.
virtual bool read(const dictionary &)
Read the blendingFactor data.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
word resultName_
Name of result field.
virtual bool write()
Write the result field.
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
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
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.
const dimensionSet dimless
Dimensionless.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
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