Loading...
Searching...
No Matches
displacementInterpolationMotionSolver.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) 2011-2016 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 "SortableList.H"
32#include "GlobalIOList.H"
33#include "Tuple2.H"
34#include "mapPolyMesh.H"
35#include "interpolateXY.H"
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42
44 (
48 );
49
51 (
54 displacement
55 );
56
57 template<>
59 (
60 "scalarVectorTable"
61 );
62}
63
64
65// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
66
67void Foam::displacementInterpolationMotionSolver::calcInterpolation()
68{
69 // Get zones and their interpolation tables for displacement
70 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
71
72 List<Pair<word>> faceZoneToTable
73 (
74 coeffDict().lookup("interpolationTables")
75 );
76
77 const faceZoneMesh& fZones = mesh().faceZones();
78
79 times_.setSize(fZones.size());
80 displacements_.setSize(fZones.size());
81
82 forAll(faceZoneToTable, i)
83 {
84 const word& zoneName = faceZoneToTable[i][0];
85 label zoneI = fZones.findZoneID(zoneName);
86
87 if (zoneI == -1)
88 {
90 << "Cannot find zone " << zoneName << endl
91 << "Valid zones are " << fZones.names()
92 << exit(FatalError);
93 }
94
95 const word& tableName = faceZoneToTable[i][1];
96
97 GlobalIOList<Tuple2<scalar, vector>> table
98 (
99 IOobject
100 (
101 tableName,
102 mesh().time().constant(),
103 "tables",
104 mesh(),
108 )
109 );
110
111 // Copy table
112 times_[zoneI].setSize(table.size());
113 displacements_[zoneI].setSize(table.size());
114
115 forAll(table, j)
116 {
117 times_[zoneI][j] = table[j].first();
118 displacements_[zoneI][j] = table[j].second();
119 }
120 }
121
122
123
124 // Sort points into bins according to position relative to faceZones
125 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
126 // Done in all three directions.
127
128 for (direction dir = 0; dir < vector::nComponents; dir++)
129 {
130 // min and max coordinates of all faceZones
131 SortableList<scalar> zoneCoordinates(2*faceZoneToTable.size());
132
133 forAll(faceZoneToTable, i)
134 {
135 const word& zoneName = faceZoneToTable[i][0];
136 const faceZone& fz = fZones[zoneName];
137
139
140 for (const label pointi : fz().meshPoints())
141 {
142 const scalar coord = points0()[pointi][dir];
143 limits.add(coord);
144 }
145
146 reduce(limits, sumOp<scalarMinMax>());
147
148 zoneCoordinates[2*i] = limits.min();
149 zoneCoordinates[2*i+1] = limits.max();
150
151 if (debug)
152 {
153 Pout<< "direction " << dir << " : "
154 << "zone " << zoneName
155 << " ranges from coordinate "
156 << limits.min() << " to " << limits.max()
157 << endl;
158 }
159 }
160 zoneCoordinates.sort();
161
162 // Slightly tweak min and max face zone so points sort within
163 zoneCoordinates[0] -= SMALL;
164 zoneCoordinates.last() += SMALL;
165
166 // Check if we have static min and max mesh bounds
167 const scalarField meshCoords(points0().component(dir));
168
169 scalar minCoord, maxCoord;
170 {
171 auto limits = gMinMax(meshCoords);
172
173 minCoord = limits.min();
174 maxCoord = limits.max();
175 }
176
177 if (debug)
178 {
179 Pout<< "direction " << dir << " : "
180 << "mesh ranges from coordinate "
181 << minCoord << " to " << maxCoord
182 << endl;
183 }
184
185 // Make copy of zoneCoordinates; include min and max of mesh
186 // if necessary. Mark min and max with zoneI=-1.
187 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
188
189 labelList& rangeZone = rangeToZone_[dir];
190 labelListList& rangePoints = rangeToPoints_[dir];
191 List<scalarField>& rangeWeights = rangeToWeights_[dir];
192
193 scalarField rangeToCoord(zoneCoordinates.size());
194 rangeZone.setSize(zoneCoordinates.size());
195 label rangeI = 0;
196
197 if (minCoord < zoneCoordinates[0])
198 {
199 label sz = rangeZone.size();
200 rangeToCoord.setSize(sz+1);
201 rangeZone.setSize(sz+1);
202 rangeToCoord[rangeI] = minCoord-SMALL;
203 rangeZone[rangeI] = -1;
204
205 if (debug)
206 {
207 Pout<< "direction " << dir << " : "
208 << "range " << rangeI << " at coordinate "
209 << rangeToCoord[rangeI] << " from min of mesh "
210 << rangeZone[rangeI] << endl;
211 }
212 rangeI = 1;
213 }
214 forAll(zoneCoordinates, i)
215 {
216 rangeToCoord[rangeI] = zoneCoordinates[i];
217 rangeZone[rangeI] = zoneCoordinates.indices()[i]/2;
218
219 if (debug)
220 {
221 Pout<< "direction " << dir << " : "
222 << "range " << rangeI << " at coordinate "
223 << rangeToCoord[rangeI]
224 << " from zone " << rangeZone[rangeI] << endl;
225 }
226 rangeI++;
227 }
228 if (maxCoord > zoneCoordinates.last())
229 {
230 label sz = rangeToCoord.size();
231 rangeToCoord.setSize(sz+1);
232 rangeZone.setSize(sz+1);
233 rangeToCoord[sz] = maxCoord+SMALL;
234 rangeZone[sz] = -1;
235
236 if (debug)
237 {
238 Pout<< "direction " << dir << " : "
239 << "range " << rangeI << " at coordinate "
240 << rangeToCoord[sz] << " from max of mesh "
241 << rangeZone[sz] << endl;
242 }
243 }
244
245
246 // Sort the points
247 // ~~~~~~~~~~~~~~~
248
249 // Count all the points inbetween rangeI and rangeI+1
250 labelList nRangePoints(rangeToCoord.size(), Zero);
251
252 forAll(meshCoords, pointi)
253 {
254 label rangeI = findLower(rangeToCoord, meshCoords[pointi]);
255
256 if (rangeI == -1 || rangeI == rangeToCoord.size()-1)
257 {
259 << "Did not find point " << points0()[pointi]
260 << " coordinate " << meshCoords[pointi]
261 << " in ranges " << rangeToCoord
262 << abort(FatalError);
263 }
264 nRangePoints[rangeI]++;
265 }
266
267 if (debug)
268 {
269 for (label rangeI = 0; rangeI < rangeToCoord.size()-1; rangeI++)
270 {
271 // Get the two zones bounding the range
272 Pout<< "direction " << dir << " : "
273 << "range from " << rangeToCoord[rangeI]
274 << " to " << rangeToCoord[rangeI+1]
275 << " contains " << nRangePoints[rangeI]
276 << " points." << endl;
277 }
278 }
279
280 // Sort
281 rangePoints.setSize(nRangePoints.size());
282 rangeWeights.setSize(nRangePoints.size());
283 forAll(rangePoints, rangeI)
284 {
285 rangePoints[rangeI].setSize(nRangePoints[rangeI]);
286 rangeWeights[rangeI].setSize(nRangePoints[rangeI]);
287 }
288 nRangePoints = 0;
289 forAll(meshCoords, pointi)
290 {
291 label rangeI = findLower(rangeToCoord, meshCoords[pointi]);
292 label& nPoints = nRangePoints[rangeI];
293 rangePoints[rangeI][nPoints] = pointi;
294 rangeWeights[rangeI][nPoints] =
295 (meshCoords[pointi]-rangeToCoord[rangeI])
296 / (rangeToCoord[rangeI+1]-rangeToCoord[rangeI]);
297 nPoints++;
299 }
300}
301
302
303// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
304
305Foam::displacementInterpolationMotionSolver::
306displacementInterpolationMotionSolver
307(
308 const polyMesh& mesh,
309 const IOdictionary& dict
310)
313{
314 calcInterpolation();
315}
316
317
318Foam::displacementInterpolationMotionSolver::
319displacementInterpolationMotionSolver
320(
321 const polyMesh& mesh,
322 const IOdictionary& dict,
323 const pointVectorField& pointDisplacement,
324 const pointIOField& points0
325)
326:
327 displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName)
328{
329 calcInterpolation();
330}
331
332
333// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
334
337{
338 if (mesh().nPoints() != points0().size())
339 {
341 << "The number of points in the mesh seems to have changed." << endl
342 << "In constant/polyMesh there are " << points0().size()
343 << " points; in the current mesh there are " << mesh().nPoints()
344 << " points." << exit(FatalError);
345 }
346
347 auto tcurPoints = tmp<pointField>::New(points0());
348 auto& curPoints = tcurPoints.ref();
349
350 // Interpolate the displacement of the face zones.
351 vectorField zoneDisp(displacements_.size(), Zero);
352 forAll(zoneDisp, zoneI)
353 {
354 if (times_[zoneI].size())
355 {
356 zoneDisp[zoneI] = interpolateXY
357 (
358 mesh().time().value(),
359 times_[zoneI],
360 displacements_[zoneI]
361 );
362 }
363 }
364 if (debug)
365 {
366 Pout<< "Zone displacements:" << zoneDisp << endl;
367 }
368
369
370 // Interpolate the point location
371 for (direction dir = 0; dir < vector::nComponents; dir++)
372 {
373 const labelList& rangeZone = rangeToZone_[dir];
374 const labelListList& rangePoints = rangeToPoints_[dir];
375 const List<scalarField>& rangeWeights = rangeToWeights_[dir];
376
377 for (label rangeI = 0; rangeI < rangeZone.size()-1; rangeI++)
378 {
379 const labelList& rPoints = rangePoints[rangeI];
380 const scalarField& rWeights = rangeWeights[rangeI];
381
382 // Get the two zones bounding the range
383 label minZoneI = rangeZone[rangeI];
384 //vector minDisp =
385 // (minZoneI == -1 ? vector::zero : zoneDisp[minZoneI]);
386 scalar minDisp = (minZoneI == -1 ? 0.0 : zoneDisp[minZoneI][dir]);
387 label maxZoneI = rangeZone[rangeI+1];
388 //vector maxDisp =
389 // (maxZoneI == -1 ? vector::zero : zoneDisp[maxZoneI]);
390 scalar maxDisp = (maxZoneI == -1 ? 0.0 : zoneDisp[maxZoneI][dir]);
391
392 forAll(rPoints, i)
393 {
394 label pointi = rPoints[i];
395 scalar w = rWeights[i];
396 //curPoints[pointi] += (1.0-w)*minDisp+w*maxDisp;
397 curPoints[pointi][dir] += (1.0-w)*minDisp+w*maxDisp;
398 }
399 }
400 }
401 return tcurPoints;
402}
403
404
405// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
The number of elements in list.
Definition DLListBase.H:194
IOList with global data (so optionally read from master).
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static constexpr direction nComponents
Number of components in this vector space.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
Definition dictionary.C:367
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Virtual base class for displacement motion solver.
pointVectorField & pointDisplacement() noexcept
Return reference to the point motion displacement field.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition polyMesh.H:671
label nPoints() const noexcept
Number of mesh points.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
label nPoints
Interpolates y values from one curve to another with a different x distribution.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
vectorIOField pointIOField
pointIOField is a vectorIOField.
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
Field< Type > interpolateXY(const scalarField &xNew, const scalarField &xOld, const Field< Type > &yOld)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
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...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))