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) 2011-2016 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"
31 #include "surfaceMesh.H"
32 #include "GeometricField.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 <
41  <
45  >
46 >
48 (
49  const GeometricField<Type, fvPatchField, volMesh>& vsf,
50  const word& name
51 ) const
52 {
53  typedef typename outerProduct<vector, Type>::type GradType;
54 
55  const fvMesh& mesh = vsf.mesh();
56 
57  tmp<GeometricField<GradType, fvPatchField, volMesh>> tlsGrad
58  (
59  new GeometricField<GradType, fvPatchField, volMesh>
60  (
61  IOobject
62  (
63  name,
64  vsf.instance(),
65  mesh,
66  IOobject::NO_READ,
67  IOobject::NO_WRITE
68  ),
69  mesh,
70  dimensioned<GradType>
71  (
72  "zero",
73  vsf.dimensions()/dimLength,
74  Zero
75  ),
76  extrapolatedCalculatedFvPatchField<GradType>::typeName
77  )
78  );
79  GeometricField<GradType, fvPatchField, volMesh>& lsGrad = tlsGrad.ref();
80 
81  // Get reference to least square vectors
82  const leastSquaresVectors& lsv = leastSquaresVectors::New(mesh);
83 
84  const surfaceVectorField& ownLs = lsv.pVectors();
85  const surfaceVectorField& neiLs = lsv.nVectors();
86 
87  const labelUList& own = mesh.owner();
88  const labelUList& nei = mesh.neighbour();
89 
90  forAll(own, facei)
91  {
92  label ownFacei = own[facei];
93  label neiFacei = nei[facei];
94 
95  Type deltaVsf = vsf[neiFacei] - vsf[ownFacei];
96 
97  lsGrad[ownFacei] += ownLs[facei]*deltaVsf;
98  lsGrad[neiFacei] -= neiLs[facei]*deltaVsf;
99  }
100 
101  // Boundary faces
102  forAll(vsf.boundaryField(), patchi)
103  {
104  const fvsPatchVectorField& patchOwnLs = ownLs.boundaryField()[patchi];
105 
106  const labelUList& faceCells =
107  vsf.boundaryField()[patchi].patch().faceCells();
108 
109  if (vsf.boundaryField()[patchi].coupled())
110  {
111  const Field<Type> neiVsf
112  (
113  vsf.boundaryField()[patchi].patchNeighbourField()
114  );
115 
116  forAll(neiVsf, patchFacei)
117  {
118  lsGrad[faceCells[patchFacei]] +=
119  patchOwnLs[patchFacei]
120  *(neiVsf[patchFacei] - vsf[faceCells[patchFacei]]);
121  }
122  }
123  else
124  {
125  const fvPatchField<Type>& patchVsf = vsf.boundaryField()[patchi];
126 
127  forAll(patchVsf, patchFacei)
128  {
129  lsGrad[faceCells[patchFacei]] +=
130  patchOwnLs[patchFacei]
131  *(patchVsf[patchFacei] - vsf[faceCells[patchFacei]]);
132  }
133  }
134  }
135 
136 
137  lsGrad.correctBoundaryConditions();
139 
140  return tlsGrad;
141 }
142 
143 
144 // ************************************************************************* //
fvsPatchField< vector > fvsPatchVectorField
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
type
Types of root.
Definition: Roots.H:52
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
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
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.
UList< label > labelUList
Definition: UList.H:64
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:53
dynamicFvMesh & mesh
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.
static const zero Zero
Definition: zero.H:91
label patchi
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
U correctBoundaryConditions()
A class for managing temporary objects.
Definition: PtrList.H:53