blendingFactor.H
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) 2016-2022 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 Class
28  Foam::functionObjects::blendingFactor
29 
30 Group
31  grpFieldFunctionObjects
32 
33 Description
34  Computes blending factor as an indicator about which
35  of the schemes is active across the domain.
36 
37  Blending factor values mean:
38  \verbatim
39  0 = scheme 0
40  1 = scheme 1
41  0-1 = a blending factor between scheme 0 and scheme 1
42  \endverbatim
43 
44  Operands:
45  \table
46  Operand | Type | Location
47  input | - | -
48  output file | dat | $FOAM_CASE/postProcessing/<FO>/<time>/<file>
49  output field | volScalarField | $FOAM_CASE/<time>/<outField>
50  \endtable
51 
52 Usage
53  Minimal example by using \c system/controlDict.functions:
54  \verbatim
55  blendingFactor1
56  {
57  // Mandatory entries
58  type blendingFactor;
59  libs (fieldFunctionObjects);
60  field <field>;
61 
62  // Optional entries
63  phi phi;
64  tolerance 0.001;
65 
66  // Inherited entries
67  ...
68  }
69  \endverbatim
70 
71  where the entries mean:
72  \table
73  Property | Description | Type | Reqd | Deflt
74  type | Type name: blendingFactor | word | yes | -
75  libs | Library name: fieldFunctionObjects | word | yes | -
76  field | Name of the operand field | word | yes | -
77  phi | Name of flux field | word | no | phi
78  tolerance | Tolerance for number of blended cells | scalar | no | 0.001
79  \endtable
80 
81  The inherited entries are elaborated in:
82  - \link functionObject.H \endlink
83  - \link fieldExpression.H \endlink
84  - \link writeFile.H \endlink
85 
86  Usage by the \c postProcess utility is not available.
87 
88 Note
89  - The \c blendingFactor function object overwrites the \c DEShybrid:Factor
90  field internally when \c blendedSchemeBase debug flag is active.
91  However, users are allowed to write out the original \c DEShybrid:Factor
92  field by executing the \c writeObjects function object before
93  any \c blendingFactor function object execution.
94 
95 SourceFiles
96  blendingFactor.C
97  blendingFactorTemplates.C
98 
99 \*---------------------------------------------------------------------------*/
100 
101 #ifndef functionObjects_blendingFactor_H
102 #define functionObjects_blendingFactor_H
103 
104 #include "fieldExpression.H"
105 #include "writeFile.H"
106 #include "convectionScheme.H"
107 
108 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
109 
110 namespace Foam
111 {
112 namespace functionObjects
113 {
114 
115 /*---------------------------------------------------------------------------*\
116  Class blendingFactor Declaration
117 \*---------------------------------------------------------------------------*/
118 
119 class blendingFactor
120 :
121  public fieldExpression,
122  public writeFile
123 {
124  // Private Data
125 
126  //- Name of flux field
127  word phiName_;
128 
129  //- Tolerance used when calculating the number of blended cells
130  scalar tolerance_;
131 
132 
133  // Private Member Functions
134 
135  //- Calculate the blending factor field
136  template<class Type>
137  void calcBlendingFactor
138  (
139  const GeometricField<Type, fvPatchField, volMesh>& field,
140  const typename fv::convectionScheme<Type>& cs
141  );
142 
143  //- Calculate the blending factor field
144  template<class Type>
145  bool calcScheme();
146 
147  //- Calculate the blending factor field and return true if successful
148  virtual bool calc();
149 
150 
151 protected:
152 
153  // Protected Member Functions
154 
155  //- Write the file header
156  virtual void writeFileHeader(Ostream& os) const;
157 
158 
159 public:
160 
161  //- Runtime type information
162  TypeName("blendingFactor");
163 
164 
165  // Constructors
167  //- Construct from Time and dictionary
169  (
170  const word& name,
171  const Time& runTime,
172  const dictionary& dict
173  );
174 
175  //- No copy construct
176  blendingFactor(const blendingFactor&) = delete;
177 
178  //- No copy assignment
179  void operator=(const blendingFactor&) = delete;
180 
181 
182  //- Destructor
183  virtual ~blendingFactor() = default;
184 
185 
186  // Member Functions
187 
188  //- Read the blendingFactor data
189  virtual bool read(const dictionary&);
190 
191  //- Write the blendingFactor
192  virtual bool write();
193 };
194 
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 } // End namespace functionObjects
199 } // End namespace Foam
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 #ifdef NoRepository
204  #include "blendingFactorTemplates.C"
205 #endif
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 #endif
210 
211 // ************************************************************************* //
dictionary dict
rDeltaTY field()
Computes blending factor as an indicator about which of the schemes is active across the domain...
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
virtual ~blendingFactor()=default
Destructor.
engineTime & runTime
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
const word & name() const noexcept
Return the name of this functionObject.
virtual bool read(const dictionary &)
Read the blendingFactor data.
A class for handling words, derived from Foam::string.
Definition: word.H:63
OBJstream os(runTime.globalPath()/outputName)
void operator=(const blendingFactor &)=delete
No copy assignment.
virtual void writeFileHeader(Ostream &os) const
Write the file header.
TypeName("blendingFactor")
Runtime type information.
Namespace for OpenFOAM.
blendingFactor(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.