Loading...
Searching...
No Matches
skewCorrectedSnGrad.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) 2019 Zeljko Tukovic, FSB Zagreb.
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 "skewCorrectedSnGrad.H"
30#include "volFields.H"
31#include "surfaceFields.H"
32#include "linear.H"
33#include "fvcGrad.H"
34#include "gaussGrad.H"
35
36// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
37
38template<class Type>
41(
43) const
44{
45 const fvMesh& mesh = this->mesh();
46
48 (
50 (
51 "snGradCorr("+vf.name()+')',
52 vf.instance(),
53 mesh,
56 ),
57 mesh,
58 vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
59 );
60 auto& ssf = tssf.ref();
61
62 ssf.setOriented();
63 ssf = Zero;
64
65
66 typedef typename
68 CmptGradType;
69
70 const labelUList& owner = mesh.owner();
71 const labelUList& neighbour = mesh.neighbour();
72
73 const vectorField& Sf = mesh.Sf().internalField();
74 const scalarField& magSf = mesh.magSf().internalField();
75
76 const vectorField& Cf = mesh.Cf().internalField();
77 const vectorField& C = mesh.C().internalField();
78
80 mesh.deltaCoeffs().internalField();
81
83 (
85 (
86 "kP",
87 vf.instance(),
88 mesh,
91 ),
92 mesh,
94 );
96
98 (
100 (
101 "kN",
102 vf.instance(),
103 mesh,
106 ),
107 mesh,
109 );
110 vectorField& kNI = kN.primitiveFieldRef();
111
112 kPI = Cf - vectorField(C, owner);
113 kPI -= Sf*(Sf&kPI)/sqr(magSf);
114
115 kNI = Cf - vectorField(C, neighbour);
116 kNI -= Sf*(Sf&kNI)/sqr(magSf);
117
118 forAll(kP.boundaryField(), patchI)
119 {
120 if (kP.boundaryField()[patchI].coupled())
121 {
122 kP.boundaryFieldRef()[patchI] =
123 mesh.boundary()[patchI].Cf()
124 - mesh.boundary()[patchI].Cn();
125
126 kP.boundaryFieldRef()[patchI] -=
127 mesh.boundary()[patchI].Sf()
128 *(
129 mesh.boundary()[patchI].Sf()
130 & kP.boundaryField()[patchI]
131 )
132 /sqr(mesh.boundary()[patchI].magSf());
133
134 kN.boundaryFieldRef()[patchI] =
135 mesh.Cf().boundaryField()[patchI]
136 - (
137 mesh.boundary()[patchI].Cn()
138 + mesh.boundary()[patchI].delta()
139 );
140
141 kN.boundaryFieldRef()[patchI] -=
142 mesh.boundary()[patchI].Sf()
143 *(
144 mesh.boundary()[patchI].Sf()
145 & kN.boundaryField()[patchI]
146 )
147 /sqr(mesh.boundary()[patchI].magSf());
148 }
149 }
150
151 for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; ++cmpt)
152 {
154 (
156 (
157 mesh,
158 mesh.gradScheme("grad(" + vf.name() + ')')
159 )()
160 .grad(vf.component(cmpt))
161 );
162
163 const Field<CmptGradType>& cmptGradI = cmptGrad.internalField();
164
165 // Skewness and non-rothogonal correction
166 {
167 ssf.primitiveFieldRef().replace
168 (
169 cmpt,
170 (
171 (kNI&Field<CmptGradType>(cmptGradI, neighbour))
172 - (kPI&Field<CmptGradType>(cmptGradI, owner))
173 )
175 );
176 }
177
178 forAll(ssf.boundaryField(), patchI)
179 {
180 if (ssf.boundaryField()[patchI].coupled())
181 {
182 ssf.boundaryFieldRef()[patchI].replace
183 (
184 cmpt,
185 (
186 (
187 kN.boundaryField()[patchI]
188 & cmptGrad.boundaryField()[patchI]
189 .patchNeighbourField()
190 )
191 - (
192 kP.boundaryField()[patchI]
193 & cmptGrad.boundaryField()[patchI]
194 .patchInternalField()
195 )
196 )
197 *mesh.deltaCoeffs().boundaryField()[patchI]
198 );
199 }
200 }
201 }
202
203 // // construct GeometricField<Type, fvsPatchField, surfaceMesh>
204 // tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tssf =
205 // linear<typename outerProduct<vector, Type>::type>(mesh).dotInterpolate
206 // (
207 // mesh.nonOrthCorrectionVectors(),
208 // gradScheme<Type>::New
209 // (
210 // mesh,
211 // mesh.gradScheme("grad(" + vf.name() + ')')
212 // )().grad(vf, "grad(" + vf.name() + ')')
213 // );
214 //
215 // tssf.ref().rename("snGradCorr(" + vf.name() + ')');
217 return tssf;
218}
219
220
221template<class Type>
224(
226) const
227{
228 const fvMesh& mesh = this->mesh();
229
231 (
233 (
234 "snGradCorr("+vf.name()+')',
235 vf.instance(),
236 mesh,
239 ),
240 mesh,
241 vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
242 );
243 auto& ssf = tssf.ref();
244 ssf.setOriented();
245
246 for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; ++cmpt)
247 {
248 ssf.replace
249 (
250 cmpt,
252 .fullGradCorrection(vf.component(cmpt))
253 );
254 }
255
256 return tssf;
257}
258
259
260// ************************************************************************* //
Graphite solid properties.
Definition C.H:49
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.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the 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
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class for gradient schemes.
Definition gradScheme.H:62
Surface gradient scheme with skewness and full explicit non-orthogonal corrections.
skewCorrectedSnGrad(const fvMesh &mesh)
Construct from mesh.
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the skewCorrectedSnGrad for the given field using the gradients of ...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fullGradCorrection(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the skewCorrectedSnGrad for the given field using the gradient of t...
const fvMesh & mesh() const
Return const reference to mesh.
static tmp< snGradScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
virtual const word & type() const =0
Runtime type information.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A class for managing temporary objects.
Definition tmp.H:75
dynamicFvMesh & mesh
Calculate the gradient of the given field.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.