All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
gaussGrad.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "gaussGrad.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 template<class Type>
33 <
35 >
37 (
38  const SurfaceField<Type>& ssf,
39  const word& name
40 )
41 {
42  typedef typename outerProduct<vector, Type>::type GradType;
43 
44  const fvMesh& mesh = ssf.mesh();
45 
46  tmp<VolField<GradType>> tgGrad
47  (
49  (
50  name,
51  mesh,
52  dimensioned<GradType>
53  (
54  "0",
55  ssf.dimensions()/dimLength,
56  Zero
57  ),
58  extrapolatedCalculatedFvPatchField<GradType>::typeName
59  )
60  );
61  VolField<GradType>& gGrad = tgGrad.ref();
62 
63  const labelUList& owner = mesh.owner();
64  const labelUList& neighbour = mesh.neighbour();
65  const vectorField& Sf = mesh.Sf();
66 
67  Field<GradType>& igGrad = gGrad;
68  const Field<Type>& issf = ssf;
69 
70  forAll(owner, facei)
71  {
72  GradType Sfssf = Sf[facei]*issf[facei];
73 
74  igGrad[owner[facei]] += Sfssf;
75  igGrad[neighbour[facei]] -= Sfssf;
76  }
77 
78  forAll(mesh.boundary(), patchi)
79  {
80  const fvPatch& p = mesh.boundary()[patchi];
81  const labelUList& pFaceCells = p.faceCells();
82  const vectorField& pSf = mesh.Sf().boundaryField()[patchi];
83  const fvsPatchField<Type>& pssf = ssf.boundaryField()[patchi];
84 
85  forAll(p, facei)
86  {
87  igGrad[pFaceCells[facei]] += pSf[facei]*pssf[facei];
88  }
89  }
90 
91  igGrad /= mesh.V();
92 
93  gGrad.correctBoundaryConditions();
94 
95  return tgGrad;
96 }
97 
98 
99 template<class Type>
100 Foam::tmp
101 <
103 >
105 (
106  const VolField<Type>& vsf,
107  const word& name
108 ) const
109 {
110  typedef typename outerProduct<vector, Type>::type GradType;
111 
112  tmp<VolField<GradType>> tgGrad
113  (
114  gradf(tinterpScheme_().interpolate(vsf), name)
115  );
116  VolField<GradType>& gGrad = tgGrad.ref();
117 
118  correctBoundaryConditions(vsf, gGrad);
119 
120  return tgGrad;
121 }
122 
123 
124 template<class Type>
126 (
127  const VolField<Type>& vsf,
129 )
130 {
132  gGradbf = gGrad.boundaryFieldRef();
133 
134  forAll(vsf.boundaryField(), patchi)
135  {
136  if (!vsf.boundaryField()[patchi].coupled())
137  {
138  const vectorField n
139  (
140  vsf.mesh().Sf().boundaryField()[patchi]
141  / vsf.mesh().magSf().boundaryField()[patchi]
142  );
143 
144  gGradbf[patchi] += n *
145  (
146  vsf.boundaryField()[patchi].snGrad()
147  - (n & gGradbf[patchi])
148  );
149  }
150  }
151 }
152 
153 
154 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< VolField< typename outerProduct< vector, Type >::type > > gradf(const SurfaceField< Type > &, const word &name)
Return the gradient of the given field.
virtual tmp< VolField< typename outerProduct< vector, Type >::type > > calcGrad(const VolField< Type > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
static void correctBoundaryConditions(const VolField< Type > &, VolField< typename outerProduct< vector, Type >::type > &)
Correct the boundary values of the gradient using the patchField.
Definition: gaussGrad.C:126
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:90
A class for managing temporary objects.
Definition: tmp.H:55
volVectorField vectorField(fieldObject, mesh)
label patchi
U correctBoundaryConditions()
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static const zero Zero
Definition: zero.H:97
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
UList< label > labelUList
Definition: UList.H:65
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
volScalarField & p