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-2013 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  (
87  "displacementInterpolationMotionSolver::"
88  "displacementInterpolationMotionSolver(const polyMesh&,"
89  "Istream&)"
90  ) << "Cannot find zone " << zoneName << endl
91  << "Valid zones are " << mesh.faceZones().names()
92  << exit(FatalError);
93  }
94 
95  const word& tableName = faceZoneToTable[i][1];
96 
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(), 0);
246 
247  forAll(meshCoords, pointI)
248  {
249  label rangeI = findLower(rangeToCoord, meshCoords[pointI]);
250 
251  if (rangeI == -1 || rangeI == rangeToCoord.size()-1)
252  {
254  (
255  "displacementInterpolationMotionSolver::"
256  "displacementInterpolationMotionSolver"
257  "(const polyMesh&, Istream&)"
258  ) << "Did not find point " << points0()[pointI]
259  << " coordinate " << meshCoords[pointI]
260  << " in ranges " << rangeToCoord
261  << abort(FatalError);
262  }
263  nRangePoints[rangeI]++;
264  }
265 
266  if (debug)
267  {
268  for (label rangeI = 0; rangeI < rangeToCoord.size()-1; rangeI++)
269  {
270  // Get the two zones bounding the range
271  Pout<< "direction " << dir << " : "
272  << "range from " << rangeToCoord[rangeI]
273  << " to " << rangeToCoord[rangeI+1]
274  << " contains " << nRangePoints[rangeI]
275  << " points." << endl;
276  }
277  }
278 
279  // Sort
280  rangePoints.setSize(nRangePoints.size());
281  rangeWeights.setSize(nRangePoints.size());
282  forAll(rangePoints, rangeI)
283  {
284  rangePoints[rangeI].setSize(nRangePoints[rangeI]);
285  rangeWeights[rangeI].setSize(nRangePoints[rangeI]);
286  }
287  nRangePoints = 0;
288  forAll(meshCoords, pointI)
289  {
290  label rangeI = findLower(rangeToCoord, meshCoords[pointI]);
291  label& nPoints = nRangePoints[rangeI];
292  rangePoints[rangeI][nPoints] = pointI;
293  rangeWeights[rangeI][nPoints] =
294  (meshCoords[pointI]-rangeToCoord[rangeI])
295  / (rangeToCoord[rangeI+1]-rangeToCoord[rangeI]);
296  nPoints++;
297  }
298  }
299 }
300 
301 
302 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
303 
306 {}
307 
308 
309 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
310 
313 {
314  if (mesh().nPoints() != points0().size())
315  {
317  (
318  "displacementInterpolationMotionSolver::curPoints() const"
319  ) << "The number of points in the mesh seems to have changed." << endl
320  << "In constant/polyMesh there are " << points0().size()
321  << " points; in the current mesh there are " << mesh().nPoints()
322  << " points." << exit(FatalError);
323  }
324 
325  tmp<pointField> tcurPoints(new pointField(points0()));
326  pointField& curPoints = tcurPoints();
327 
328  // Interpolate the displacement of the face zones.
329  vectorField zoneDisp(displacements_.size(), vector::zero);
330  forAll(zoneDisp, zoneI)
331  {
332  if (times_[zoneI].size())
333  {
334  zoneDisp[zoneI] = interpolateXY
335  (
336  mesh().time().value(),
337  times_[zoneI],
338  displacements_[zoneI]
339  );
340  }
341  }
342  if (debug)
343  {
344  Pout<< "Zone displacements:" << zoneDisp << endl;
345  }
346 
347 
348  // Interpolate the point location
349  for (direction dir = 0; dir < vector::nComponents; dir++)
350  {
351  const labelList& rangeZone = rangeToZone_[dir];
352  const labelListList& rangePoints = rangeToPoints_[dir];
353  const List<scalarField>& rangeWeights = rangeToWeights_[dir];
354 
355  for (label rangeI = 0; rangeI < rangeZone.size()-1; rangeI++)
356  {
357  const labelList& rPoints = rangePoints[rangeI];
358  const scalarField& rWeights = rangeWeights[rangeI];
359 
360  // Get the two zones bounding the range
361  label minZoneI = rangeZone[rangeI];
362  //vector minDisp =
363  // (minZoneI == -1 ? vector::zero : zoneDisp[minZoneI]);
364  scalar minDisp = (minZoneI == -1 ? 0.0 : zoneDisp[minZoneI][dir]);
365  label maxZoneI = rangeZone[rangeI+1];
366  //vector maxDisp =
367  // (maxZoneI == -1 ? vector::zero : zoneDisp[maxZoneI]);
368  scalar maxDisp = (maxZoneI == -1 ? 0.0 : zoneDisp[maxZoneI][dir]);
369 
370  forAll(rPoints, i)
371  {
372  label pointI = rPoints[i];
373  scalar w = rWeights[i];
374  //curPoints[pointI] += (1.0-w)*minDisp+w*maxDisp;
375  curPoints[pointI][dir] += (1.0-w)*minDisp+w*maxDisp;
376  }
377  }
378  }
379  return tcurPoints;
380 }
381 
382 
383 // ************************************************************************* //
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
unsigned char direction
Definition: direction.H:43
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:142
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
A class for handling words, derived from string.
Definition: word.H:59
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dynamicFvMesh & mesh
T & first()
Return the first element of the list.
Definition: UListI.H:117
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Namespace for OpenFOAM.
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
Number of components in this vector space.
Definition: VectorSpace.H:88
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
stressControl lookup("compactNormalStress") >> compactNormalStress
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)
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:65
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Macros for easy insertion into run-time selection tables.
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,.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Type gMin(const FieldField< Field, Type > &f)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Field< Type > interpolateXY(const scalarField &xNew, const scalarField &xOld, const Field< Type > &yOld)
Definition: interpolateXY.C:38
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
error FatalError
static const Vector zero
Definition: Vector.H:80
A List of objects of type <T> with automated input and output.
Definition: IOList.H:50
Interpolates y values from one curve to another with a different x distribution.
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:354
const Time & time() const
Return time.
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:269
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label nPoints
Type gMax(const FieldField< Field, Type > &f)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
A class for managing temporary objects.
Definition: PtrList.H:118
label nPoints() const
defineTypeNameAndDebug(combustionModel, 0)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
Virtual base class for displacement motion solver.