outletStabilised.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 | Copyright (C) 2011-2016 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::outletStabilised
26 
27 Description
28  Outlet-stabilised interpolation scheme which applies upwind differencing
29  to the faces of the cells adjacent to outlets.
30 
31  This is particularly useful to stabilise the velocity at entrainment
32  boundaries for LES cases using linear or other centred differencing
33  schemes.
34 
35 SourceFiles
36  outletStabilised.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef outletStabilised_H
41 #define outletStabilised_H
42 
44 #include "skewCorrectionVectors.H"
45 #include "linear.H"
46 #include "gaussGrad.H"
48 #include "mixedFvPatchField.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class outletStabilised Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class Type>
61 class outletStabilised
62 :
63  public surfaceInterpolationScheme<Type>
64 {
65  // Private member data
66 
67  const surfaceScalarField& faceFlux_;
69 
70 
71  // Private Member Functions
72 
73  //- Disallow default bitwise copy construct
75 
76  //- Disallow default bitwise assignment
77  void operator=(const outletStabilised&);
78 
79 
80 public:
81 
82  //- Runtime type information
83  TypeName("outletStabilised");
84 
85 
86  // Constructors
87 
88  //- Construct from mesh and Istream
90  (
91  const fvMesh& mesh,
92  Istream& is
93  )
94  :
95  surfaceInterpolationScheme<Type>(mesh),
96  faceFlux_
97  (
98  mesh.lookupObject<surfaceScalarField>
99  (
100  word(is)
101  )
102  ),
103  tScheme_
104  (
105  surfaceInterpolationScheme<Type>::New(mesh, faceFlux_, is)
106  )
107  {}
108 
109 
110  //- Construct from mesh, faceFlux and Istream
112  (
113  const fvMesh& mesh,
114  const surfaceScalarField& faceFlux,
115  Istream& is
116  )
117  :
118  surfaceInterpolationScheme<Type>(mesh),
119  faceFlux_(faceFlux),
120  tScheme_
121  (
122  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
123  )
124  {}
125 
126 
127  // Member Functions
128 
129  //- Return the interpolation weighting factors
131  (
133  ) const
134  {
135  tmp<surfaceScalarField> tw = tScheme_().weights(vf);
136  surfaceScalarField& w = tw.ref();
137 
138  const fvMesh& mesh_ = this->mesh();
139  const cellList& cells = mesh_.cells();
140 
142  {
143  if
144  (
146  (vf.boundaryField()[patchi])
149  (vf.boundaryField()[patchi])
150  )
151  {
152  const labelList& pFaceCells =
153  mesh_.boundary()[patchi].faceCells();
154 
155  forAll(pFaceCells, pFacei)
156  {
157  const cell& pFaceCell = cells[pFaceCells[pFacei]];
158 
159  forAll(pFaceCell, fi)
160  {
161  label facei = pFaceCell[fi];
162 
163  if (mesh_.isInternalFace(facei))
164  {
165  // Apply upwind differencing
166  w[facei] = pos(faceFlux_[facei]);
167  }
168  }
169  }
170  }
171  }
172 
173  return tw;
174  }
175 
176  //- Return true if this scheme uses an explicit correction
177  virtual bool corrected() const
178  {
179  return tScheme_().corrected();
180  }
181 
182  //- Return the explicit correction to the face-interpolate
183  // set to zero on the near-boundary faces where upwinf is applied
186  (
188  ) const
189  {
190  if (tScheme_().corrected())
191  {
193  tScheme_().correction(vf);
194 
196  tcorr.ref();
197 
198  const fvMesh& mesh_ = this->mesh();
199  const cellList& cells = mesh_.cells();
200 
202  {
203  if
204  (
206  (vf.boundaryField()[patchi])
208  (vf.boundaryField()[patchi])
209  )
210  {
211  const labelList& pFaceCells =
212  mesh_.boundary()[patchi].faceCells();
213 
214  forAll(pFaceCells, pFacei)
215  {
216  const cell& pFaceCell = cells[pFaceCells[pFacei]];
217 
218  forAll(pFaceCell, fi)
219  {
220  label facei = pFaceCell[fi];
221 
222  if (mesh_.isInternalFace(facei))
223  {
224  // Remove correction
225  corr[facei] = Zero;
226  }
227  }
228  }
229  }
230  }
231 
232  return tcorr;
233  }
234  else
235  {
237  (
238  NULL
239  );
240  }
241  }
242 };
243 
244 
245 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
246 
247 } // End namespace Foam
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 #endif
252 
253 // ************************************************************************* //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const fvMesh & mesh() const
Return mesh reference.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
This boundary condition provides a base class for &#39;mixed&#39; type boundary conditions, i.e. conditions that mix fixed value and patch-normal gradient conditions.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
const cellList & cells() const
TypeName("outletStabilised")
Runtime type information.
Base class for direction-mixed boundary conditions.
dimensionedScalar pos(const dimensionedScalar &ds)
Outlet-stabilised interpolation scheme which applies upwind differencing to the faces of the cells ad...
const cellShapeList & cells
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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.
static const zero Zero
Definition: zero.H:91
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
label patchi
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
Abstract base class for surface interpolation schemes.
Namespace for OpenFOAM.