leastSquaresGrad.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 -------------------------------------------------------------------------------
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 "leastSquaresGrad.H"
30 #include "leastSquaresVectors.H"
31 #include "gaussGrad.H"
32 #include "fvMesh.H"
33 #include "volMesh.H"
34 #include "surfaceMesh.H"
35 #include "GeometricField.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 template<class Type>
42 <
44  <
48  >
49 >
51 (
52  const GeometricField<Type, fvPatchField, volMesh>& vsf,
53  const word& name
54 ) const
55 {
56  typedef typename outerProduct<vector, Type>::type GradType;
57  typedef GeometricField<GradType, fvPatchField, volMesh> GradFieldType;
58 
59  const fvMesh& mesh = vsf.mesh();
60 
61  tmp<GradFieldType> tlsGrad
62  (
63  new GradFieldType
64  (
65  IOobject
66  (
67  name,
68  vsf.instance(),
69  mesh,
70  IOobject::NO_READ,
71  IOobject::NO_WRITE
72  ),
73  mesh,
74  dimensioned<GradType>(vsf.dimensions()/dimLength, Zero),
76  )
77  );
78  GradFieldType& lsGrad = tlsGrad.ref();
79 
80  // Get reference to least square vectors
81  const leastSquaresVectors& lsv = leastSquaresVectors::New(mesh);
82 
83  const surfaceVectorField& ownLs = lsv.pVectors();
84  const surfaceVectorField& neiLs = lsv.nVectors();
85 
86  const labelUList& own = mesh.owner();
87  const labelUList& nei = mesh.neighbour();
88 
89  forAll(own, facei)
90  {
91  const label ownFacei = own[facei];
92  const label neiFacei = nei[facei];
93 
94  const Type deltaVsf(vsf[neiFacei] - vsf[ownFacei]);
95 
96  lsGrad[ownFacei] += ownLs[facei]*deltaVsf;
97  lsGrad[neiFacei] -= neiLs[facei]*deltaVsf;
98  }
99 
100  // Boundary faces
101  forAll(vsf.boundaryField(), patchi)
102  {
103  const fvsPatchVectorField& patchOwnLs = ownLs.boundaryField()[patchi];
104 
105  const labelUList& faceCells =
106  vsf.boundaryField()[patchi].patch().faceCells();
107 
108  if (vsf.boundaryField()[patchi].coupled())
109  {
110  const Field<Type> neiVsf
111  (
112  vsf.boundaryField()[patchi].patchNeighbourField()
113  );
114 
115  forAll(neiVsf, patchFacei)
116  {
117  lsGrad[faceCells[patchFacei]] +=
118  patchOwnLs[patchFacei]
119  *(neiVsf[patchFacei] - vsf[faceCells[patchFacei]]);
120  }
121  }
122  else
123  {
124  const fvPatchField<Type>& patchVsf = vsf.boundaryField()[patchi];
125 
126  forAll(patchVsf, patchFacei)
127  {
128  lsGrad[faceCells[patchFacei]] +=
129  patchOwnLs[patchFacei]
130  *(patchVsf[patchFacei] - vsf[faceCells[patchFacei]]);
131  }
132  }
133  }
134 
135 
136  lsGrad.correctBoundaryConditions();
138 
139  return tlsGrad;
140 }
141 
142 
143 // ************************************************************************* //
fvsPatchField< vector > fvsPatchVectorField
type
Types of root.
Definition: Roots.H:52
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
cellMask correctBoundaryConditions()
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:118
Generic GeometricField class.
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
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:45
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const word extrapolatedCalculatedType
A combined zero-gradient and calculated patch field type.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127