OpenFOAM: API Guide
v2112
The open source CFD toolbox
deferredCorrection.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) 2017 OpenCFD Ltd.
9
-------------------------------------------------------------------------------
10
License
11
This file is part of OpenFOAM.
12
13
OpenFOAM is free software: you can redistribute it and/or modify it
14
under the terms of the GNU General Public License as published by
15
the Free Software Foundation, either version 3 of the License, or
16
(at your option) any later version.
17
18
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21
for more details.
22
23
You should have received a copy of the GNU General Public License
24
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26
\*---------------------------------------------------------------------------*/
27
28
#include "
deferredCorrection.H
"
29
30
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31
32
template
<
class
Type>
33
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>
>
34
Foam::deferredCorrection<Type>::correction
35
(
36
const
GeometricField<Type, fvPatchField, volMesh>
& vf
37
)
const
38
{
39
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>
> tsfCorr
40
(
41
new
GeometricField<Type, fvsPatchField, surfaceMesh>
42
(
43
"deferredCorrection::correction("
+ vf.
name
() +
')'
,
44
tbaseScheme_().
interpolate
(vf)
45
)
46
);
47
GeometricField<Type, fvsPatchField, surfaceMesh>
& sfCorr = tsfCorr.
ref
();
48
49
// Interpolate using the upwind weights to avoid circular reference to
50
// [this] explicit correction
51
sfCorr -=
upwind<Type>::interpolate
(vf,
upwind<Type>::weights
());
52
/*
53
auto& sfCorrBf = sfCorr.boundaryFieldRef();
54
for (auto& pf : sfCorrBf)
55
{
56
if (!pf.coupled())
57
{
58
pf = pTraits<Type>::zero;
59
}
60
}
61
*/
62
return
tsfCorr;
63
}
64
65
66
namespace
Foam
67
{
68
//makeSurfaceInterpolationScheme(deferredCorrection)
69
makeSurfaceInterpolationTypeScheme
(
deferredCorrection
, scalar)
70
makeSurfaceInterpolationTypeScheme
(
deferredCorrection
,
vector
)
71
}
72
73
// ************************************************************************* //
Foam::GeometricField
Generic GeometricField class.
Definition:
GeometricField.H:80
Foam::IOobject::name
const word & name() const noexcept
Return the object name.
Definition:
IOobjectI.H:65
Foam::Vector< scalar >
Foam::deferredCorrection
Deferred correction interpolation scheme wrapper around a run-time selectable base scheme.
Definition:
deferredCorrection.H:75
Foam::deferredCorrection::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
Definition:
deferredCorrection.C:35
Foam::tmp
A class for managing temporary objects.
Definition:
tmp.H:65
Foam::tmp::ref
T & ref() const
Definition:
tmpI.H:227
Foam::upwind
Upwind differencing scheme class.
Definition:
upwind.H:59
interpolate
mesh interpolate(rAU)
deferredCorrection.H
Foam
Namespace for OpenFOAM.
Definition:
atmBoundaryLayer.C:34
makeSurfaceInterpolationTypeScheme
#define makeSurfaceInterpolationTypeScheme(SS, Type)
Definition:
surfaceInterpolationScheme.H:286
src
finiteVolume
interpolation
surfaceInterpolation
schemes
deferredCorrection
deferredCorrection.C
Generated by
1.9.5
OPENFOAM® is a registered
trademark
of OpenCFD Ltd.