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-2022 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 >
41 (
42  const VolField<Type>& vtf,
43  const word& name
44 ) const
45 {
46  typedef typename outerProduct<vector, Type>::type GradType;
47 
48  const fvMesh& mesh = vtf.mesh();
49 
50  // Get reference to least square vectors
51  const LeastSquaresVectors<Stencil>& lsv = LeastSquaresVectors<Stencil>::New
52  (
53  mesh
54  );
55 
56  tmp<VolField<GradType>> tlsGrad
57  (
59  (
60  name,
61  mesh,
62  dimensioned<GradType>
63  (
64  "zero",
65  vtf.dimensions()/dimLength,
66  Zero
67  ),
68  extrapolatedCalculatedFvPatchField<GradType>::typeName
69  )
70  );
71  VolField<GradType>& lsGrad = tlsGrad.ref();
72  Field<GradType>& lsGradIf = lsGrad;
73 
74  const extendedCentredCellToCellStencil& stencil = lsv.stencil();
75  const List<List<label>>& stencilAddr = stencil.stencil();
76  const List<List<vector>>& lsvs = lsv.vectors();
77 
78  // Construct flat version of vtf
79  // including all values referred to by the stencil
80  List<Type> flatVtf(stencil.map().constructSize(), Zero);
81 
82  // Insert internal values
83  forAll(vtf, celli)
84  {
85  flatVtf[celli] = vtf[celli];
86  }
87 
88  // Insert boundary values
89  forAll(vtf.boundaryField(), patchi)
90  {
91  const fvPatchField<Type>& ptf = vtf.boundaryField()[patchi];
92 
93  label nCompact =
94  ptf.patch().start()
95  - mesh.nInternalFaces()
96  + mesh.nCells();
97 
98  forAll(ptf, i)
99  {
100  flatVtf[nCompact++] = ptf[i];
101  }
102  }
103 
104  // Do all swapping to complete flatVtf
105  stencil.map().distribute(flatVtf);
106 
107  // Accumulate the cell-centred gradient from the
108  // weighted least-squares vectors and the flattened field values
109  forAll(stencilAddr, celli)
110  {
111  const labelList& compactCells = stencilAddr[celli];
112  const List<vector>& lsvc = lsvs[celli];
113 
114  forAll(compactCells, i)
115  {
116  lsGradIf[celli] += lsvc[i]*flatVtf[compactCells[i]];
117  }
118  }
119 
120  // Correct the boundary conditions
121  lsGrad.correctBoundaryConditions();
123 
124  return tlsGrad;
125 }
126 
127 
128 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
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.
A class for managing temporary objects.
Definition: tmp.H:55
label patchi
U correctBoundaryConditions()
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static const zero Zero
Definition: zero.H:97
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488