outletStabilised.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "outletStabilised.H"
28 #include "mixedFvPatchField.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<class Type>
36 (
38 ) const
39 {
40  tmp<surfaceScalarField> tw = tScheme_().weights(vf);
41  surfaceScalarField& w = tw.ref();
42 
43  const fvMesh& mesh_ = this->mesh();
44  const cellList& cells = mesh_.cells();
45 
47  {
48  if
49  (
51  (vf.boundaryField()[patchi])
54  (vf.boundaryField()[patchi])
55  )
56  {
57  const labelList& pFaceCells =
58  mesh_.boundary()[patchi].faceCells();
59 
60  forAll(pFaceCells, pFacei)
61  {
62  const cell& pFaceCell = cells[pFaceCells[pFacei]];
63 
64  forAll(pFaceCell, fi)
65  {
66  label facei = pFaceCell[fi];
67 
68  if (mesh_.isInternalFace(facei))
69  {
70  // Apply upwind interpolation
71  w[facei] = pos0(faceFlux_[facei]);
72  }
73  }
74  }
75  }
76  }
77 
78  return tw;
79 }
80 
81 
82 template<class Type>
85 (
87 ) const
88 {
89  if (tScheme_().corrected())
90  {
92  tScheme_().correction(vf);
93 
95  tcorr.ref();
96 
97  const fvMesh& mesh_ = this->mesh();
98  const cellList& cells = mesh_.cells();
99 
101  {
102  if
103  (
105  (vf.boundaryField()[patchi])
107  (vf.boundaryField()[patchi])
108  )
109  {
110  const labelList& pFaceCells =
111  mesh_.boundary()[patchi].faceCells();
112 
113  forAll(pFaceCells, pFacei)
114  {
115  const cell& pFaceCell = cells[pFaceCells[pFacei]];
116 
117  forAll(pFaceCell, fi)
118  {
119  label facei = pFaceCell[fi];
120 
121  if (mesh_.isInternalFace(facei))
122  {
123  // Remove correction
124  corr[facei] = Zero;
125  }
126  }
127  }
128  }
129  }
130 
131  return tcorr;
132  }
133  else
134  {
136  (
137  nullptr
138  );
139  }
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 namespace Foam
146 {
148 }
149 
150 // ************************************************************************* //
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:434
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
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.
const cellList & cells() const
Generic GeometricField class.
fvMesh & mesh
Base class for direction-mixed boundary conditions.
Outlet-stabilised interpolation scheme which applies upwind interpolation to the faces of the cells a...
const cellShapeList & cells
makeSurfaceInterpolationScheme(cellCoBlended)
static const zero Zero
Definition: zero.H:97
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
dimensionedScalar pos0(const dimensionedScalar &ds)
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
label patchi
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Namespace for OpenFOAM.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800