leastSquaresFaGrad.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 "leastSquaresFaGrad.H"
30 #include "leastSquaresFaVectors.H"
31 #include "gaussFaGrad.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace fa
41 {
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 template<class Type>
45 tmp
46 <
47  GeometricField
48  <
50  faPatchField,
51  areaMesh
52  >
53 >
55 (
56  const GeometricField<Type, faPatchField, areaMesh>& vsf,
57  const word& name
58 ) const
59 {
60  typedef typename outerProduct<vector, Type>::type GradType;
61 
62  const faMesh& mesh = vsf.mesh();
63 
64  tmp<GeometricField<GradType, faPatchField, areaMesh>> tlsGrad
65  (
66  new GeometricField<GradType, faPatchField, areaMesh>
67  (
68  IOobject
69  (
70  name,
71  vsf.instance(),
72  vsf.db(),
75  ),
76  mesh,
77  dimensioned<GradType>(vsf.dimensions()/dimLength, Zero),
79  )
80  );
81  GeometricField<GradType, faPatchField, areaMesh>& lsGrad = tlsGrad.ref();
82 
83  // Get reference to least square vectors
84  const leastSquaresFaVectors& lsv = leastSquaresFaVectors::New(mesh);
85 
86  const edgeVectorField& ownLs = lsv.pVectors();
87  const edgeVectorField& neiLs = lsv.nVectors();
88 
89  const labelUList& own = mesh.owner();
90  const labelUList& nei = mesh.neighbour();
91 
92  forAll(own, edgei)
93  {
94  const label ownEdgeI = own[edgei];
95  const label neiEdgeI = nei[edgei];
96 
97  const Type deltaVsf(vsf[neiEdgeI] - vsf[ownEdgeI]);
98 
99  lsGrad[ownEdgeI] += ownLs[edgei]*deltaVsf;
100  lsGrad[neiEdgeI] -= neiLs[edgei]*deltaVsf;
101  }
102 
103  // Boundary edges
104  forAll(vsf.boundaryField(), patchi)
105  {
106  const faPatchField<Type>& bf = vsf.boundaryField()[patchi];
107 
108  tmp<Field<Type>> tvsfp(bf);
109 
110  if (bf.coupled())
111  {
112  tvsfp = bf.patchNeighbourField();
113  }
114  const auto& vsfp = tvsfp();
115 
116  const faePatchVectorField& ownLsp = ownLs.boundaryField()[patchi];
117  const labelUList& edgeFaces =
118  lsGrad.boundaryField()[patchi].patch().edgeFaces();
119 
120  forAll(vsfp, pEdgei)
121  {
122  lsGrad[edgeFaces[pEdgei]] +=
123  ownLsp[pEdgei]*(vsfp[pEdgei] - vsf[edgeFaces[pEdgei]]);
124  }
125  }
126 
127  // Remove component of gradient normal to surface (area)
128  lsGrad.correctBoundaryConditions();
129 
131 
132  return tlsGrad;
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
138 } // End namespace fa
139 
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 
142 } // End namespace Foam
143 
144 // ************************************************************************* //
static const leastSquaresFaVectors & New(const faMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
Definition: MeshObject.C:53
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:118
Ignore writing from objectRegistry::writeObject()
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: faPatchField.H:191
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
GeometricField< vector, faePatchField, edgeMesh > edgeVectorField
Definition: edgeFieldsFwd.H:60
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Nothing to be read.
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
faePatchField< vector > faePatchVectorField
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.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127