Loading...
Searching...
No Matches
patchTransformedInterpolation.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) 2015 OpenFOAM Foundation
9 Copyright (C) 2015 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
31#include "pointFields.H"
32#include "symmTensor2D.H"
33#include "tensor2D.H"
34#include "syncTools.H"
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42
44 (
48 );
49}
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
53Foam::labelList Foam::patchTransformedInterpolation::getPatches
54(
55 Istream& entry
56) const
57{
58 wordList patchNames(entry);
59
60 labelList patches(patchNames.size(), -1);
61
62 forAll(patchNames, patchI)
63 {
64 patches[patchI] =
66 (
67 patchNames[patchI]
68 );
69
70 if (patches[patchI] == -1)
71 {
73 << "patch \"" << patchNames[patchI]
74 << "\" not found" << exit(FatalError);
75 }
76 }
77
78 return patches;
79}
80
81// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82
84(
85 const fvMesh& mesh,
86 Istream& entry
87)
88:
89 motionInterpolation(mesh, entry),
90 patches_(getPatches(entry))
91{}
92
93
94// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97{}
98
99
100// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
101
103(
104 const volScalarField&,
106) const
107{
109}
110
111
113(
114 const volVectorField& cellDisplacement,
115 pointVectorField& pointDisplacement
116) const
117{
118 const pointField& points(mesh().points());
119 const label nPoints(points.size());
120
122 (
123 cellDisplacement,
124 pointDisplacement
125 );
126
127 pointDisplacement.correctBoundaryConditions();
128
129 vectorField pointRotation(nPoints, Zero);
130 scalarField pointExpansion(nPoints, Zero);
131
132 labelList pointDisplacementNSum(nPoints, Zero);
133 vectorField pointDisplacementSum(nPoints, Zero);
134
135 for (const label patchi : patches_)
136 {
137 const polyPatch& patch = mesh().boundaryMesh()[patchi];
138
139 forAll(patch, pFaceI)
140 {
141 const face& f(patch[pFaceI]);
142
143 const label cellI(patch.faceCells()[pFaceI]);
144 const cell& c(mesh().cells()[cellI]);
145 const labelList cPoints(c.labels(mesh().faces()));
146
147 // Consider movement around the face centre
148 const point xOrigin(patch.faceCentres()[pFaceI]);
149
150 // Mean translation
151 const vector uMean(f.average(points, pointDisplacement));
152
153 // Calculate rotation and expansion for each point
154 forAll(f, fPointI)
155 {
156 const label pointI(f[fPointI]);
157 const vector& x(points[pointI]);
158 const vector r(x - xOrigin);
159 const vector u(pointDisplacement[pointI] - uMean);
160
161 pointRotation[pointI] = 2*(r ^ u)/magSqr(r);
162 pointExpansion[pointI] = (r & u)/magSqr(r);
163 }
164
165 // Mean rotation and expansion
166 const vector omegaMean(f.average(points, pointRotation));
167 const scalar divMean(f.average(points, pointExpansion));
168
169 // Apply mean solid body motion to all cell points
170 forAll(cPoints, cPointI)
171 {
172 const label pointI(cPoints[cPointI]);
173 const vector& x(points[pointI]);
174 const vector r(x - xOrigin);
175
176 pointDisplacementNSum[pointI] += 1;
177 pointDisplacementSum[pointI] +=
178 uMean + (omegaMean ^ r) + (divMean*r);
179 }
180 }
181 }
182
184 (
185 mesh(),
186 pointDisplacementNSum,
188 label(0)
189 );
190
192 (
193 mesh(),
194 pointDisplacementSum,
197 );
198
199 forAll(points, pointI)
200 {
201 if (pointDisplacementNSum[pointI])
202 {
203 pointDisplacement[pointI] =
204 pointDisplacementSum[pointI]/pointDisplacementNSum[pointI];
205 }
206 }
207
208 // Correct the faces
209 pointDisplacement.correctBoundaryConditions();
210}
211
212
213// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void correctBoundaryConditions()
Correct boundary field.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
static FOAM_NO_DANGLING_REFERENCE const volPointInterpolation & New(const fvMesh &mesh, Args &&... args)
static const Form zero
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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 of patch-adjacent ...
patchTransformedInterpolation(const fvMesh &mesh, Istream &entry)
Construct from an fvMesh and an Istream.
virtual void interpolate(const volScalarField &, pointScalarField &) const
Interpolate the given scalar cell displacement.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
const cellShapeList & cells
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
wordList patchNames(nPatches)
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299