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,
64 IOobject::NO_READ,
65 IOobject::NO_WRITE
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
137 correctBoundaryConditions(vsf, gGrad);
138
139 return tgGrad;
140}
141
142
143template<class Type>
145(
148 <
150 >& gGrad
151)
152{
153 auto& gGradbf = gGrad.boundaryFieldRef();
154
155 forAll(vsf.boundaryField(), patchi)
156 {
157 if (!vsf.boundaryField()[patchi].coupled())
158 {
159 const vectorField n
160 (
161 vsf.mesh().Sf().boundaryField()[patchi]
162 / vsf.mesh().magSf().boundaryField()[patchi]
163 );
164
165 gGradbf[patchi] += n *
166 (
167 vsf.boundaryField()[patchi].snGrad()
168 - (n & gGradbf[patchi])
169 );
170 }
171 }
172}
173
174
175// ************************************************************************* //
label n
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
Generic dimensioned Type class.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
Basic second-order gradient scheme using face-interpolation and Gauss' theorem.
Definition: gaussGrad.H:66
static tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > gradf(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const word &name)
Definition: gaussGrad.C:45
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:54
type
Volume classification types.
Definition: volumeType.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
mesh interpolate(rAU)
dynamicFvMesh & mesh
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
cellMask correctBoundaryConditions()
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333