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 >
49>
51(
53 const word& name
54) const
55{
56 // The fourth-order gradient is calculated in two passes. First,
57 // the standard least-square gradient is assembled. Then, the
58 // fourth-order correction is added to the second-order accurate
59 // gradient to complete the accuracy.
60
61 typedef typename outerProduct<vector, Type>::type GradType;
63
64 const fvMesh& mesh = vsf.mesh();
65
66 // Assemble the second-order least-square gradient
67 // Calculate the second-order least-square gradient
68 tmp<GradFieldType> tsecondfGrad
70 (
71 vsf,
72 "leastSquaresGrad(" + vsf.name() + ")"
73 );
74 const GradFieldType& secondfGrad =
75 tsecondfGrad();
76
78 (
79 new GradFieldType
80 (
82 (
83 name,
84 vsf.instance(),
85 mesh,
86 IOobject::NO_READ,
87 IOobject::NO_WRITE
88 ),
89 secondfGrad
90 )
91 );
92 GradFieldType& fGrad = tfGrad.ref();
93
94 const vectorField& C = mesh.C();
95
96 const surfaceScalarField& lambda = mesh.weights();
97
98 // Get reference to least square vectors
99 const leastSquaresVectors& lsv = leastSquaresVectors::New(mesh);
100 const surfaceVectorField& ownLs = lsv.pVectors();
101 const surfaceVectorField& neiLs = lsv.nVectors();
102
103 // owner/neighbour addressing
104 const labelUList& own = mesh.owner();
105 const labelUList& nei = mesh.neighbour();
106
107 // Assemble the fourth-order gradient
108
109 // Internal faces
110 forAll(own, facei)
111 {
112 Type dDotGradDelta = 0.5*
113 (
114 (C[nei[facei]] - C[own[facei]])
115 & (secondfGrad[nei[facei]] - secondfGrad[own[facei]])
116 );
117
118 fGrad[own[facei]] -= lambda[facei]*ownLs[facei]*dDotGradDelta;
119 fGrad[nei[facei]] -= (1.0 - lambda[facei])*neiLs[facei]*dDotGradDelta;
120 }
121
122 // Boundary faces
123 forAll(vsf.boundaryField(), patchi)
124 {
125 if (secondfGrad.boundaryField()[patchi].coupled())
126 {
127 const fvsPatchVectorField& patchOwnLs =
128 ownLs.boundaryField()[patchi];
129
130 const scalarField& lambdap = lambda.boundaryField()[patchi];
131
132 const fvPatch& p = vsf.boundaryField()[patchi].patch();
133
134 const labelUList& faceCells = p.faceCells();
135
136 // Build the d-vectors
137 const vectorField pd(p.delta());
138
139 const Field<GradType> neighbourSecondfGrad
140 (
141 secondfGrad.boundaryField()[patchi].patchNeighbourField()
142 );
143
144 forAll(faceCells, patchFacei)
145 {
146 fGrad[faceCells[patchFacei]] -=
147 0.5*lambdap[patchFacei]*patchOwnLs[patchFacei]
148 *(
149 pd[patchFacei]
150 & (
151 neighbourSecondfGrad[patchFacei]
152 - secondfGrad[faceCells[patchFacei]]
153 )
154 );
155 }
156 }
157 }
158
159 fGrad.correctBoundaryConditions();
161
162 return tfGrad;
163}
164
165
166// ************************************************************************* //
Graphite solid properties.
Definition: C.H:53
const Mesh & mesh() const
Return mesh.
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 word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Second-order gradient scheme using least-squares.
Definition: fourthGrad.H:63
Basic second-order gradient scheme using face-interpolation and Gauss' theorem.
Definition: gaussGrad.H:66
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvPatchField, volMesh > &, const word &name) const
Definition: gradScheme.C:88
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
volScalarField & p
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333