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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 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"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 template<class Type, class Stencil>
40 <
42  <
46  >
47 >
49 (
50  const GeometricField<Type, fvPatchField, volMesh>& vtf,
51  const word& name
52 ) const
53 {
54  typedef typename outerProduct<vector, Type>::type GradType;
55  typedef GeometricField<GradType, fvPatchField, volMesh> GradFieldType;
56 
57  const fvMesh& mesh = vtf.mesh();
58 
59  // Get reference to least square vectors
60  const LeastSquaresVectors<Stencil>& lsv = LeastSquaresVectors<Stencil>::New
61  (
62  mesh
63  );
64 
65  tmp<GradFieldType> tlsGrad
66  (
67  new GradFieldType
68  (
69  IOobject
70  (
71  name,
72  vtf.instance(),
73  mesh,
74  IOobject::NO_READ,
75  IOobject::NO_WRITE
76  ),
77  mesh,
78  dimensioned<GradType>(vtf.dimensions()/dimLength, Zero),
79  extrapolatedCalculatedFvPatchField<GradType>::typeName
80  )
81  );
82  GradFieldType& lsGrad = tlsGrad.ref();
83  Field<GradType>& lsGradIf = lsGrad;
84 
85  const extendedCentredCellToCellStencil& stencil = lsv.stencil();
86  const List<List<label>>& stencilAddr = stencil.stencil();
87  const List<List<vector>>& lsvs = lsv.vectors();
88 
89  // Construct flat version of vtf
90  // including all values referred to by the stencil
91  List<Type> flatVtf(stencil.map().constructSize(), Zero);
92 
93  // Insert internal values
94  forAll(vtf, celli)
95  {
96  flatVtf[celli] = vtf[celli];
97  }
98 
99  // Insert boundary values
100  forAll(vtf.boundaryField(), patchi)
101  {
102  const fvPatchField<Type>& ptf = vtf.boundaryField()[patchi];
103 
104  label nCompact =
105  ptf.patch().start()
106  - mesh.nInternalFaces()
107  + mesh.nCells();
108 
109  forAll(ptf, i)
110  {
111  flatVtf[nCompact++] = ptf[i];
112  }
113  }
114 
115  // Do all swapping to complete flatVtf
116  stencil.map().distribute(flatVtf);
117 
118  // Accumulate the cell-centred gradient from the
119  // weighted least-squares vectors and the flattened field values
120  forAll(stencilAddr, celli)
121  {
122  const labelList& compactCells = stencilAddr[celli];
123  const List<vector>& lsvc = lsvs[celli];
124 
125  forAll(compactCells, i)
126  {
127  lsGradIf[celli] += lsvc[i]*flatVtf[compactCells[i]];
128  }
129  }
130 
131  // Correct the boundary conditions
132  lsGrad.correctBoundaryConditions();
134 
135  return tlsGrad;
136 }
137 
138 
139 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::fv::LeastSquaresGrad::calcGrad
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
extrapolatedCalculatedFvPatchField.H
Foam::volMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:51
Foam::outerProduct::type
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:114
volMesh.H
LeastSquaresVectors.H
LeastSquaresGrad.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
correctBoundaryConditions
cellMask correctBoundaryConditions()
gaussGrad.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
fvMesh.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53