gaussFaGrad.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "gaussFaGrad.H"
30 #include "facGrad.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fa
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class Type>
45 tmp
46 <
47  GeometricField
48  <
50  faPatchField,
51  areaMesh
52  >
53 >
55 (
56  const GeometricField<Type, faePatchField, edgeMesh>& ssf,
57  const word& name
58 )
59 {
60  const areaVectorField &n = ssf.mesh().faceAreaNormals();
61  typedef typename outerProduct<vector,Type>::type GradType;
62 
63  tmp<GeometricField<GradType, faPatchField, areaMesh>> tgGrad =
64  fac::edgeIntegrate(ssf.mesh().Sf()*ssf);
65 
66  GeometricField<GradType, faPatchField, areaMesh>& gGrad = tgGrad.ref();
67 
68  gGrad -= n*(n & gGrad);
69  gGrad.correctBoundaryConditions();
70 
71  return tgGrad;
72 }
73 
74 
75 template<class Type>
76 tmp
77 <
78  GeometricField
79  <
81  faPatchField,
82  areaMesh
83  >
84 >
86 (
87  const GeometricField<Type, faPatchField, areaMesh>& vsf,
88  const word& name
89 ) const
90 {
91  typedef typename outerProduct<vector, Type>::type GradType;
92  typedef GeometricField<GradType, faPatchField, areaMesh> GradFieldType;
93 
94  tmp<GradFieldType> tgGrad
95  (
97  (
98  vsf.mesh().Le()
99  *tinterpScheme_().interpolate(vsf)
100  )
101  );
102  GradFieldType& gGrad = tgGrad.ref();
103 
104  gGrad.correctBoundaryConditions();
105  gGrad.rename(name);
106 
107  correctBoundaryConditions(vsf, gGrad);
109  return tgGrad;
110 }
111 
112 
113 template<class Type>
115 (
118  <
120  >& gGrad
121 )
122 {
123  forAll(vsf.boundaryField(), patchI)
124  {
125  if (!vsf.boundaryField()[patchI].coupled())
126  {
127  const vectorField m
128  (
129  vsf.mesh().Le().boundaryField()[patchI]/
130  vsf.mesh().magLe().boundaryField()[patchI]
131  );
132 
133  gGrad.boundaryFieldRef()[patchI] += m*
134  (
135  vsf.boundaryField()[patchI].snGrad()
136  - (m & gGrad.boundaryField()[patchI])
137  );
138  }
139  }
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 } // End namespace fa
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 } // End namespace Foam
150 
151 // ************************************************************************* //
cellMask correctBoundaryConditions()
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:118
Generic GeometricField class.
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
Definition: areaFieldsFwd.H:46
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Calculate the gradient of the given field.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const Mesh & mesh() const noexcept
Return mesh.
tmp< GeometricField< Type, faPatchField, areaMesh > > edgeIntegrate(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh > > calcGrad(const GeometricField< Type, faPatchField, areaMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad for optional caching.
Mesh data needed to do the Finite Area discretisation.
Definition: areaFaMesh.H:47
static tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh > > gradf(const GeometricField< Type, faePatchField, edgeMesh > &, const word &name)
Return the gradient of the given field calculated using Gauss&#39; theorem on the given surface field...
static void correctBoundaryConditions(const GeometricField< Type, faPatchField, areaMesh > &, GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh > &)
Correct the boundary values of the gradient using the patchField.
Definition: gaussFaGrad.C:108
label n
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:76
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.