localBlended.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::localBlended
26 
27 Description
28  Two-scheme localBlended interpolation scheme.
29 
30 SourceFiles
31  localBlended.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef localBlended_H
36 #define localBlended_H
37 
39 #include "blendedSchemeBase.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class localBlended Declaration
48 \*---------------------------------------------------------------------------*/
49 
50 template<class Type>
51 class localBlended
52 :
53  public surfaceInterpolationScheme<Type>,
54  public blendedSchemeBase<Type>
55 {
56  // Private Member Functions
57 
58  //- Scheme 1
60 
61  //- Scheme 2
63 
64 
65 public:
66 
67  //- Runtime type information
68  TypeName("localBlended");
69 
70 
71  // Constructors
72 
73  //- Construct from mesh and Istream.
74  // The name of the flux field is read from the Istream and looked-up
75  // from the mesh objectRegistry
77  (
78  const fvMesh& mesh,
79  Istream& is
80  )
81  :
82  surfaceInterpolationScheme<Type>(mesh),
83  tScheme1_
84  (
85  surfaceInterpolationScheme<Type>::New(mesh, is)
86  ),
87  tScheme2_
88  (
89  surfaceInterpolationScheme<Type>::New(mesh, is)
90  )
91  {}
92 
93  //- Construct from mesh, faceFlux and Istream
95  (
96  const fvMesh& mesh,
97  const surfaceScalarField& faceFlux,
98  Istream& is
99  )
100  :
101  surfaceInterpolationScheme<Type>(mesh),
102  tScheme1_
103  (
104  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
105  ),
106  tScheme2_
107  (
108  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
109  )
110  {}
111 
112  //- Disallow default bitwise copy construction
113  localBlended(const localBlended&) = delete;
114 
115 
116  //- Destructor
117  virtual ~localBlended()
118  {}
119 
120 
121  // Member Functions
122 
123  //- Return the face-based blending factor
125  (
127  ) const
128  {
129  return
130  this->mesh().objectRegistry::template
131  lookupObject<const surfaceScalarField>
132  (
133  word(vf.name() + "BlendingFactor")
134  );
135  }
136 
137  //- Return the interpolation weighting factors
139  (
141  ) const
142  {
144  this->mesh().objectRegistry::template
145  lookupObject<const surfaceScalarField>
146  (
147  word(vf.name() + "BlendingFactor")
148  );
149 
150  return
151  blendingFactor*tScheme1_().weights(vf)
152  + (scalar(1) - blendingFactor)*tScheme2_().weights(vf);
153  }
154 
155  //- Return the face-interpolate of the given cell field
156  // with explicit correction
159  {
161  (
162  this->mesh().objectRegistry::template
163  lookupObject<const surfaceScalarField>
164  (
165  word(vf.name() + "BlendingFactor")
166  )
167  );
168 
169  return
170  blendingFactor*tScheme1_().interpolate(vf)
171  + (scalar(1) - blendingFactor)*tScheme2_().interpolate(vf);
172  }
173 
174 
175  //- Return true if this scheme uses an explicit correction
176  virtual bool corrected() const
177  {
178  return tScheme1_().corrected() || tScheme2_().corrected();
179  }
180 
181 
182  //- Return the explicit correction to the face-interpolate
183  // for the given field
186  (
188  ) const
189  {
191  this->mesh().objectRegistry::template
192  lookupObject<const surfaceScalarField>
193  (
194  word(vf.name() + "BlendingFactor")
195  );
196 
197  if (tScheme1_().corrected())
198  {
199  if (tScheme2_().corrected())
200  {
201  return
202  (
203  blendingFactor
204  * tScheme1_().correction(vf)
205  + (scalar(1) - blendingFactor)
206  * tScheme2_().correction(vf)
207  );
208  }
209  else
210  {
211  return
212  (
213  blendingFactor
214  * tScheme1_().correction(vf)
215  );
216  }
217  }
218  else if (tScheme2_().corrected())
219  {
220  return
221  (
222  (scalar(1) - blendingFactor)
223  * tScheme2_().correction(vf)
224  );
225  }
226  else
227  {
229  (
230  nullptr
231  );
232  }
233  }
234 
235 
236  // Member Operators
237 
238  //- Disallow default bitwise assignment
239  void operator=(const localBlended&) = delete;
240 };
241 
242 
243 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 
245 } // End namespace Foam
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 #endif
250 
251 // ************************************************************************* //
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
Definition: localBlended.H:124
const word & name() const
Return name.
Definition: IOobject.H:315
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
localBlended(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
Definition: localBlended.H:76
const fvMesh & mesh() const
Return mesh reference.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Definition: localBlended.H:157
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Definition: localBlended.H:138
void operator=(const localBlended &)=delete
Disallow default bitwise assignment.
virtual ~localBlended()
Destructor.
Definition: localBlended.H:116
A class for handling words, derived from string.
Definition: word.H:59
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Base class for blended schemes to provide access to the blending factor surface field.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: localBlended.H:175
Two-scheme localBlended interpolation scheme.
Definition: localBlended.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: localBlended.H:185
A class for managing temporary objects.
Definition: PtrList.H:53
TypeName("localBlended")
Runtime type information.
Abstract base class for surface interpolation schemes.
Namespace for OpenFOAM.