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 | Copyright (C) 2013 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  (
63  new GeometricField<GradType, fvPatchField, volMesh>
64  (
65  IOobject
66  (
67  name,
68  vtf.instance(),
69  mesh,
70  IOobject::NO_READ,
71  IOobject::NO_WRITE
72  ),
73  mesh,
74  dimensioned<GradType>
75  (
76  "zero",
77  vtf.dimensions()/dimLength,
78  pTraits<GradType>::zero
79  ),
80  zeroGradientFvPatchField<GradType>::typeName
81  )
82  );
83  GeometricField<GradType, fvPatchField, volMesh>& lsGrad = tlsGrad();
84  Field<GradType>& lsGradIf = lsGrad;
85 
86  const extendedCentredCellToCellStencil& stencil = lsv.stencil();
87  const List<List<label> >& stencilAddr = stencil.stencil();
88  const List<List<vector> >& lsvs = lsv.vectors();
89 
90  // Construct flat version of vtf
91  // including all values referred to by the stencil
92  List<Type> flatVtf(stencil.map().constructSize(), pTraits<Type>::zero);
93 
94  // Insert internal values
95  forAll(vtf, celli)
96  {
97  flatVtf[celli] = vtf[celli];
98  }
99 
100  // Insert boundary values
101  forAll(vtf.boundaryField(), patchi)
102  {
103  const fvPatchField<Type>& ptf = vtf.boundaryField()[patchi];
104 
105  label nCompact =
106  ptf.patch().start()
107  - mesh.nInternalFaces()
108  + mesh.nCells();
109 
110  forAll(ptf, i)
111  {
112  flatVtf[nCompact++] = ptf[i];
113  }
114  }
115 
116  // Do all swapping to complete flatVtf
117  stencil.map().distribute(flatVtf);
118 
119  // Accumulate the cell-centred gradient from the
120  // weighted least-squares vectors and the flattened field values
121  forAll(stencilAddr, celli)
122  {
123  const labelList& compactCells = stencilAddr[celli];
124  const List<vector>& lsvc = lsvs[celli];
125 
126  forAll(compactCells, i)
127  {
128  lsGradIf[celli] += lsvc[i]*flatVtf[compactCells[i]];
129  }
130  }
131 
132  // Correct the boundary conditions
133  lsGrad.correctBoundaryConditions();
135 
136  return tlsGrad;
137 }
138 
139 
140 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
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
Return the gradient of the given field to the gradScheme::grad.
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::volMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:47
volMesh.H
Foam::outerProduct::type
typeOfRank< typename pTraits< arg1 >::cmptType, int(pTraits< arg1 >::rank)+int(pTraits< arg2 >::rank) >::type type
Definition: products.H:72
LeastSquaresVectors.H
LeastSquaresGrad.H
zeroGradientFvPatchField.H
gaussGrad.H
Foam::label
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
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fvMesh.H
correctBoundaryConditions
U correctBoundaryConditions()
patchi
label patchi
Definition: getPatchFieldScalar.H:1
List
Definition: Test.C:19
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47