linearUpwind.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-2019 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 "linearUpwind.H"
27 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 template<class Type>
34 (
36 ) const
37 {
38  const fvMesh& mesh = this->mesh();
39 
41  (
43  (
44  "linearUpwind::correction(" + vf.name() + ')',
45  mesh,
47  )
48  );
49 
51 
52  const surfaceScalarField& faceFlux = this->faceFlux_;
53 
54  const labelList& owner = mesh.owner();
55  const labelList& neighbour = mesh.neighbour();
56 
57  const volVectorField& C = mesh.C();
58  const surfaceVectorField& Cf = mesh.Cf();
59 
60  tmp<fv::gradScheme<scalar>> gradScheme_
61  (
63  (
64  mesh,
65  mesh.gradScheme(gradSchemeName_)
66  )
67  );
68 
69  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
70  {
71  tmp<volVectorField> tgradVf =
72  gradScheme_().grad(vf.component(cmpt), gradSchemeName_);
73 
74  const volVectorField& gradVf = tgradVf();
75 
76  forAll(faceFlux, facei)
77  {
78  const label celli =
79  (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
80 
81  setComponent(sfCorr[facei], cmpt) =
82  (Cf[facei] - C[celli]) & gradVf[celli];
83  }
84 
86  Boundary& bSfCorr = sfCorr.boundaryFieldRef();
87 
88  forAll(bSfCorr, patchi)
89  {
90  fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
91 
92  if (pSfCorr.coupled())
93  {
94  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
95  const vectorField& pCf = Cf.boundaryField()[patchi];
96  const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
97 
98  const vectorField pGradVfNei
99  (
100  gradVf.boundaryField()[patchi].patchNeighbourField()
101  );
102 
103  // Build the d-vectors
104  const vectorField pd
105  (
106  Cf.boundaryField()[patchi].patch().delta()
107  );
108 
109  forAll(pOwner, facei)
110  {
111  label own = pOwner[facei];
112 
113  if (pFaceFlux[facei] > 0)
114  {
115  setComponent(pSfCorr[facei], cmpt) =
116  (pCf[facei] - C[own])
117  & gradVf[own];
118  }
119  else
120  {
121  setComponent(pSfCorr[facei], cmpt) =
122  (pCf[facei] - pd[facei] - C[own])
123  & pGradVfNei[facei];
124  }
125  }
126  }
127  }
128  }
129 
130  return tsfCorr;
131 }
132 
133 
134 template<>
137 (
138  const volVectorField& vf
139 ) const
140 {
141  const fvMesh& mesh = this->mesh();
142 
144  (
146  (
147  "linearUpwind::correction(" + vf.name() + ')',
148  mesh,
150  )
151  );
152 
153  surfaceVectorField& sfCorr = tsfCorr.ref();
154 
155  const surfaceScalarField& faceFlux = this->faceFlux_;
156 
157  const labelList& owner = mesh.owner();
158  const labelList& neighbour = mesh.neighbour();
159 
160  const volVectorField& C = mesh.C();
161  const surfaceVectorField& Cf = mesh.Cf();
162 
163  tmp<fv::gradScheme<vector>> gradScheme_
164  (
166  (
167  mesh,
168  mesh.gradScheme(gradSchemeName_)
169  )
170  );
171 
172  tmp<volTensorField> tgradVf = gradScheme_().grad(vf, gradSchemeName_);
173  const volTensorField& gradVf = tgradVf();
174 
175  forAll(faceFlux, facei)
176  {
177  const label celli =
178  (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
179  sfCorr[facei] = (Cf[facei] - C[celli]) & gradVf[celli];
180  }
181 
182 
183  typename surfaceVectorField::Boundary& bSfCorr = sfCorr.boundaryFieldRef();
184 
185  forAll(bSfCorr, patchi)
186  {
187  fvsPatchVectorField& pSfCorr = bSfCorr[patchi];
188 
189  if (pSfCorr.coupled())
190  {
191  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
192  const vectorField& pCf = Cf.boundaryField()[patchi];
193  const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
194 
195  const tensorField pGradVfNei
196  (
197  gradVf.boundaryField()[patchi].patchNeighbourField()
198  );
199 
200  // Build the d-vectors
201  vectorField pd(Cf.boundaryField()[patchi].patch().delta());
202 
203  forAll(pOwner, facei)
204  {
205  label own = pOwner[facei];
206 
207  if (pFaceFlux[facei] > 0)
208  {
209  pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own];
210  }
211  else
212  {
213  pSfCorr[facei] =
214  (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei];
215  }
216  }
217  }
218  }
219 
220  return tsfCorr;
221 }
222 
223 
224 namespace Foam
225 {
227 }
228 
229 // ************************************************************************* //
#define makelimitedSurfaceInterpolationScheme(SS)
Graphite solid properties.
Definition: C.H:48
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 word & name() const
Return name.
Definition: IOobject.H:303
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
uint8_t direction
Definition: direction.H:45
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Generic GeometricField class.
Generic dimensioned Type class.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:284
const dimensionSet & dimensions() const
Return dimensions.
dynamicFvMesh & mesh
virtual bool coupled() const
Return true if this patch field is coupled.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
static const zero Zero
Definition: zero.H:97
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
ITstream & gradScheme(const word &name) const
Definition: fvSchemes.C:440
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:278
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
Abstract base class for gradient schemes.
Definition: gradScheme.H:60
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
linearUpwind interpolation scheme class derived from upwind and returns upwind weighting factors and ...
Definition: linearUpwind.H:52
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
Definition: linearUpwind.C:34
const volVectorField & C() const
Return cell centres as volVectorField.
A class for managing temporary objects.
Definition: PtrList.H:53
label & setComponent(label &l, const direction)
Definition: label.H:86
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540