Loading...
Searching...
No Matches
linearUpwindV.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 "linearUpwindV.H"
29#include "fvMesh.H"
30#include "volFields.H"
31#include "surfaceFields.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35template<class Type>
38(
40) const
41{
42 const fvMesh& mesh = this->mesh();
43
45 (
47 (
49 (
50 "linearUpwindV::correction(" + vf.name() + ')',
51 mesh.time().timeName(),
52 mesh,
56 ),
57 mesh,
59 )
60 );
61
63
64 const surfaceScalarField& faceFlux = this->faceFlux_;
65 const surfaceScalarField& w = mesh.weights();
66
67 const auto& own = mesh.owner();
68 const auto& nei = mesh.neighbour();
69
70 const volVectorField& C = mesh.C();
71 const surfaceVectorField& Cf = mesh.Cf();
72
73 tmp
74 <
76 <
80 >
81 > tgradVf = gradScheme_().grad(vf, gradSchemeName_);
82
83 const GeometricField
84 <
88 >& gradVf = tgradVf();
89
90 forAll(faceFlux, facei)
91 {
92 vector maxCorr;
93
94 if (faceFlux[facei] > 0.0)
95 {
96 maxCorr =
97 (1.0 - w[facei])*(vf[nei[facei]] - vf[own[facei]]);
98
99 sfCorr[facei] =
100 (Cf[facei] - C[own[facei]]) & gradVf[own[facei]];
101 }
102 else
103 {
104 maxCorr =
105 w[facei]*(vf[own[facei]] - vf[nei[facei]]);
106
107 sfCorr[facei] =
108 (Cf[facei] - C[nei[facei]]) & gradVf[nei[facei]];
109 }
110
111 scalar sfCorrs = magSqr(sfCorr[facei]);
112 scalar maxCorrs = sfCorr[facei] & maxCorr;
113
114 if (sfCorrs > 0)
115 {
116 if (maxCorrs < 0)
117 {
118 sfCorr[facei] = Zero;
119 }
120 else if (sfCorrs > maxCorrs)
121 {
122 sfCorr[facei] *= maxCorrs/(sfCorrs + VSMALL);
123 }
124 }
125 }
126
127
129 Boundary& bSfCorr = sfCorr.boundaryFieldRef();
130
131 forAll(bSfCorr, patchi)
132 {
133 fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
134
135 if (pSfCorr.coupled())
136 {
137 const labelUList& pOwner =
138 mesh.boundary()[patchi].faceCells();
139
140 const vectorField& pCf = Cf.boundaryField()[patchi];
141 const scalarField& pW = w.boundaryField()[patchi];
142
143 const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
144
146 (
147 gradVf.boundaryField()[patchi].patchNeighbourField()
148 );
149
150 const Field<Type> pVfNei
151 (
152 vf.boundaryField()[patchi].patchNeighbourField()
153 );
154
155 // Build the d-vectors
156 vectorField pd(Cf.boundaryField()[patchi].patch().delta());
157
158 forAll(pOwner, facei)
159 {
160 label own = pOwner[facei];
161
162 vector maxCorr;
163
164 if (pFaceFlux[facei] > 0)
165 {
166 pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own];
167
168 maxCorr = (1.0 - pW[facei])*(pVfNei[facei] - vf[own]);
169 }
170 else
171 {
172 pSfCorr[facei] =
173 (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei];
174
175 maxCorr = pW[facei]*(vf[own] - pVfNei[facei]);
176 }
177
178 scalar pSfCorrs = magSqr(pSfCorr[facei]);
179 scalar maxCorrs = pSfCorr[facei] & maxCorr;
180
181 if (pSfCorrs > 0)
182 {
183 if (maxCorrs < 0)
184 {
185 pSfCorr[facei] = Zero;
186 }
187 else if (pSfCorrs > maxCorrs)
188 {
189 pSfCorr[facei] *= maxCorrs/(pSfCorrs + VSMALL);
190 }
191 }
192 }
193 }
194 }
195
196 return tsfCorr;
197}
198
199
200namespace Foam
201{
203}
204
205
206// ************************************************************************* //
Graphite solid properties.
Definition C.H:49
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual bool coupled() const
True if the patch field is coupled.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
linearUpwindV interpolation scheme class derived from upwind and returns upwind weighting factors but...
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
dynamicFvMesh & mesh
#define makelimitedSurfaceInterpolationTypeScheme(SS, Type)
Namespace for OpenFOAM.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.