linearUpwind.C
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-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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\*---------------------------------------------------------------------------*/
27
28#include "linearUpwind.H"
29#include "fvMesh.H"
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33template<class Type>
36(
38) const
39{
40 const fvMesh& mesh = this->mesh();
41
43 (
45 (
47 (
48 "linearUpwind::correction(" + vf.name() + ')',
49 mesh.time().timeName(),
50 mesh,
51 IOobject::NO_READ,
52 IOobject::NO_WRITE,
53 false
54 ),
55 mesh,
56 dimensioned<Type>(vf.name(), vf.dimensions(), Zero)
57 )
58 );
59
61
62 const surfaceScalarField& faceFlux = this->faceFlux_;
63
64 const labelList& owner = mesh.owner();
65 const labelList& neighbour = mesh.neighbour();
66
67 const volVectorField& C = mesh.C();
68 const surfaceVectorField& Cf = mesh.Cf();
69
71 (
73 (
74 mesh,
75 mesh.gradScheme(gradSchemeName_)
76 )
77 );
78
79 for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
80 {
81 tmp<volVectorField> tgradVf =
82 gradScheme_().grad(vf.component(cmpt), gradSchemeName_);
83
84 const volVectorField& gradVf = tgradVf();
85
86 forAll(faceFlux, facei)
87 {
88 const label celli =
89 (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
90
91 setComponent(sfCorr[facei], cmpt) =
92 (Cf[facei] - C[celli]) & gradVf[celli];
93 }
94
96 Boundary& bSfCorr = sfCorr.boundaryFieldRef();
97
98 forAll(bSfCorr, patchi)
99 {
100 fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
101
102 if (pSfCorr.coupled())
103 {
104 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
105 const vectorField& pCf = Cf.boundaryField()[patchi];
106 const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
107
108 const vectorField pGradVfNei
109 (
110 gradVf.boundaryField()[patchi].patchNeighbourField()
111 );
112
113 // Build the d-vectors
114 const vectorField pd
115 (
116 Cf.boundaryField()[patchi].patch().delta()
117 );
118
119 forAll(pOwner, facei)
120 {
121 label own = pOwner[facei];
122
123 if (pFaceFlux[facei] > 0)
124 {
125 setComponent(pSfCorr[facei], cmpt) =
126 (pCf[facei] - C[own])
127 & gradVf[own];
128 }
129 else
130 {
131 setComponent(pSfCorr[facei], cmpt) =
132 (pCf[facei] - pd[facei] - C[own])
133 & pGradVfNei[facei];
134 }
135 }
136 }
137 }
138 }
139
140 return tsfCorr;
141}
142
143
144template<>
147(
148 const volVectorField& vf
149) const
150{
151 const fvMesh& mesh = this->mesh();
152
154 (
156 (
158 (
159 "linearUpwind::correction(" + vf.name() + ')',
160 mesh.time().timeName(),
161 mesh,
162 IOobject::NO_READ,
163 IOobject::NO_WRITE,
164 false
165 ),
166 mesh,
167 dimensioned<vector>(vf.name(), vf.dimensions(), Zero)
168 )
169 );
170
171 surfaceVectorField& sfCorr = tsfCorr.ref();
172
173 const surfaceScalarField& faceFlux = this->faceFlux_;
174
175 const labelList& owner = mesh.owner();
176 const labelList& neighbour = mesh.neighbour();
177
178 const volVectorField& C = mesh.C();
179 const surfaceVectorField& Cf = mesh.Cf();
180
181 tmp<fv::gradScheme<vector>> gradScheme_
182 (
184 (
185 mesh,
186 mesh.gradScheme(gradSchemeName_)
187 )
188 );
189
190 tmp<volTensorField> tgradVf = gradScheme_().grad(vf, gradSchemeName_);
191 const volTensorField& gradVf = tgradVf();
192
193 forAll(faceFlux, facei)
194 {
195 const label celli =
196 (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
197 sfCorr[facei] = (Cf[facei] - C[celli]) & gradVf[celli];
198 }
199
200
201 typename surfaceVectorField::Boundary& bSfCorr = sfCorr.boundaryFieldRef();
202
203 forAll(bSfCorr, patchi)
204 {
205 fvsPatchVectorField& pSfCorr = bSfCorr[patchi];
206
207 if (pSfCorr.coupled())
208 {
209 const labelUList& pOwner = mesh.boundary()[patchi].faceCells();
210 const vectorField& pCf = Cf.boundaryField()[patchi];
211 const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
212
213 const tensorField pGradVfNei
214 (
215 gradVf.boundaryField()[patchi].patchNeighbourField()
216 );
217
218 // Build the d-vectors
219 vectorField pd(Cf.boundaryField()[patchi].patch().delta());
220
221 forAll(pOwner, facei)
222 {
223 label own = pOwner[facei];
224
225 if (pFaceFlux[facei] > 0)
226 {
227 pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own];
228 }
229 else
230 {
231 pSfCorr[facei] =
232 (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei];
233 }
234 }
235 }
236 }
237
238 return tsfCorr;
239}
240
241
242namespace Foam
243{
245}
246
247// ************************************************************************* //
Graphite solid properties.
Definition: C.H:53
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for gradient schemes.
Definition: gradScheme.H:66
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
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:60
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
Definition: linearUpwind.C:36
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
dynamicFvMesh & mesh
#define makelimitedSurfaceInterpolationScheme(SS)
Namespace for OpenFOAM.
label & setComponent(label &val, const direction) noexcept
Non-const access to integer-type (has no components)
Definition: label.H:123
uint8_t direction
Definition: direction.H:56
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333