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-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"
31 #include "surfaceMesh.H"
32 #include "GeometricField.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 <
41 >
43 (
44  const VolField<Type>& vsf,
45  const word& name
46 ) const
47 {
48  typedef typename outerProduct<vector, Type>::type GradType;
49 
50  const fvMesh& mesh = vsf.mesh();
51 
52  tmp<VolField<GradType>> tlsGrad
53  (
55  (
56  name,
57  mesh,
58  dimensioned<GradType>
59  (
60  "zero",
61  vsf.dimensions()/dimLength,
62  Zero
63  ),
64  extrapolatedCalculatedFvPatchField<GradType>::typeName
65  )
66  );
67  VolField<GradType>& lsGrad = tlsGrad.ref();
68 
69  // Get reference to least square vectors
70  const leastSquaresVectors& lsv = leastSquaresVectors::New(mesh);
71 
72  const surfaceVectorField& ownLs = lsv.pVectors();
73  const surfaceVectorField& neiLs = lsv.nVectors();
74 
75  const labelUList& own = mesh.owner();
76  const labelUList& nei = mesh.neighbour();
77 
78  forAll(own, facei)
79  {
80  label ownFacei = own[facei];
81  label neiFacei = nei[facei];
82 
83  Type deltaVsf = vsf[neiFacei] - vsf[ownFacei];
84 
85  lsGrad[ownFacei] += ownLs[facei]*deltaVsf;
86  lsGrad[neiFacei] -= neiLs[facei]*deltaVsf;
87  }
88 
89  // Boundary faces
90  forAll(vsf.boundaryField(), patchi)
91  {
92  const fvsPatchVectorField& patchOwnLs = ownLs.boundaryField()[patchi];
93 
94  const labelUList& faceCells =
95  vsf.boundaryField()[patchi].patch().faceCells();
96 
97  if (vsf.boundaryField()[patchi].coupled())
98  {
99  const Field<Type> neiVsf
100  (
101  vsf.boundaryField()[patchi].patchNeighbourField()
102  );
103 
104  forAll(neiVsf, patchFacei)
105  {
106  lsGrad[faceCells[patchFacei]] +=
107  patchOwnLs[patchFacei]
108  *(neiVsf[patchFacei] - vsf[faceCells[patchFacei]]);
109  }
110  }
111  else
112  {
113  const fvPatchField<Type>& patchVsf = vsf.boundaryField()[patchi];
114 
115  forAll(patchVsf, patchFacei)
116  {
117  lsGrad[faceCells[patchFacei]] +=
118  patchOwnLs[patchFacei]
119  *(patchVsf[patchFacei] - vsf[faceCells[patchFacei]]);
120  }
121  }
122  }
123 
124 
125  lsGrad.correctBoundaryConditions();
127 
128  return tlsGrad;
129 }
130 
131 
132 // ************************************************************************* //
#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
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
SurfaceField< vector > surfaceVectorField
UList< label > labelUList
Definition: UList.H:65
fvsPatchField< vector > fvsPatchVectorField
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488