Loading...
Searching...
No Matches
linearUpwindNormal.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) 2007-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
31#include "fvMesh.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35template<class Type>
38(
40) const
41{
42 const fvMesh& mesh = this->mesh();
43
45 (
47 (
49 (
50 "linearUpwind::correction(" + vf.name() + ')',
51 mesh.time().timeName(),
52 mesh,
56 ),
57 mesh,
59 )
60 );
61
63
64 const surfaceScalarField& faceFlux = this->faceFlux_;
65
66 const auto& owner = mesh.owner();
67 const auto& neighbour = mesh.neighbour();
68
69 const volVectorField& C = mesh.C();
70 const surfaceVectorField& Cf = mesh.Cf();
71
72 tmp
73 <
75 <
79 >
80 > tgradVf = gradScheme_().grad(vf, gradSchemeName_);
81
83 <
87 >& gradVf = tgradVf();
88 gradVf /= mag(gradVf) + 1.e-12;
89
90 forAll(faceFlux, facei)
91 {
92 label celli = (faceFlux[facei] > 0) ? owner[facei] : neighbour[facei];
93 sfCorr[facei] = (Cf[facei] - C[celli]) & gradVf[celli];
94 }
95
96
98 GeometricBoundaryField& bSfCorr = sfCorr.boundaryField();
99
100 forAll(bSfCorr, patchi)
101 {
102 fvsPatchField<Type>& pSfCorr = bSfCorr[patchi];
103
104 if (pSfCorr.coupled())
105 {
106 const labelUList& pOwner =
107 mesh.boundary()[patchi].faceCells();
108
109 const vectorField& pCf = Cf.boundaryField()[patchi];
110
111 const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi];
112
114 (
115 gradVf.boundaryField()[patchi].patchNeighbourField()
116 );
117
118 // Build the d-vectors
119 vectorField pd(Cf.boundaryField()[patchi].patch().delta());
120
121 forAll(pOwner, facei)
122 {
123 label own = pOwner[facei];
124
125 if (pFaceFlux[facei] > 0)
126 {
127 pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own];
128 }
129 else
130 {
131 pSfCorr[facei] =
132 (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei];
133 }
134 }
135 }
136 }
137
138 return tsfCorr;
139}
140
141
142namespace Foam
143{
144 //makelimitedSurfaceInterpolationScheme(linearUpwindNormal)
147}
148
149// ************************************************************************* //
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.
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...
linearUpwindNormal interpolation scheme class derived from upwind and returns upwind weighting factor...
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
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.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299