All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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) 2011-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"
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  (
60  (
61  name,
62  mesh,
63  dimensioned<GradType>
64  (
65  "zero",
66  vsf.dimensions()/dimLength,
67  Zero
68  ),
69  extrapolatedCalculatedFvPatchField<GradType>::typeName
70  )
71  );
72  GeometricField<GradType, fvPatchField, volMesh>& lsGrad = tlsGrad.ref();
73 
74  // Get reference to least square vectors
75  const leastSquaresVectors& lsv = leastSquaresVectors::New(mesh);
76 
77  const surfaceVectorField& ownLs = lsv.pVectors();
78  const surfaceVectorField& neiLs = lsv.nVectors();
79 
80  const labelUList& own = mesh.owner();
81  const labelUList& nei = mesh.neighbour();
82 
83  forAll(own, facei)
84  {
85  label ownFacei = own[facei];
86  label neiFacei = nei[facei];
87 
88  Type deltaVsf = vsf[neiFacei] - vsf[ownFacei];
89 
90  lsGrad[ownFacei] += ownLs[facei]*deltaVsf;
91  lsGrad[neiFacei] -= neiLs[facei]*deltaVsf;
92  }
93 
94  // Boundary faces
95  forAll(vsf.boundaryField(), patchi)
96  {
97  const fvsPatchVectorField& patchOwnLs = ownLs.boundaryField()[patchi];
98 
99  const labelUList& faceCells =
100  vsf.boundaryField()[patchi].patch().faceCells();
101 
102  if (vsf.boundaryField()[patchi].coupled())
103  {
104  const Field<Type> neiVsf
105  (
106  vsf.boundaryField()[patchi].patchNeighbourField()
107  );
108 
109  forAll(neiVsf, patchFacei)
110  {
111  lsGrad[faceCells[patchFacei]] +=
112  patchOwnLs[patchFacei]
113  *(neiVsf[patchFacei] - vsf[faceCells[patchFacei]]);
114  }
115  }
116  else
117  {
118  const fvPatchField<Type>& patchVsf = vsf.boundaryField()[patchi];
119 
120  forAll(patchVsf, patchFacei)
121  {
122  lsGrad[faceCells[patchFacei]] +=
123  patchOwnLs[patchFacei]
124  *(patchVsf[patchFacei] - vsf[faceCells[patchFacei]]);
125  }
126  }
127  }
128 
129 
130  lsGrad.correctBoundaryConditions();
132 
133  return tlsGrad;
134 }
135 
136 
137 // ************************************************************************* //
fvsPatchField< vector > fvsPatchVectorField
#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()
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:65
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:53
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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: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