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-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 "linearUpwind.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 template<>
33 (
34  const volVectorField& vf
35 ) const
36 {
37  const fvMesh& mesh = this->mesh();
38 
40  (
42  (
43  "linearUpwind::correction(" + vf.name() + ')',
44  mesh,
46  )
47  );
48 
49  surfaceVectorField& sfCorr = tsfCorr.ref();
50 
51  const surfaceScalarField& faceFlux = this->faceFlux_;
52 
53  const labelList& owner = mesh.owner();
54  const labelList& neighbour = mesh.neighbour();
55 
56  const volVectorField& C = mesh.C();
57  const surfaceVectorField& Cf = mesh.Cf();
58 
59  tmp<fv::gradScheme<vector>> gradScheme_
60  (
62  (
63  mesh,
64  mesh.schemes().grad(gradSchemeName_)
65  )
66  );
67 
68  tmp<volTensorField> tgradVf = gradScheme_().grad(vf, gradSchemeName_);
69  const volTensorField& gradVf = tgradVf();
70 
71  forAll(faceFlux, facei)
72  {
73  const label celli =
74  (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
75  sfCorr[facei] = (Cf[facei] - C[celli]) & gradVf[celli];
76  }
77 
78 
79  typename surfaceVectorField::Boundary& bSfCorr = sfCorr.boundaryFieldRef();
80 
81  forAll(bSfCorr, patchi)
82  {
83  fvsPatchVectorField& pSfCorr = bSfCorr[patchi];
84 
85  if (pSfCorr.coupled())
86  {
87  const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
88  const vectorField& pCf = Cf.boundaryField()[patchi];
89  const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
90 
91  const tensorField pGradVfNei
92  (
93  gradVf.boundaryField()[patchi].patchNeighbourField()
94  );
95 
96  // Build the d-vectors
97  vectorField pd(Cf.boundaryField()[patchi].patch().delta());
98 
99  forAll(pOwner, facei)
100  {
101  label own = pOwner[facei];
102 
103  if (pFaceFlux[facei] > 0)
104  {
105  pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own];
106  }
107  else
108  {
109  pSfCorr[facei] =
110  (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei];
111  }
112  }
113  }
114  }
115 
116  return tsfCorr;
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 namespace Foam
123 {
125 }
126 
127 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Graphite solid properties.
Definition: C.H:51
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const word & name() const
Return name.
Definition: IOobject.H:310
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const volVectorField & C() const
Return cell centres.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:451
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
const surfaceVectorField & Cf() const
Return face centres.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:457
ITstream & grad(const word &name) const
Definition: fvSchemes.C:485
Abstract base class for gradient schemes.
Definition: gradScheme.H:63
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
virtual bool coupled() const
Return true if this patch field is coupled.
linearUpwind interpolation scheme class derived from upwind and returns upwind weighting factors and ...
Definition: linearUpwind.H:55
virtual tmp< SurfaceField< Type > > correction(const VolField< Type > &) const
Return the explicit correction to the face-interpolate.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
label patchi
#define makelimitedSurfaceInterpolationScheme(SS)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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