Loading...
Searching...
No Matches
gaussGrad.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 Copyright (C) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "gaussGrad.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34template<class Type>
36<
38 <
42 >
43>
45(
47 const word& name
48)
49{
50 typedef typename outerProduct<vector, Type>::type GradType;
52
53 const fvMesh& mesh = ssf.mesh();
54
56 (
57 new GradFieldType
58 (
60 (
61 name,
62 ssf.instance(),
63 mesh,
66 ),
67 mesh,
70 )
71 );
72 GradFieldType& gGrad = tgGrad.ref();
73
74 const labelUList& owner = mesh.owner();
75 const labelUList& neighbour = mesh.neighbour();
76 const vectorField& Sf = mesh.Sf();
77
78 Field<GradType>& igGrad = gGrad;
79 const Field<Type>& issf = ssf;
80
81 forAll(owner, facei)
82 {
83 const GradType Sfssf = Sf[facei]*issf[facei];
84
85 igGrad[owner[facei]] += Sfssf;
86 igGrad[neighbour[facei]] -= Sfssf;
87 }
88
89 forAll(mesh.boundary(), patchi)
90 {
91 const labelUList& pFaceCells =
92 mesh.boundary()[patchi].faceCells();
93
94 const vectorField& pSf = mesh.Sf().boundaryField()[patchi];
95
96 const fvsPatchField<Type>& pssf = ssf.boundaryField()[patchi];
97
98 forAll(mesh.boundary()[patchi], facei)
99 {
100 igGrad[pFaceCells[facei]] += pSf[facei]*pssf[facei];
101 }
102 }
103
104 igGrad /= mesh.V();
105
106 gGrad.correctBoundaryConditions();
107
108 return tgGrad;
109}
110
111
112template<class Type>
114<
116 <
120 >
121>
123(
125 const word& name
126) const
127{
128 typedef typename outerProduct<vector, Type>::type GradType;
130
131 tmp<GradFieldType> tgGrad
132 (
133 gradf(tinterpScheme_().interpolate(vsf), name)
134 );
135 GradFieldType& gGrad = tgGrad.ref();
136
138
139 return tgGrad;
140}
141
142
143template<class Type>
145(
147 <
150 volMesh
151 >& gGrad,
153) const
154{
155 typedef typename outerProduct<vector, Type>::type GradType;
156
157 DebugPout<< "gaussGrad<Type>::calcGrad on " << vsf.name()
158 << " into " << gGrad.name() << endl;
159
160 const fvMesh& mesh = vsf.mesh();
161 const labelUList& owner = mesh.owner();
162 const labelUList& neighbour = mesh.neighbour();
163 const vectorField& Sf = mesh.Sf();
164
165 const auto tssf(tinterpScheme_().interpolate(vsf));
166 const auto& ssf = tssf();
167
169
170 Field<GradType>& igGrad = gGrad;
171 const Field<Type>& issf = ssf;
172
173 forAll(owner, facei)
174 {
175 const GradType Sfssf = Sf[facei]*issf[facei];
176
177 igGrad[owner[facei]] += Sfssf;
178 igGrad[neighbour[facei]] -= Sfssf;
179 }
180
181 forAll(mesh.boundary(), patchi)
182 {
183 const labelUList& pFaceCells =
184 mesh.boundary()[patchi].faceCells();
185
186 const vectorField& pSf = mesh.Sf().boundaryField()[patchi];
187
188 const fvsPatchField<Type>& pssf = ssf.boundaryField()[patchi];
189
190 forAll(mesh.boundary()[patchi], facei)
191 {
192 igGrad[pFaceCells[facei]] += pSf[facei]*pssf[facei];
193 }
194 }
195
196 igGrad /= mesh.V();
197
198 gGrad.correctBoundaryConditions();
199
200 correctBoundaryConditions(vsf, gGrad);
201}
202
203
204template<class Type>
206(
209 <
211 >& gGrad
212)
213{
214 auto& gGradbf = gGrad.boundaryFieldRef();
215
216 forAll(vsf.boundaryField(), patchi)
217 {
218 if (!vsf.boundaryField()[patchi].coupled())
219 {
220 const vectorField n
221 (
222 vsf.mesh().Sf().boundaryField()[patchi]
223 / vsf.mesh().magSf().boundaryField()[patchi]
224 );
225
226 gGradbf[patchi] += n *
227 (
228 vsf.boundaryField()[patchi].snGrad()
229 - (n & gGradbf[patchi])
230 );
231 }
232 }
233}
234
235
236// ************************************************************************* //
label n
const Mesh & mesh() const noexcept
Return const reference to mesh.
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_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
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad for optional caching.
Definition gaussGrad.C:116
static void correctBoundaryConditions(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > &)
Correct the boundary values of the gradient using the patchField snGrad functions.
Definition gaussGrad.C:199
static tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > gradf(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const word &name)
Return the gradient of the given field calculated using Gauss' theorem on the given surface field.
Definition gaussGrad.C:38
const fvMesh & mesh() const
Return const reference to mesh.
Definition gradScheme.H:138
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
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
A class for handling words, derived from Foam::string.
Definition word.H:66
mesh interpolate(rAU)
dynamicFvMesh & mesh
#define DebugPout
Report an information message using Foam::Pout.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
cellMask correctBoundaryConditions()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299