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 -------------------------------------------------------------------------------
11 License
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"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(displacementInterpolationMotionSolver, 0);
42 
44  (
45  motionSolver,
46  displacementInterpolationMotionSolver,
47  dictionary
48  );
49 
51  (
52  displacementMotionSolver,
53  displacementInterpolationMotionSolver,
54  displacement
55  );
56 
57  template<>
58  const word GlobalIOList<Tuple2<scalar, vector>>::typeName
59  (
60  "scalarVectorTable"
61  );
62 }
63 
64 
65 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
66 
67 void 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(),
107  false
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 
138  scalar minCoord = VGREAT;
139  scalar maxCoord = -VGREAT;
140 
141  forAll(fz().meshPoints(), localI)
142  {
143  label pointi = fz().meshPoints()[localI];
144  const scalar coord = points0()[pointi][dir];
145  minCoord = min(minCoord, coord);
146  maxCoord = max(maxCoord, coord);
147  }
148 
149  zoneCoordinates[2*i] = returnReduce(minCoord, minOp<scalar>());
150  zoneCoordinates[2*i+1] = returnReduce(maxCoord, maxOp<scalar>());
151 
152  if (debug)
153  {
154  Pout<< "direction " << dir << " : "
155  << "zone " << zoneName
156  << " ranges from coordinate " << zoneCoordinates[2*i]
157  << " to " << zoneCoordinates[2*i+1]
158  << endl;
159  }
160  }
161  zoneCoordinates.sort();
162 
163  // Slightly tweak min and max face zone so points sort within
164  zoneCoordinates[0] -= SMALL;
165  zoneCoordinates.last() += SMALL;
166 
167  // Check if we have static min and max mesh bounds
168  const scalarField meshCoords(points0().component(dir));
169 
170  scalar minCoord = gMin(meshCoords);
171  scalar maxCoord = gMax(meshCoords);
172 
173  if (debug)
174  {
175  Pout<< "direction " << dir << " : "
176  << "mesh ranges from coordinate " << minCoord << " to "
177  << maxCoord << endl;
178  }
179 
180  // Make copy of zoneCoordinates; include min and max of mesh
181  // if necessary. Mark min and max with zoneI=-1.
182  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
183 
184  labelList& rangeZone = rangeToZone_[dir];
185  labelListList& rangePoints = rangeToPoints_[dir];
186  List<scalarField>& rangeWeights = rangeToWeights_[dir];
187 
188  scalarField rangeToCoord(zoneCoordinates.size());
189  rangeZone.setSize(zoneCoordinates.size());
190  label rangeI = 0;
191 
192  if (minCoord < zoneCoordinates[0])
193  {
194  label sz = rangeZone.size();
195  rangeToCoord.setSize(sz+1);
196  rangeZone.setSize(sz+1);
197  rangeToCoord[rangeI] = minCoord-SMALL;
198  rangeZone[rangeI] = -1;
199 
200  if (debug)
201  {
202  Pout<< "direction " << dir << " : "
203  << "range " << rangeI << " at coordinate "
204  << rangeToCoord[rangeI] << " from min of mesh "
205  << rangeZone[rangeI] << endl;
206  }
207  rangeI = 1;
208  }
209  forAll(zoneCoordinates, i)
210  {
211  rangeToCoord[rangeI] = zoneCoordinates[i];
212  rangeZone[rangeI] = zoneCoordinates.indices()[i]/2;
213 
214  if (debug)
215  {
216  Pout<< "direction " << dir << " : "
217  << "range " << rangeI << " at coordinate "
218  << rangeToCoord[rangeI]
219  << " from zone " << rangeZone[rangeI] << endl;
220  }
221  rangeI++;
222  }
223  if (maxCoord > zoneCoordinates.last())
224  {
225  label sz = rangeToCoord.size();
226  rangeToCoord.setSize(sz+1);
227  rangeZone.setSize(sz+1);
228  rangeToCoord[sz] = maxCoord+SMALL;
229  rangeZone[sz] = -1;
230 
231  if (debug)
232  {
233  Pout<< "direction " << dir << " : "
234  << "range " << rangeI << " at coordinate "
235  << rangeToCoord[sz] << " from max of mesh "
236  << rangeZone[sz] << endl;
237  }
238  }
239 
240 
241  // Sort the points
242  // ~~~~~~~~~~~~~~~
243 
244  // Count all the points inbetween rangeI and rangeI+1
245  labelList nRangePoints(rangeToCoord.size(), Zero);
246 
247  forAll(meshCoords, pointi)
248  {
249  label rangeI = findLower(rangeToCoord, meshCoords[pointi]);
250 
251  if (rangeI == -1 || rangeI == rangeToCoord.size()-1)
252  {
254  << "Did not find point " << points0()[pointi]
255  << " coordinate " << meshCoords[pointi]
256  << " in ranges " << rangeToCoord
257  << abort(FatalError);
258  }
259  nRangePoints[rangeI]++;
260  }
261 
262  if (debug)
263  {
264  for (label rangeI = 0; rangeI < rangeToCoord.size()-1; rangeI++)
265  {
266  // Get the two zones bounding the range
267  Pout<< "direction " << dir << " : "
268  << "range from " << rangeToCoord[rangeI]
269  << " to " << rangeToCoord[rangeI+1]
270  << " contains " << nRangePoints[rangeI]
271  << " points." << endl;
272  }
273  }
274 
275  // Sort
276  rangePoints.setSize(nRangePoints.size());
277  rangeWeights.setSize(nRangePoints.size());
278  forAll(rangePoints, rangeI)
279  {
280  rangePoints[rangeI].setSize(nRangePoints[rangeI]);
281  rangeWeights[rangeI].setSize(nRangePoints[rangeI]);
282  }
283  nRangePoints = 0;
284  forAll(meshCoords, pointi)
285  {
286  label rangeI = findLower(rangeToCoord, meshCoords[pointi]);
287  label& nPoints = nRangePoints[rangeI];
288  rangePoints[rangeI][nPoints] = pointi;
289  rangeWeights[rangeI][nPoints] =
290  (meshCoords[pointi]-rangeToCoord[rangeI])
291  / (rangeToCoord[rangeI+1]-rangeToCoord[rangeI]);
292  nPoints++;
293  }
294  }
295 }
296 
297 
298 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
299 
300 Foam::displacementInterpolationMotionSolver::
301 displacementInterpolationMotionSolver
302 (
303  const polyMesh& mesh,
304  const IOdictionary& dict
305 )
306 :
308 {
309  calcInterpolation();
310 }
311 
312 
313 Foam::displacementInterpolationMotionSolver::
314 displacementInterpolationMotionSolver
315 (
316  const polyMesh& mesh,
317  const IOdictionary& dict,
318  const pointVectorField& pointDisplacement,
319  const pointIOField& points0
320 )
321 :
322  displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName)
323 {
324  calcInterpolation();
325 }
326 
327 
328 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
329 
332 {
333  if (mesh().nPoints() != points0().size())
334  {
336  << "The number of points in the mesh seems to have changed." << endl
337  << "In constant/polyMesh there are " << points0().size()
338  << " points; in the current mesh there are " << mesh().nPoints()
339  << " points." << exit(FatalError);
340  }
341 
342  tmp<pointField> tcurPoints(new pointField(points0()));
343  pointField& curPoints = tcurPoints.ref();
344 
345  // Interpolate the displacement of the face zones.
346  vectorField zoneDisp(displacements_.size(), Zero);
347  forAll(zoneDisp, zoneI)
348  {
349  if (times_[zoneI].size())
350  {
351  zoneDisp[zoneI] = interpolateXY
352  (
353  mesh().time().value(),
354  times_[zoneI],
355  displacements_[zoneI]
356  );
357  }
358  }
359  if (debug)
360  {
361  Pout<< "Zone displacements:" << zoneDisp << endl;
362  }
363 
364 
365  // Interpolate the point location
366  for (direction dir = 0; dir < vector::nComponents; dir++)
367  {
368  const labelList& rangeZone = rangeToZone_[dir];
369  const labelListList& rangePoints = rangeToPoints_[dir];
370  const List<scalarField>& rangeWeights = rangeToWeights_[dir];
371 
372  for (label rangeI = 0; rangeI < rangeZone.size()-1; rangeI++)
373  {
374  const labelList& rPoints = rangePoints[rangeI];
375  const scalarField& rWeights = rangeWeights[rangeI];
376 
377  // Get the two zones bounding the range
378  label minZoneI = rangeZone[rangeI];
379  //vector minDisp =
380  // (minZoneI == -1 ? vector::zero : zoneDisp[minZoneI]);
381  scalar minDisp = (minZoneI == -1 ? 0.0 : zoneDisp[minZoneI][dir]);
382  label maxZoneI = rangeZone[rangeI+1];
383  //vector maxDisp =
384  // (maxZoneI == -1 ? vector::zero : zoneDisp[maxZoneI]);
385  scalar maxDisp = (maxZoneI == -1 ? 0.0 : zoneDisp[maxZoneI][dir]);
386 
387  forAll(rPoints, i)
388  {
389  label pointi = rPoints[i];
390  scalar w = rWeights[i];
391  //curPoints[pointi] += (1.0-w)*minDisp+w*maxDisp;
392  curPoints[pointi][dir] += (1.0-w)*minDisp+w*maxDisp;
393  }
394  }
395  }
396  return tcurPoints;
397 }
398 
399 
400 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
GlobalIOList.H
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Tuple2.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
mapPolyMesh.H
Foam::IOobject::IOobject
IOobject(const IOobject &)=default
Copy construct.
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::IOobject::time
const Time & time() const
Return Time associated with the objectRegistry.
Definition: IOobject.C:493
Foam::interpolateXY
Field< Type > interpolateXY(const scalarField &xNew, const scalarField &xOld, const Field< Type > &yOld)
Definition: interpolateXY.C:40
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::primitiveMesh::nPoints
label nPoints() const noexcept
Number of mesh points.
Definition: primitiveMeshI.H:37
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::faceZoneMesh
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
Definition: faceZoneMeshFwd.H:44
SortableList.H
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< vector >
Foam::findLower
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::displacementInterpolationMotionSolver::curPoints
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Definition: displacementInterpolationMotionSolver.C:331
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:486
interpolateXY.H
Interpolates y values from one curve to another with a different x distribution.
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::motionSolver::coeffDict
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: motionSolver.H:150
Foam::displacementMotionSolver
Virtual base class for displacement motion solver.
Definition: displacementMotionSolver.H:53
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::points0MotionSolver::points0
pointField & points0()
Return reference to the reference field.
Definition: points0MotionSolver.H:117
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
points0
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, false)))
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::List< label >
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::motionSolver::mesh
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:144
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::GeometricField< vector, pointPatchField, pointMesh >
constant
constant condensation/saturation model.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
displacementInterpolationMotionSolver.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::IOobject::MUST_READ
Definition: IOobject.H:185