/
faceCorrectedSnGrad.C
167 lines (133 loc) · 4.95 KB
/
faceCorrectedSnGrad.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "faceCorrectedSnGrad.H"
#include "volPointInterpolation.H"
#include "triangle.H"
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
template<class Type>
Foam::fv::faceCorrectedSnGrad<Type>::~faceCorrectedSnGrad()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>>
Foam::fv::faceCorrectedSnGrad<Type>::fullGradCorrection
(
const GeometricField<Type, fvPatchField, volMesh>& vf
) const
{
const fvMesh& mesh = this->mesh();
GeometricField<Type, pointPatchField, pointMesh> pvf
(
volPointInterpolation::New(mesh).interpolate(vf)
);
// construct GeometricField<Type, fvsPatchField, surfaceMesh>
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsfCorr
(
new GeometricField<Type, fvsPatchField, surfaceMesh>
(
IOobject
(
"snGradCorr("+vf.name()+')',
vf.instance(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
)
);
Field<Type>& sfCorr = tsfCorr.ref().primitiveFieldRef();
const pointField& points = mesh.points();
const faceList& faces = mesh.faces();
const vectorField& Sf = mesh.Sf();
const vectorField& C = mesh.C();
const scalarField& magSf = mesh.magSf();
const labelList& owner = mesh.owner();
const labelList& neighbour = mesh.neighbour();
forAll(sfCorr, facei)
{
typename outerProduct<vector, Type>::type fgrad
(
outerProduct<vector, Type>::type::zero
);
const face& fi = faces[facei];
vector nf(Sf[facei]/magSf[facei]);
for (label pi=0; pi<fi.size(); pi++)
{
// Next point index
label pj = (pi+1)%fi.size();
// Edge normal in plane of face
vector edgen(nf^(points[fi[pj]] - points[fi[pi]]));
// Edge centre field value
Type pvfe(0.5*(pvf[fi[pj]] + pvf[fi[pi]]));
// Integrate face gradient
fgrad += edgen*pvfe;
}
// Finalize face-gradient by dividing by face area
fgrad /= magSf[facei];
// Calculate correction vector
vector dCorr(C[neighbour[facei]] - C[owner[facei]]);
dCorr /= (nf & dCorr);
// if (mag(dCorr) > 2) dCorr *= 2/mag(dCorr);
sfCorr[facei] = dCorr&fgrad;
}
tsfCorr.ref().boundaryFieldRef() = Zero;
return tsfCorr;
}
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>>
Foam::fv::faceCorrectedSnGrad<Type>::correction
(
const GeometricField<Type, fvPatchField, volMesh>& vf
) const
{
const fvMesh& mesh = this->mesh();
// construct GeometricField<Type, fvsPatchField, surfaceMesh>
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tssf
(
new GeometricField<Type, fvsPatchField, surfaceMesh>
(
IOobject
(
"snGradCorr("+vf.name()+')',
vf.instance(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
vf.dimensions()*mesh.nonOrthDeltaCoeffs().dimensions()
)
);
GeometricField<Type, fvsPatchField, surfaceMesh>& ssf = tssf.ref();
for (direction cmpt = 0; cmpt < pTraits<Type>::nComponents; cmpt++)
{
ssf.replace
(
cmpt,
faceCorrectedSnGrad<typename pTraits<Type>::cmptType>(mesh)
.fullGradCorrection(vf.component(cmpt))
);
}
return tssf;
}
// ************************************************************************* //