patchCorrectedInterpolation.H
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) 2015 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::patchCorrectedInterpolation
28
29Description
30 Interpolation of cell-based displacements to the points with additional
31 correction for interpolation inconsistency on patches.
32
33 The default interpolation method interpolates from the cells to all points
34 except those on boundaries with value boundary conditions. The discrepancy
35 across the boundary cell can cause shearing and inversion if the cells are
36 of very high aspect ratio.
37
38 This method applies the default interpolation to *all* points, including
39 those on value boundaries. The difference between the interpolated value on
40 the boundary and the desired boundary condition is then propagated into the
41 mesh with a wave. Contributions from different patches are inverse-distance
42 weighted, and the result is added to the default interpolation. The result
43 of this is that thin boundary cells are maintained much more accurately for
44 non-uniform patch displacements.
45
46 The user must specify the patch groups from which to propagate the motion.
47 Ideally, these groups will be opposing; i.e. one group with the aerofoil,
48 and one with the far field:
49
50 \verbatim
51 interpolation patchCorrected
52 (
53 (aerofoilUpper aerofoilLower)
54 (farField)
55 );
56 \endverbatim
57
58SourceFiles
59 patchCorrectedInterpolation.C
60
61\*---------------------------------------------------------------------------*/
62
63#ifndef patchCorrectedInterpolation_H
64#define patchCorrectedInterpolation_H
65
66#include "motionInterpolation.H"
67
68// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69
70namespace Foam
71{
72
73/*---------------------------------------------------------------------------*\
74 Class patchCorrectedInterpolation Declaration
75\*---------------------------------------------------------------------------*/
78:
80{
81 // Private data
82
83 //- Patch groups from which to propagate the displacement
84 const labelListList patchGroups_;
85
86
87 // Private member functions
88
89 //- Get patch groups from the input stream
90 labelListList getPatchGroups(Istream& entry) const;
91
92 //- Interpolate the given cell displacement
93 template <class Type>
94 void interpolateType
95 (
98 ) const;
99
100 //- Interpolate patch data by inverse distance weighting
101 template <class Type>
102 void interpolateDataFromPatchGroups
103 (
105 ) const;
106
107 //- Propagate data from a number of patches into the field
108 template <class Type>
109 void propagateDataFromPatchGroup
110 (
111 const label,
114 ) const;
115
116
117public:
118
119 //- Runtime type information
120 TypeName("patchCorrected");
121
122
123 // Constructors
124
125 //- Construct from an fvMesh and an Istream
127 (
128 const fvMesh& mesh,
130 );
131
132
133 //- Destructor
135
136
137 // Member Functions
138
139 //- Interpolate the given scalar cell displacement
140 virtual void interpolate
141 (
142 const volScalarField&,
144 ) const;
145
146 //- Interpolate the given vector cell displacement
147 virtual void interpolate
148 (
149 const volVectorField&,
151 ) const;
152};
153
154
155// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156
157} // End namespace Foam
158
159// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160
162
163// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164
165#endif
166
167// ************************************************************************* //
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
const fvMesh & mesh() const
Return const-reference to the mesh.
Interpolation of cell-based displacements to the points with additional correction for interpolation ...
virtual void interpolate(const volScalarField &, pointScalarField &) const
Interpolate the given scalar cell displacement.
TypeName("patchCorrected")
Runtime type information.
Namespace for OpenFOAM.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73