Loading...
Searching...
No Matches
fourthGrad.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) 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 "fourthGrad.H"
30#include "leastSquaresGrad.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 >& fGrad,
50 const GeometricField
51 <
55 >& secondfGrad
56) const
57{
58 // The fourth-order gradient is calculated in two passes. First,
59 // the standard least-square gradient is assembled. Then, the
60 // fourth-order correction is added to the second-order accurate
61 // gradient to complete the accuracy.
62
63 typedef typename outerProduct<vector, Type>::type GradType;
64
65 const fvMesh& mesh = vsf.mesh();
66
67 // Initialise to gradient
68 fGrad = secondfGrad;
69
70 const vectorField& C = mesh.C();
71
72 const surfaceScalarField& lambda = mesh.weights();
73
74 // Get reference to least square vectors
76 const surfaceVectorField& ownLs = lsv.pVectors();
77 const surfaceVectorField& neiLs = lsv.nVectors();
78
79 // owner/neighbour addressing
80 const labelUList& own = mesh.owner();
81 const labelUList& nei = mesh.neighbour();
82
83 // Assemble the fourth-order gradient
84
85 // Internal faces
86 forAll(own, facei)
87 {
88 Type dDotGradDelta = 0.5*
89 (
90 (C[nei[facei]] - C[own[facei]])
91 & (secondfGrad[nei[facei]] - secondfGrad[own[facei]])
92 );
93
94 fGrad[own[facei]] -= lambda[facei]*ownLs[facei]*dDotGradDelta;
95 fGrad[nei[facei]] -= (1.0 - lambda[facei])*neiLs[facei]*dDotGradDelta;
96 }
97
98 // Boundary faces
99 forAll(vsf.boundaryField(), patchi)
100 {
101 if (secondfGrad.boundaryField()[patchi].coupled())
102 {
103 const fvsPatchVectorField& patchOwnLs =
104 ownLs.boundaryField()[patchi];
105
106 const scalarField& lambdap = lambda.boundaryField()[patchi];
107
108 const fvPatch& p = vsf.boundaryField()[patchi].patch();
109
110 const labelUList& faceCells = p.faceCells();
111
112 // Build the d-vectors
113 const vectorField pd(p.delta());
114
115 const Field<GradType> neighbourSecondfGrad
116 (
117 secondfGrad.boundaryField()[patchi].patchNeighbourField()
118 );
119
120 forAll(faceCells, patchFacei)
121 {
122 fGrad[faceCells[patchFacei]] -=
123 0.5*lambdap[patchFacei]*patchOwnLs[patchFacei]
124 *(
125 pd[patchFacei]
126 & (
127 neighbourSecondfGrad[patchFacei]
128 - secondfGrad[faceCells[patchFacei]]
129 )
130 );
131 }
132 }
133 }
134
135 fGrad.correctBoundaryConditions();
137}
138
139
140template<class Type>
142<
144 <
148 >
149>
151(
153 const word& name
154) const
155{
156 // The fourth-order gradient is calculated in two passes. First,
157 // the standard least-square gradient is assembled. Then, the
158 // fourth-order correction is added to the second-order accurate
159 // gradient to complete the accuracy.
160
161 typedef typename outerProduct<vector, Type>::type GradType;
163
164 const fvMesh& mesh = vsf.mesh();
165
166 // Assemble the second-order least-square gradient
167 // Calculate the second-order least-square gradient
168 tmp<GradFieldType> tsecondfGrad
170 (
171 vsf,
172 "leastSquaresGrad(" + vsf.name() + ")"
173 );
174 const GradFieldType& secondfGrad = tsecondfGrad();
175
176 tmp<GradFieldType> tfGrad
177 (
178 new GradFieldType
179 (
181 (
182 name,
183 vsf.instance(),
184 mesh,
187 ),
188 mesh,
189 vsf.dimensions()/dimLength,
191 )
192 );
193
194 calcGrad(tfGrad.ref(), vsf, secondfGrad);
195
196 return tfGrad;
197}
198
199
200template<class Type>
202(
204 <
207 volMesh
208 >& fGrad,
210) const
211{
212 // The fourth-order gradient is calculated in two passes. First,
213 // the standard least-square gradient is assembled. Then, the
214 // fourth-order correction is added to the second-order accurate
215 // gradient to complete the accuracy.
216
217 typedef typename outerProduct<vector, Type>::type GradType;
219
220 const fvMesh& mesh = vsf.mesh();
221
222 // Assemble the second-order least-square gradient
223 // Calculate the second-order least-square gradient
224 tmp<GradFieldType> tsecondfGrad
226 (
227 vsf,
228 "leastSquaresGrad(" + vsf.name() + ")"
229 );
230
231 calcGrad(fGrad, vsf, tsecondfGrad());
232}
233
234
235// ************************************************************************* //
Graphite solid properties.
Definition C.H:49
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
static FOAM_NO_DANGLING_REFERENCE const leastSquaresVectors & New(const fvMesh &mesh, Args &&... args)
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
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 for optional caching.
Definition fourthGrad.C:144
static void correctBoundaryConditions(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > &)
Correct the boundary values of the gradient using the patchField snGrad functions.
Definition gaussGrad.C:199
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvPatchField, volMesh > &, const word &name) const
Calculate and return the grad of the given field which may have been cached.
Definition gradScheme.C:81
const fvMesh & mesh() const
Return const reference to mesh.
Definition gradScheme.H:138
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.
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvsPatchField< vector > fvsPatchVectorField
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299