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 -------------------------------------------------------------------------------
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 "blendingFactor.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(blendingFactor, 0);
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 
58 bool 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  (
89  IOobject
90  (
92  time_.timeName(),
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 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
word resultName_
Name of result field.
Computes blending factor as an indicator about which of the schemes is active across the domain...
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition: writeFile.C:339
virtual bool write()
Write the blendingFactor.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Tab [isspace].
Definition: token.H:129
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual bool read(const dictionary &dict)
Read the fieldExpression data.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Abstract base-class for Time/database function objects.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual bool read(const dictionary &)
Read the blendingFactor data.
A class for handling words, derived from Foam::string.
Definition: word.H:63
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
virtual bool write()
Write the result field.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
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
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc...
virtual void writeFileHeader(Ostream &os) const
Write the file header.
Nothing to be read.
#define Log
Definition: PDRblock.C:28
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
Base class for writing single files from the function objects.
Definition: writeFile.H:112
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
Namespace for OpenFOAM.
blendingFactor(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329
const Time & time_
Reference to the time database.