displacementInterpolationMotionSolver.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "SortableList.H"
29 #include "IOList.H"
30 #include "Tuple2.H"
31 #include "mapPolyMesh.H"
32 #include "interpolateXY.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(displacementInterpolationMotionSolver, 0);
39 
41  (
42  motionSolver,
43  displacementInterpolationMotionSolver,
44  dictionary
45  );
46 
47  template<>
48  const word IOList<Tuple2<scalar, vector>>::typeName("scalarVectorTable");
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
56 (
57  const polyMesh& mesh,
58  const dictionary& dict
59 )
60 :
61  points0MotionSolver(mesh, dict, typeName)
62 {
63  // Get zones and their interpolation tables for displacement
64  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
65 
66  List<Pair<word>> faceZoneToTable
67  (
68  coeffDict().lookup("interpolationTables")
69  );
70 
71  const faceZoneMesh& fZones = mesh.faceZones();
72 
73  times_.setSize(fZones.size());
74  displacements_.setSize(fZones.size());
75 
76  forAll(faceZoneToTable, i)
77  {
78  const word& zoneName = faceZoneToTable[i][0];
79  label zoneI = fZones.findZoneID(zoneName);
80 
81  if (zoneI == -1)
82  {
84  << "Cannot find zone " << zoneName << endl
85  << "Valid zones are " << mesh.faceZones().names()
86  << exit(FatalError);
87  }
88 
89  const word& tableName = faceZoneToTable[i][1];
90 
92  (
93  IOobject
94  (
95  tableName,
96  mesh.time().constant(),
97  "tables",
98  mesh,
101  false
102  )
103  );
104 
105  // Copy table
106  times_[zoneI].setSize(table.size());
107  displacements_[zoneI].setSize(table.size());
108 
109  forAll(table, j)
110  {
111  times_[zoneI][j] = table[j].first();
112  displacements_[zoneI][j] = table[j].second();
113  }
114  }
115 
116 
117 
118  // Sort points into bins according to position relative to faceZones
119  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
120  // Done in all three directions.
121 
122  for (direction dir = 0; dir < vector::nComponents; dir++)
123  {
124  // min and max coordinates of all faceZones
125  SortableList<scalar> zoneCoordinates(2*faceZoneToTable.size());
126 
127  forAll(faceZoneToTable, i)
128  {
129  const word& zoneName = faceZoneToTable[i][0];
130  const faceZone& fz = fZones[zoneName];
131 
132  scalar minCoord = vGreat;
133  scalar maxCoord = -vGreat;
134 
135  forAll(fz().meshPoints(), localI)
136  {
137  label pointi = fz().meshPoints()[localI];
138  const scalar coord = points0()[pointi][dir];
139  minCoord = min(minCoord, coord);
140  maxCoord = max(maxCoord, coord);
141  }
142 
143  zoneCoordinates[2*i] = returnReduce(minCoord, minOp<scalar>());
144  zoneCoordinates[2*i+1] = returnReduce(maxCoord, maxOp<scalar>());
145 
146  if (debug)
147  {
148  Pout<< "direction " << dir << " : "
149  << "zone " << zoneName
150  << " ranges from coordinate " << zoneCoordinates[2*i]
151  << " to " << zoneCoordinates[2*i+1]
152  << endl;
153  }
154  }
155  zoneCoordinates.sort();
156 
157  // Slightly tweak min and max face zone so points sort within
158  zoneCoordinates[0] -= small;
159  zoneCoordinates.last() += small;
160 
161  // Check if we have static min and max mesh bounds
162  const scalarField meshCoords(points0().component(dir));
163 
164  scalar minCoord = gMin(meshCoords);
165  scalar maxCoord = gMax(meshCoords);
166 
167  if (debug)
168  {
169  Pout<< "direction " << dir << " : "
170  << "mesh ranges from coordinate " << minCoord << " to "
171  << maxCoord << endl;
172  }
173 
174  // Make copy of zoneCoordinates; include min and max of mesh
175  // if necessary. Mark min and max with zoneI=-1.
176  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
177 
178  labelList& rangeZone = rangeToZone_[dir];
179  labelListList& rangePoints = rangeToPoints_[dir];
180  List<scalarField>& rangeWeights = rangeToWeights_[dir];
181 
182  scalarField rangeToCoord(zoneCoordinates.size());
183  rangeZone.setSize(zoneCoordinates.size());
184  label rangeI = 0;
185 
186  if (minCoord < zoneCoordinates[0])
187  {
188  label sz = rangeZone.size();
189  rangeToCoord.setSize(sz+1);
190  rangeZone.setSize(sz+1);
191  rangeToCoord[rangeI] = minCoord-small;
192  rangeZone[rangeI] = -1;
193 
194  if (debug)
195  {
196  Pout<< "direction " << dir << " : "
197  << "range " << rangeI << " at coordinate "
198  << rangeToCoord[rangeI] << " from min of mesh "
199  << rangeZone[rangeI] << endl;
200  }
201  rangeI = 1;
202  }
203  forAll(zoneCoordinates, i)
204  {
205  rangeToCoord[rangeI] = zoneCoordinates[i];
206  rangeZone[rangeI] = zoneCoordinates.indices()[i]/2;
207 
208  if (debug)
209  {
210  Pout<< "direction " << dir << " : "
211  << "range " << rangeI << " at coordinate "
212  << rangeToCoord[rangeI]
213  << " from zone " << rangeZone[rangeI] << endl;
214  }
215  rangeI++;
216  }
217  if (maxCoord > zoneCoordinates.last())
218  {
219  label sz = rangeToCoord.size();
220  rangeToCoord.setSize(sz+1);
221  rangeZone.setSize(sz+1);
222  rangeToCoord[sz] = maxCoord+small;
223  rangeZone[sz] = -1;
224 
225  if (debug)
226  {
227  Pout<< "direction " << dir << " : "
228  << "range " << rangeI << " at coordinate "
229  << rangeToCoord[sz] << " from max of mesh "
230  << rangeZone[sz] << endl;
231  }
232  }
233 
234 
235  // Sort the points
236  // ~~~~~~~~~~~~~~~
237 
238  // Count all the points in between rangeI and rangeI+1
239  labelList nRangePoints(rangeToCoord.size(), 0);
240 
241  forAll(meshCoords, pointi)
242  {
243  label rangeI = findLower(rangeToCoord, meshCoords[pointi]);
244 
245  if (rangeI == -1 || rangeI == rangeToCoord.size()-1)
246  {
248  << "Did not find point " << points0()[pointi]
249  << " coordinate " << meshCoords[pointi]
250  << " in ranges " << rangeToCoord
251  << abort(FatalError);
252  }
253  nRangePoints[rangeI]++;
254  }
255 
256  if (debug)
257  {
258  for (label rangeI = 0; rangeI < rangeToCoord.size()-1; rangeI++)
259  {
260  // Get the two zones bounding the range
261  Pout<< "direction " << dir << " : "
262  << "range from " << rangeToCoord[rangeI]
263  << " to " << rangeToCoord[rangeI+1]
264  << " contains " << nRangePoints[rangeI]
265  << " points." << endl;
266  }
267  }
268 
269  // Sort
270  rangePoints.setSize(nRangePoints.size());
271  rangeWeights.setSize(nRangePoints.size());
272  forAll(rangePoints, rangeI)
273  {
274  rangePoints[rangeI].setSize(nRangePoints[rangeI]);
275  rangeWeights[rangeI].setSize(nRangePoints[rangeI]);
276  }
277  nRangePoints = 0;
278  forAll(meshCoords, pointi)
279  {
280  label rangeI = findLower(rangeToCoord, meshCoords[pointi]);
281  label& nPoints = nRangePoints[rangeI];
282  rangePoints[rangeI][nPoints] = pointi;
283  rangeWeights[rangeI][nPoints] =
284  (meshCoords[pointi]-rangeToCoord[rangeI])
285  / (rangeToCoord[rangeI+1]-rangeToCoord[rangeI]);
286  nPoints++;
287  }
288  }
289 }
290 
291 
292 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
293 
296 {}
297 
298 
299 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
300 
303 {
304  if (mesh().nPoints() != points0().size())
305  {
307  << "The number of points in the mesh seems to have changed." << endl
308  << "In constant/polyMesh there are " << points0().size()
309  << " points; in the current mesh there are " << mesh().nPoints()
310  << " points." << exit(FatalError);
311  }
312 
313  tmp<pointField> tcurPoints(new pointField(points0()));
314  pointField& curPoints = tcurPoints.ref();
315 
316  // Interpolate the displacement of the face zones.
317  vectorField zoneDisp(displacements_.size(), Zero);
318  forAll(zoneDisp, zoneI)
319  {
320  if (times_[zoneI].size())
321  {
322  zoneDisp[zoneI] = interpolateXY
323  (
324  mesh().time().value(),
325  times_[zoneI],
326  displacements_[zoneI]
327  );
328  }
329  }
330  if (debug)
331  {
332  Pout<< "Zone displacements:" << zoneDisp << endl;
333  }
334 
335 
336  // Interpolate the point location
337  for (direction dir = 0; dir < vector::nComponents; dir++)
338  {
339  const labelList& rangeZone = rangeToZone_[dir];
340  const labelListList& rangePoints = rangeToPoints_[dir];
341  const List<scalarField>& rangeWeights = rangeToWeights_[dir];
342 
343  for (label rangeI = 0; rangeI < rangeZone.size()-1; rangeI++)
344  {
345  const labelList& rPoints = rangePoints[rangeI];
346  const scalarField& rWeights = rangeWeights[rangeI];
347 
348  // Get the two zones bounding the range
349  label minZoneI = rangeZone[rangeI];
350  // vector minDisp =
351  // (minZoneI == -1 ? vector::zero : zoneDisp[minZoneI]);
352  scalar minDisp = (minZoneI == -1 ? 0.0 : zoneDisp[minZoneI][dir]);
353  label maxZoneI = rangeZone[rangeI+1];
354  // vector maxDisp =
355  // (maxZoneI == -1 ? vector::zero : zoneDisp[maxZoneI]);
356  scalar maxDisp = (maxZoneI == -1 ? 0.0 : zoneDisp[maxZoneI][dir]);
357 
358  forAll(rPoints, i)
359  {
360  label pointi = rPoints[i];
361  scalar w = rWeights[i];
362  // curPoints[pointi] += (1.0-w)*minDisp+w*maxDisp;
363  curPoints[pointi][dir] += (1.0-w)*minDisp+w*maxDisp;
364  }
365  }
366  }
367  return tcurPoints;
368 }
369 
370 
371 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:476
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Type gMin(const FieldField< Field, Type > &f)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:80
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
uint8_t direction
Definition: direction.H:45
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A List of objects of type <T> with automated input and output.
Definition: IOList.H:50
T & first()
Return the first element of the list.
Definition: UListI.H:114
Macros for easy insertion into run-time selection tables.
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Interpolates y values from one curve to another with a different x distribution.
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
static const zero Zero
Definition: zero.H:97
Field< Type > interpolateXY(const scalarField &xNew, const scalarField &xOld, const Field< Type > &yOld)
Definition: interpolateXY.C:38
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Type gMax(const FieldField< Field, Type > &f)
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Virtual base class for displacement motion solvers.
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:256
label findLower(const ListType &, typename ListType::const_reference, const label stary, const BinaryOp &bop)
Find last element < given value in sorted list and return index,.
label nPoints() const
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
displacementInterpolationMotionSolver(const polyMesh &, const dictionary &dict)
Construct from polyMesh and dictionary.
Namespace for OpenFOAM.