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 | Copyright (C) 2011-2016 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 Foam::displacementInterpolationMotionSolver::
57 displacementInterpolationMotionSolver
58 (
59  const polyMesh& mesh,
60  const IOdictionary& dict
61 )
62 :
63  displacementMotionSolver(mesh, dict, typeName)
64 {
65  // Get zones and their interpolation tables for displacement
66  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
67 
68  List<Pair<word>> faceZoneToTable
69  (
70  coeffDict().lookup("interpolationTables")
71  );
72 
73  const faceZoneMesh& fZones = mesh.faceZones();
74 
75  times_.setSize(fZones.size());
76  displacements_.setSize(fZones.size());
77 
78  forAll(faceZoneToTable, i)
79  {
80  const word& zoneName = faceZoneToTable[i][0];
81  label zoneI = fZones.findZoneID(zoneName);
82 
83  if (zoneI == -1)
84  {
86  << "Cannot find zone " << zoneName << endl
87  << "Valid zones are " << mesh.faceZones().names()
88  << exit(FatalError);
89  }
90 
91  const word& tableName = faceZoneToTable[i][1];
92 
94  (
95  IOobject
96  (
97  tableName,
98  mesh.time().constant(),
99  "tables",
100  mesh,
103  false
104  )
105  );
106 
107  // Copy table
108  times_[zoneI].setSize(table.size());
109  displacements_[zoneI].setSize(table.size());
110 
111  forAll(table, j)
112  {
113  times_[zoneI][j] = table[j].first();
114  displacements_[zoneI][j] = table[j].second();
115  }
116  }
117 
118 
119 
120  // Sort points into bins according to position relative to faceZones
121  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
122  // Done in all three directions.
123 
124  for (direction dir = 0; dir < vector::nComponents; dir++)
125  {
126  // min and max coordinates of all faceZones
127  SortableList<scalar> zoneCoordinates(2*faceZoneToTable.size());
128 
129  forAll(faceZoneToTable, i)
130  {
131  const word& zoneName = faceZoneToTable[i][0];
132  const faceZone& fz = fZones[zoneName];
133 
134  scalar minCoord = VGREAT;
135  scalar maxCoord = -VGREAT;
136 
137  forAll(fz().meshPoints(), localI)
138  {
139  label pointi = fz().meshPoints()[localI];
140  const scalar coord = points0()[pointi][dir];
141  minCoord = min(minCoord, coord);
142  maxCoord = max(maxCoord, coord);
143  }
144 
145  zoneCoordinates[2*i] = returnReduce(minCoord, minOp<scalar>());
146  zoneCoordinates[2*i+1] = returnReduce(maxCoord, maxOp<scalar>());
147 
148  if (debug)
149  {
150  Pout<< "direction " << dir << " : "
151  << "zone " << zoneName
152  << " ranges from coordinate " << zoneCoordinates[2*i]
153  << " to " << zoneCoordinates[2*i+1]
154  << endl;
155  }
156  }
157  zoneCoordinates.sort();
158 
159  // Slightly tweak min and max face zone so points sort within
160  zoneCoordinates[0] -= SMALL;
161  zoneCoordinates.last() += SMALL;
162 
163  // Check if we have static min and max mesh bounds
164  const scalarField meshCoords(points0().component(dir));
165 
166  scalar minCoord = gMin(meshCoords);
167  scalar maxCoord = gMax(meshCoords);
168 
169  if (debug)
170  {
171  Pout<< "direction " << dir << " : "
172  << "mesh ranges from coordinate " << minCoord << " to "
173  << maxCoord << endl;
174  }
175 
176  // Make copy of zoneCoordinates; include min and max of mesh
177  // if necessary. Mark min and max with zoneI=-1.
178  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
179 
180  labelList& rangeZone = rangeToZone_[dir];
181  labelListList& rangePoints = rangeToPoints_[dir];
182  List<scalarField>& rangeWeights = rangeToWeights_[dir];
183 
184  scalarField rangeToCoord(zoneCoordinates.size());
185  rangeZone.setSize(zoneCoordinates.size());
186  label rangeI = 0;
187 
188  if (minCoord < zoneCoordinates[0])
189  {
190  label sz = rangeZone.size();
191  rangeToCoord.setSize(sz+1);
192  rangeZone.setSize(sz+1);
193  rangeToCoord[rangeI] = minCoord-SMALL;
194  rangeZone[rangeI] = -1;
195 
196  if (debug)
197  {
198  Pout<< "direction " << dir << " : "
199  << "range " << rangeI << " at coordinate "
200  << rangeToCoord[rangeI] << " from min of mesh "
201  << rangeZone[rangeI] << endl;
202  }
203  rangeI = 1;
204  }
205  forAll(zoneCoordinates, i)
206  {
207  rangeToCoord[rangeI] = zoneCoordinates[i];
208  rangeZone[rangeI] = zoneCoordinates.indices()[i]/2;
209 
210  if (debug)
211  {
212  Pout<< "direction " << dir << " : "
213  << "range " << rangeI << " at coordinate "
214  << rangeToCoord[rangeI]
215  << " from zone " << rangeZone[rangeI] << endl;
216  }
217  rangeI++;
218  }
219  if (maxCoord > zoneCoordinates.last())
220  {
221  label sz = rangeToCoord.size();
222  rangeToCoord.setSize(sz+1);
223  rangeZone.setSize(sz+1);
224  rangeToCoord[sz] = maxCoord+SMALL;
225  rangeZone[sz] = -1;
226 
227  if (debug)
228  {
229  Pout<< "direction " << dir << " : "
230  << "range " << rangeI << " at coordinate "
231  << rangeToCoord[sz] << " from max of mesh "
232  << rangeZone[sz] << endl;
233  }
234  }
235 
236 
237  // Sort the points
238  // ~~~~~~~~~~~~~~~
239 
240  // Count all the points inbetween rangeI and rangeI+1
241  labelList nRangePoints(rangeToCoord.size(), 0);
242 
243  forAll(meshCoords, pointi)
244  {
245  label rangeI = findLower(rangeToCoord, meshCoords[pointi]);
246 
247  if (rangeI == -1 || rangeI == rangeToCoord.size()-1)
248  {
250  << "Did not find point " << points0()[pointi]
251  << " coordinate " << meshCoords[pointi]
252  << " in ranges " << rangeToCoord
253  << abort(FatalError);
254  }
255  nRangePoints[rangeI]++;
256  }
257 
258  if (debug)
259  {
260  for (label rangeI = 0; rangeI < rangeToCoord.size()-1; rangeI++)
261  {
262  // Get the two zones bounding the range
263  Pout<< "direction " << dir << " : "
264  << "range from " << rangeToCoord[rangeI]
265  << " to " << rangeToCoord[rangeI+1]
266  << " contains " << nRangePoints[rangeI]
267  << " points." << endl;
268  }
269  }
270 
271  // Sort
272  rangePoints.setSize(nRangePoints.size());
273  rangeWeights.setSize(nRangePoints.size());
274  forAll(rangePoints, rangeI)
275  {
276  rangePoints[rangeI].setSize(nRangePoints[rangeI]);
277  rangeWeights[rangeI].setSize(nRangePoints[rangeI]);
278  }
279  nRangePoints = 0;
280  forAll(meshCoords, pointi)
281  {
282  label rangeI = findLower(rangeToCoord, meshCoords[pointi]);
283  label& nPoints = nRangePoints[rangeI];
284  rangePoints[rangeI][nPoints] = pointi;
285  rangeWeights[rangeI][nPoints] =
286  (meshCoords[pointi]-rangeToCoord[rangeI])
287  / (rangeToCoord[rangeI+1]-rangeToCoord[rangeI]);
288  nPoints++;
289  }
290  }
291 }
292 
293 
294 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
295 
298 {}
299 
300 
301 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
302 
305 {
306  if (mesh().nPoints() != points0().size())
307  {
309  << "The number of points in the mesh seems to have changed." << endl
310  << "In constant/polyMesh there are " << points0().size()
311  << " points; in the current mesh there are " << mesh().nPoints()
312  << " points." << exit(FatalError);
313  }
314 
315  tmp<pointField> tcurPoints(new pointField(points0()));
316  pointField& curPoints = tcurPoints.ref();
317 
318  // Interpolate the displacement of the face zones.
319  vectorField zoneDisp(displacements_.size(), Zero);
320  forAll(zoneDisp, zoneI)
321  {
322  if (times_[zoneI].size())
323  {
324  zoneDisp[zoneI] = interpolateXY
325  (
326  mesh().time().value(),
327  times_[zoneI],
328  displacements_[zoneI]
329  );
330  }
331  }
332  if (debug)
333  {
334  Pout<< "Zone displacements:" << zoneDisp << endl;
335  }
336 
337 
338  // Interpolate the point location
339  for (direction dir = 0; dir < vector::nComponents; dir++)
340  {
341  const labelList& rangeZone = rangeToZone_[dir];
342  const labelListList& rangePoints = rangeToPoints_[dir];
343  const List<scalarField>& rangeWeights = rangeToWeights_[dir];
344 
345  for (label rangeI = 0; rangeI < rangeZone.size()-1; rangeI++)
346  {
347  const labelList& rPoints = rangePoints[rangeI];
348  const scalarField& rWeights = rangeWeights[rangeI];
349 
350  // Get the two zones bounding the range
351  label minZoneI = rangeZone[rangeI];
352  //vector minDisp =
353  // (minZoneI == -1 ? vector::zero : zoneDisp[minZoneI]);
354  scalar minDisp = (minZoneI == -1 ? 0.0 : zoneDisp[minZoneI][dir]);
355  label maxZoneI = rangeZone[rangeI+1];
356  //vector maxDisp =
357  // (maxZoneI == -1 ? vector::zero : zoneDisp[maxZoneI]);
358  scalar maxDisp = (maxZoneI == -1 ? 0.0 : zoneDisp[maxZoneI][dir]);
359 
360  forAll(rPoints, i)
361  {
362  label pointi = rPoints[i];
363  scalar w = rWeights[i];
364  //curPoints[pointi] += (1.0-w)*minDisp+w*maxDisp;
365  curPoints[pointi][dir] += (1.0-w)*minDisp+w*maxDisp;
366  }
367  }
368  }
369  return tcurPoints;
370 }
371 
372 
373 // ************************************************************************* //
const Time & time() const
Return time.
Virtual base class for displacement motion solver.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
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:68
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:341
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:96
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:256
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:91
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
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void setSize(const label)
Reset size of List.
Definition: List.C:295
label nPoints() const
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,.
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
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:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29