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-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 
48  (
49  displacementMotionSolver,
50  displacementInterpolationMotionSolver,
51  displacement
52  );
53 
54  template<>
55  const word IOList<Tuple2<scalar, vector> >::typeName("scalarVectorTable");
56 }
57 
58 
59 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
60 
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 " << fZones.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 inbetween 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
293 
296 (
297  const polyMesh& mesh,
298  const IOdictionary& dict
299 )
300 :
302 {
303  calcInterpolation();
304 }
305 
306 
309 (
310  const polyMesh& mesh,
311  const IOdictionary& dict,
312  const pointVectorField& pointDisplacement,
313  const pointIOField& points0
314 )
315 :
316  displacementMotionSolver(mesh, dict, pointDisplacement, points0, typeName)
317 {
318  calcInterpolation();
319 }
320 
321 
322 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
323 
326 {}
327 
328 
329 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
330 
333 {
334  if (mesh().nPoints() != points0().size())
335  {
337  << "The number of points in the mesh seems to have changed." << endl
338  << "In constant/polyMesh there are " << points0().size()
339  << " points; in the current mesh there are " << mesh().nPoints()
340  << " points." << exit(FatalError);
341  }
342 
343  tmp<pointField> tcurPoints(new pointField(points0()));
344  pointField& curPoints = tcurPoints();
345 
346  // Interpolate the displacement of the face zones.
347  vectorField zoneDisp(displacements_.size(), vector::zero);
348  forAll(zoneDisp, zoneI)
349  {
350  if (times_[zoneI].size())
351  {
352  zoneDisp[zoneI] = interpolateXY
353  (
354  mesh().time().value(),
355  times_[zoneI],
356  displacements_[zoneI]
357  );
358  }
359  }
360  if (debug)
361  {
362  Pout<< "Zone displacements:" << zoneDisp << endl;
363  }
364 
365 
366  // Interpolate the point location
367  for (direction dir = 0; dir < vector::nComponents; dir++)
368  {
369  const labelList& rangeZone = rangeToZone_[dir];
370  const labelListList& rangePoints = rangeToPoints_[dir];
371  const List<scalarField>& rangeWeights = rangeToWeights_[dir];
372 
373  for (label rangeI = 0; rangeI < rangeZone.size()-1; rangeI++)
374  {
375  const labelList& rPoints = rangePoints[rangeI];
376  const scalarField& rWeights = rangeWeights[rangeI];
377 
378  // Get the two zones bounding the range
379  label minZoneI = rangeZone[rangeI];
380  //vector minDisp =
381  // (minZoneI == -1 ? vector::zero : zoneDisp[minZoneI]);
382  scalar minDisp = (minZoneI == -1 ? 0.0 : zoneDisp[minZoneI][dir]);
383  label maxZoneI = rangeZone[rangeI+1];
384  //vector maxDisp =
385  // (maxZoneI == -1 ? vector::zero : zoneDisp[maxZoneI]);
386  scalar maxDisp = (maxZoneI == -1 ? 0.0 : zoneDisp[maxZoneI][dir]);
387 
388  forAll(rPoints, i)
389  {
390  label pointI = rPoints[i];
391  scalar w = rWeights[i];
392  //curPoints[pointI] += (1.0-w)*minDisp+w*maxDisp;
393  curPoints[pointI][dir] += (1.0-w)*minDisp+w*maxDisp;
394  }
395  }
396  }
397  return tcurPoints;
398 }
399 
400 
401 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::maxOp
Definition: ops.H:172
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:41
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOField< vector >
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::SortableList::sort
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:95
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Tuple2.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
mapPolyMesh.H
Foam::minOp
Definition: ops.H:173
Foam::displacementInterpolationMotionSolver::displacements_
List< vectorField > displacements_
Interpolation table. From faceZone to displacements.
Definition: displacementInterpolationMotionSolver.H:70
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::displacementInterpolationMotionSolver::times_
List< scalarField > times_
Interpolation table. From faceZone to times.
Definition: displacementInterpolationMotionSolver.H:67
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
@ nComponents
Number of components in this vector space.
Definition: VectorSpace.H:88
Foam::displacementInterpolationMotionSolver::rangeToPoints_
FixedList< labelListList, 3 > rangeToPoints_
Per direction, per range the points that are in it.
Definition: displacementInterpolationMotionSolver.H:79
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::IOobject::time
const Time & time() const
Return time.
Definition: IOobject.C:245
Foam::interpolateXY
Field< Type > interpolateXY(const scalarField &xNew, const scalarField &xOld, const Field< Type > &yOld)
Definition: interpolateXY.C:38
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
constant
Constant dispersed-phase particle diameter model.
Foam::displacementInterpolationMotionSolver::rangeToWeights_
FixedList< List< scalarField >, 3 > rangeToWeights_
Per direction, per range the weight of the points relative.
Definition: displacementInterpolationMotionSolver.H:83
IOList.H
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
SortableList.H
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::displacementInterpolationMotionSolver::curPoints
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Definition: displacementInterpolationMotionSolver.C:332
Foam::ZoneMesh< faceZone, polyMesh >
interpolateXY.H
Interpolates y values from one curve to another with a different x distribution.
Foam::findLower
label findLower(const ListType &, typename ListType::const_reference, const label start, const BinaryOp &bop)
Find last element < given value in sorted list and return index,.
Foam::displacementInterpolationMotionSolver::calcInterpolation
void calcInterpolation()
Read settings.
Definition: displacementInterpolationMotionSolver.C:61
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:131
Foam::displacementMotionSolver
Virtual base class for displacement motion solver.
Definition: displacementMotionSolver.H:55
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:348
Foam::ZoneMesh::names
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:263
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::displacementMotionSolver::points0
pointField & points0()
Return reference to the reference field.
Definition: displacementMotionSolver.H:151
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::displacementInterpolationMotionSolver::rangeToZone_
FixedList< labelList, 3 > rangeToZone_
Per direction, per range the index of the lower.
Definition: displacementInterpolationMotionSolver.H:76
Foam::displacementInterpolationMotionSolver::~displacementInterpolationMotionSolver
~displacementInterpolationMotionSolver()
Destructor.
Definition: displacementInterpolationMotionSolver.C:325
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::displacementInterpolationMotionSolver::displacementInterpolationMotionSolver
displacementInterpolationMotionSolver(const displacementInterpolationMotionSolver &)
Disallow default bitwise copy construct.
Foam::List
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
Foam::IOList
A List of objects of type <T> with automated input and output.
Definition: IOList.H:50
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::motionSolver::mesh
const polyMesh & mesh() const
Return reference to mesh.
Definition: motionSolver.H:125
Foam::SortableList::indices
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:90
Foam::List< T >::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
displacementInterpolationMotionSolver.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562