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-2018 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  IOobject
45  (
46  "linearUpwind::correction(" + vf.name() + ')',
47  mesh.time().timeName(),
48  mesh,
49  IOobject::NO_READ,
50  IOobject::NO_WRITE,
51  false
52  ),
53  mesh,
55  )
56  );
57 
59 
60  const surfaceScalarField& faceFlux = this->faceFlux_;
61 
62  const labelList& owner = mesh.owner();
63  const labelList& neighbour = mesh.neighbour();
64 
65  const volVectorField& C = mesh.C();
66  const surfaceVectorField& Cf = mesh.Cf();
67 
68  tmp<fv::gradScheme<scalar>> gradScheme_
69  (
71  (
72  mesh,
73  mesh.gradScheme(gradSchemeName_)
74  )
75  );
76 
77  for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
78  {
79  tmp<volVectorField> tgradVf =
80  gradScheme_().grad(vf.component(cmpt), gradSchemeName_);
81 
82  const volVectorField& gradVf = tgradVf();
83 
84  forAll(faceFlux, facei)
85  {
86  const label celli =
87  (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
88 
89  setComponent(sfCorr[facei], cmpt) =
90  (Cf[facei] - C[celli]) & gradVf[celli];
91  }
92 
94  Boundary& bSfCorr = sfCorr.boundaryFieldRef();
95 
96  forAll(bSfCorr, patchi)
97  {
98  fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
99 
100  if (pSfCorr.coupled())
101  {
102  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
103  const vectorField& pCf = Cf.boundaryField()[patchi];
104  const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
105 
106  const vectorField pGradVfNei
107  (
108  gradVf.boundaryField()[patchi].patchNeighbourField()
109  );
110 
111  // Build the d-vectors
112  const vectorField pd
113  (
114  Cf.boundaryField()[patchi].patch().delta()
115  );
116 
117  forAll(pOwner, facei)
118  {
119  label own = pOwner[facei];
120 
121  if (pFaceFlux[facei] > 0)
122  {
123  setComponent(pSfCorr[facei], cmpt) =
124  (pCf[facei] - C[own])
125  & gradVf[own];
126  }
127  else
128  {
129  setComponent(pSfCorr[facei], cmpt) =
130  (pCf[facei] - pd[facei] - C[own])
131  & pGradVfNei[facei];
132  }
133  }
134  }
135  }
136  }
137 
138  return tsfCorr;
139 }
140 
141 
142 template<>
145 (
146  const volVectorField& vf
147 ) const
148 {
149  const fvMesh& mesh = this->mesh();
150 
152  (
154  (
155  "linearUpwind::correction(" + vf.name() + ')',
156  mesh,
158  )
159  );
160 
161  surfaceVectorField& sfCorr = tsfCorr.ref();
162 
163  const surfaceScalarField& faceFlux = this->faceFlux_;
164 
165  const labelList& owner = mesh.owner();
166  const labelList& neighbour = mesh.neighbour();
167 
168  const volVectorField& C = mesh.C();
169  const surfaceVectorField& Cf = mesh.Cf();
170 
171  tmp<fv::gradScheme<vector>> gradScheme_
172  (
174  (
175  mesh,
176  mesh.gradScheme(gradSchemeName_)
177  )
178  );
179 
180  tmp<volTensorField> tgradVf = gradScheme_().grad(vf, gradSchemeName_);
181  const volTensorField& gradVf = tgradVf();
182 
183  forAll(faceFlux, facei)
184  {
185  const label celli =
186  (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
187  sfCorr[facei] = (Cf[facei] - C[celli]) & gradVf[celli];
188  }
189 
190 
191  typename surfaceVectorField::Boundary& bSfCorr = sfCorr.boundaryFieldRef();
192 
193  forAll(bSfCorr, patchi)
194  {
195  fvsPatchVectorField& pSfCorr = bSfCorr[patchi];
196 
197  if (pSfCorr.coupled())
198  {
199  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
200  const vectorField& pCf = Cf.boundaryField()[patchi];
201  const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
202 
203  const tensorField pGradVfNei
204  (
205  gradVf.boundaryField()[patchi].patchNeighbourField()
206  );
207 
208  // Build the d-vectors
209  vectorField pd(Cf.boundaryField()[patchi].patch().delta());
210 
211  forAll(pOwner, facei)
212  {
213  label own = pOwner[facei];
214 
215  if (pFaceFlux[facei] > 0)
216  {
217  pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own];
218  }
219  else
220  {
221  pSfCorr[facei] =
222  (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei];
223  }
224  }
225  }
226  }
227 
228  return tsfCorr;
229 }
230 
231 
232 namespace Foam
233 {
235 }
236 
237 // ************************************************************************* //
#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:295
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.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
Generic dimensioned Type class.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
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.
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:79
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
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