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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2018-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "leastSquaresGrad.H"
30#include "leastSquaresVectors.H"
31#include "gaussGrad.H"
32#include "fvMesh.H"
33#include "volMesh.H"
34#include "surfaceMesh.H"
35#include "GeometricField.H"
37
38// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39
40template<class Type>
42<
44 <
48 >
49>
51(
53 const word& name
54) const
55{
56 typedef typename outerProduct<vector, Type>::type GradType;
58
59 const fvMesh& mesh = vsf.mesh();
60
61 tmp<GradFieldType> tlsGrad
62 (
63 new GradFieldType
64 (
66 (
67 name,
68 vsf.instance(),
69 mesh,
70 IOobject::NO_READ,
71 IOobject::NO_WRITE
72 ),
73 mesh,
76 )
77 );
78 GradFieldType& lsGrad = tlsGrad.ref();
79
80 // Get reference to least square vectors
81 const leastSquaresVectors& lsv = leastSquaresVectors::New(mesh);
82
83 const surfaceVectorField& ownLs = lsv.pVectors();
84 const surfaceVectorField& neiLs = lsv.nVectors();
85
86 const labelUList& own = mesh.owner();
87 const labelUList& nei = mesh.neighbour();
88
89 forAll(own, facei)
90 {
91 const label ownFacei = own[facei];
92 const label neiFacei = nei[facei];
93
94 const Type deltaVsf(vsf[neiFacei] - vsf[ownFacei]);
95
96 lsGrad[ownFacei] += ownLs[facei]*deltaVsf;
97 lsGrad[neiFacei] -= neiLs[facei]*deltaVsf;
98 }
99
100 // Boundary faces
101 forAll(vsf.boundaryField(), patchi)
102 {
103 const fvsPatchVectorField& patchOwnLs = ownLs.boundaryField()[patchi];
104
105 const labelUList& faceCells =
106 vsf.boundaryField()[patchi].patch().faceCells();
107
108 if (vsf.boundaryField()[patchi].coupled())
109 {
110 const Field<Type> neiVsf
111 (
112 vsf.boundaryField()[patchi].patchNeighbourField()
113 );
114
115 forAll(neiVsf, patchFacei)
116 {
117 lsGrad[faceCells[patchFacei]] +=
118 patchOwnLs[patchFacei]
119 *(neiVsf[patchFacei] - vsf[faceCells[patchFacei]]);
120 }
121 }
122 else
123 {
124 const fvPatchField<Type>& patchVsf = vsf.boundaryField()[patchi];
125
126 forAll(patchVsf, patchFacei)
127 {
128 lsGrad[faceCells[patchFacei]] +=
129 patchOwnLs[patchFacei]
130 *(patchVsf[patchFacei] - vsf[faceCells[patchFacei]]);
131 }
132 }
133 }
134
135
136 lsGrad.correctBoundaryConditions();
138
139 return tlsGrad;
140}
141
142
143// ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
Generic dimensioned Type class.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
Basic second-order gradient scheme using face-interpolation and Gauss' theorem.
Definition: gaussGrad.H:66
Second-order gradient scheme using least-squares.
Least-squares gradient scheme vectors.
const surfaceVectorField & pVectors() const
Return const reference to owner least square vectors.
const surfaceVectorField & nVectors() const
Return const reference to neighbour least square vectors.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:54
type
Volume classification types.
Definition: volumeType.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333