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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::outletStabilised
28 
29 Group
30  grpFvSurfaceInterpolationSchemes
31 
32 Description
33  Outlet-stabilised interpolation scheme which applies upwind differencing
34  to the faces of the cells adjacent to outlets.
35 
36  This is particularly useful to stabilise the velocity at entrainment
37  boundaries for LES cases using linear or other centred differencing
38  schemes.
39 
40 SourceFiles
41  outletStabilised.C
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #ifndef outletStabilised_H
46 #define outletStabilised_H
47 
49 #include "skewCorrectionVectors.H"
50 #include "linear.H"
51 #include "gaussGrad.H"
53 #include "mixedFvPatchField.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 /*---------------------------------------------------------------------------*\
62  Class outletStabilised Declaration
63 \*---------------------------------------------------------------------------*/
64 
65 template<class Type>
66 class outletStabilised
67 :
68  public surfaceInterpolationScheme<Type>
69 {
70  // Private member data
71 
72  const surfaceScalarField& faceFlux_;
74 
75 
76  // Private Member Functions
77 
78  //- No copy construct
79  outletStabilised(const outletStabilised&) = delete;
80 
81  //- No copy assignment
82  void operator=(const outletStabilised&) = delete;
83 
84 
85 public:
86 
87  //- Runtime type information
88  TypeName("outletStabilised");
89 
90 
91  // Constructors
92 
93  //- Construct from mesh and Istream
95  (
96  const fvMesh& mesh,
97  Istream& is
98  )
99  :
101  faceFlux_
102  (
103  mesh.lookupObject<surfaceScalarField>
104  (
105  word(is)
106  )
107  ),
108  tScheme_
109  (
110  surfaceInterpolationScheme<Type>::New(mesh, faceFlux_, is)
111  )
112  {}
113 
114 
115  //- Construct from mesh, faceFlux and Istream
117  (
118  const fvMesh& mesh,
119  const surfaceScalarField& faceFlux,
120  Istream& is
121  )
122  :
124  faceFlux_(faceFlux),
125  tScheme_
126  (
127  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
128  )
129  {}
130 
131 
132  // Member Functions
133 
134  //- Return the interpolation weighting factors
136  (
138  ) const
139  {
140  tmp<surfaceScalarField> tw = tScheme_().weights(vf);
141  surfaceScalarField& w = tw.ref();
142 
143  const fvMesh& mesh_ = this->mesh();
144  const cellList& cells = mesh_.cells();
145 
146  forAll(vf.boundaryField(), patchi)
147  {
148  if
149  (
151  (vf.boundaryField()[patchi])
152  || isA<mixedFvPatchField<Type>>(vf.boundaryField()[patchi])
154  (vf.boundaryField()[patchi])
155  )
156  {
157  const labelList& pFaceCells =
158  mesh_.boundary()[patchi].faceCells();
159 
160  forAll(pFaceCells, pFacei)
161  {
162  const cell& pFaceCell = cells[pFaceCells[pFacei]];
163 
164  forAll(pFaceCell, fi)
165  {
166  label facei = pFaceCell[fi];
167 
168  if (mesh_.isInternalFace(facei))
169  {
170  // Apply upwind differencing
171  w[facei] = pos0(faceFlux_[facei]);
172  }
173  }
174  }
175  }
176  }
177 
178  return tw;
179  }
180 
181  //- Return true if this scheme uses an explicit correction
182  virtual bool corrected() const
183  {
184  return tScheme_().corrected();
185  }
186 
187  //- Return the explicit correction to the face-interpolate
188  // set to zero on the near-boundary faces where upwind is applied
191  (
193  ) const
194  {
195  if (tScheme_().corrected())
196  {
198  tScheme_().correction(vf);
199 
201  tcorr.ref();
202 
203  const fvMesh& mesh_ = this->mesh();
204  const cellList& cells = mesh_.cells();
205 
206  forAll(vf.boundaryField(), patchi)
207  {
208  if
209  (
211  (vf.boundaryField()[patchi])
213  (vf.boundaryField()[patchi])
214  )
215  {
216  const labelList& pFaceCells =
217  mesh_.boundary()[patchi].faceCells();
218 
219  forAll(pFaceCells, pFacei)
220  {
221  const cell& pFaceCell = cells[pFaceCells[pFacei]];
222 
223  forAll(pFaceCell, fi)
224  {
225  label facei = pFaceCell[fi];
226 
227  if (mesh_.isInternalFace(facei))
228  {
229  // Remove correction
230  corr[facei] = Zero;
231  }
232  }
233  }
234  }
235  }
236 
237  return tcorr;
238  }
239  else
240  {
242  (
243  nullptr
244  );
245  }
246  }
247 };
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace Foam
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 #endif
257 
258 // ************************************************************************* //
skewCorrectionVectors.H
Foam::surfaceInterpolationScheme::New
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Definition: surfaceInterpolationScheme.C:40
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
directionMixedFvPatchField.H
Foam::outletStabilised::TypeName
TypeName("outletStabilised")
Runtime type information.
mixedFvPatchField.H
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
zeroGradientFvPatchField.H
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:64
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
gaussGrad.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::outletStabilised::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: outletStabilised.H:190
Foam::outletStabilised::weights
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Definition: outletStabilised.H:135
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::outletStabilised
Outlet-stabilised interpolation scheme which applies upwind differencing to the faces of the cells ad...
Definition: outletStabilised.H:65
Foam::mixedFvPatchField
This boundary condition provides a base class for 'mixed' type boundary conditions,...
Definition: mixedFvPatchField.H:123
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
Foam::isA
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:103
Foam::fvBoundaryMesh::faceCells
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
Definition: fvBoundaryMesh.C:134
Foam::outletStabilised::corrected
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: outletStabilised.H:181
Foam::List< cell >
Foam::surfaceInterpolationScheme
Abstract base class for surface interpolation schemes.
Definition: surfaceInterpolationScheme.H:57
Foam::surfaceInterpolationScheme::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: surfaceInterpolationScheme.H:144
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::directionMixedFvPatchField
Base class for direction-mixed boundary conditions.
Definition: directionMixedFvPatchField.H:54
surfaceInterpolationScheme.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62