Loading...
Searching...
No Matches
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"
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
59{
60 bool processed = false;
61
62 processed = processed || calcScheme<scalar>();
63 processed = processed || calcScheme<vector>();
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);
86
87 auto indicatorPtr = tmp<volScalarField>::New
88 (
90 (
93 mesh_,
97 ),
98 mesh_,
101 );
103 store(resultName_, indicatorPtr);
104}
105
106
107// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108
110{
112 {
113 phiName_ = dict.getOrDefault<word>("phi", "phi");
114
115 tolerance_ =
116 dict.getCheckOrDefault
117 (
118 "tolerance",
119 0.001,
120 [](scalar tol){ return (tol > 0 && tol < 1); }
121 );
122
123 return true;
124 }
125
126 return false;
127}
128
129
131{
133 {
134 const volScalarField& indicator =
135 lookupObject<volScalarField>(resultName_);
136
137 // Generate scheme statistics
138 label nCellsScheme1 = 0;
139 label nCellsScheme2 = 0;
140 label nCellsBlended = 0;
141 for (const auto i : indicator)
142 {
143 if (i < tolerance_)
144 {
145 nCellsScheme1++;
146 }
147 else if (i > (1 - tolerance_))
148 {
149 nCellsScheme2++;
150 }
151 else
152 {
153 nCellsBlended++;
154 }
155 }
156
157 reduce(nCellsScheme1, sumOp<label>());
158 reduce(nCellsScheme2, sumOp<label>());
159 reduce(nCellsBlended, sumOp<label>());
160
161 Log << " scheme 1 cells : " << nCellsScheme1 << nl
162 << " scheme 2 cells : " << nCellsScheme2 << nl
163 << " blended cells : " << nCellsBlended << nl
164 << endl;
165
166 writeCurrentTime(file());
167
168 file()
169 << token::TAB << nCellsScheme1
170 << token::TAB << nCellsScheme2
171 << token::TAB << nCellsBlended
172 << endl;
173 }
174
175 return true;
176}
177
178
179// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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 bool read(const dictionary &dict)
Read the function-object dictionary.
virtual void writeFileHeader(Ostream &os) const
Write the file header.
blendingFactor(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
virtual bool write()
Write the function-object results.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
word resultName_
Name of result field.
fieldExpression(const word &name, const Time &runTime, const dictionary &dict, const word &fieldName=word::null, const word &resultName=word::null)
Construct from name, Time and dictionary.
virtual bool write()
Write the function-object results.
virtual bool calc()=0
Calculate the components of the field and return true if successful.
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
const fvMesh & mesh_
Reference to the fvMesh.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
const objectRegistry & obr_
Reference to the region objectRegistry.
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:113
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition writeFile.C:344
virtual bool read(const dictionary &dict)
Read.
Definition writeFile.C:240
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition writeFile.C:354
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
@ TAB
Tab [isspace].
Definition token.H:142
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
engineTime & runTime
OBJstream os(runTime.globalPath()/outputName)
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict