LeastSquaresGrad.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) 2013-2019 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 "LeastSquaresGrad.H"
27 #include "LeastSquaresVectors.H"
28 #include "gaussGrad.H"
29 #include "fvMesh.H"
30 #include "volMesh.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 template<class Type, class Stencil>
37 <
39  <
43  >
44 >
46 (
47  const GeometricField<Type, fvPatchField, volMesh>& vtf,
48  const word& name
49 ) const
50 {
51  typedef typename outerProduct<vector, Type>::type GradType;
52 
53  const fvMesh& mesh = vtf.mesh();
54 
55  // Get reference to least square vectors
56  const LeastSquaresVectors<Stencil>& lsv = LeastSquaresVectors<Stencil>::New
57  (
58  mesh
59  );
60 
61  tmp<GeometricField<GradType, fvPatchField, volMesh>> tlsGrad
62  (
64  (
65  name,
66  mesh,
67  dimensioned<GradType>
68  (
69  "zero",
70  vtf.dimensions()/dimLength,
71  Zero
72  ),
73  extrapolatedCalculatedFvPatchField<GradType>::typeName
74  )
75  );
76  GeometricField<GradType, fvPatchField, volMesh>& lsGrad = tlsGrad.ref();
77  Field<GradType>& lsGradIf = lsGrad;
78 
79  const extendedCentredCellToCellStencil& stencil = lsv.stencil();
80  const List<List<label>>& stencilAddr = stencil.stencil();
81  const List<List<vector>>& lsvs = lsv.vectors();
82 
83  // Construct flat version of vtf
84  // including all values referred to by the stencil
85  List<Type> flatVtf(stencil.map().constructSize(), Zero);
86 
87  // Insert internal values
88  forAll(vtf, celli)
89  {
90  flatVtf[celli] = vtf[celli];
91  }
92 
93  // Insert boundary values
94  forAll(vtf.boundaryField(), patchi)
95  {
96  const fvPatchField<Type>& ptf = vtf.boundaryField()[patchi];
97 
98  label nCompact =
99  ptf.patch().start()
100  - mesh.nInternalFaces()
101  + mesh.nCells();
102 
103  forAll(ptf, i)
104  {
105  flatVtf[nCompact++] = ptf[i];
106  }
107  }
108 
109  // Do all swapping to complete flatVtf
110  stencil.map().distribute(flatVtf);
111 
112  // Accumulate the cell-centred gradient from the
113  // weighted least-squares vectors and the flattened field values
114  forAll(stencilAddr, celli)
115  {
116  const labelList& compactCells = stencilAddr[celli];
117  const List<vector>& lsvc = lsvs[celli];
118 
119  forAll(compactCells, i)
120  {
121  lsGradIf[celli] += lsvc[i]*flatVtf[compactCells[i]];
122  }
123  }
124 
125  // Correct the boundary conditions
126  lsGrad.correctBoundaryConditions();
128 
129  return tlsGrad;
130 }
131 
132 
133 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
rDeltaT correctBoundaryConditions()
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:90
Generic GeometricField class.
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:53
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.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
label patchi
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
A class for managing temporary objects.
Definition: PtrList.H:53