67void Foam::displacementInterpolationMotionSolver::calcInterpolation()
72 List<Pair<word>> faceZoneToTable
80 displacements_.
setSize(fZones.size());
84 const word& zoneName = faceZoneToTable[i][0];
85 label zoneI = fZones.findZoneID(zoneName);
90 <<
"Cannot find zone " << zoneName <<
endl
91 <<
"Valid zones are " << fZones.names()
95 const word& tableName = faceZoneToTable[i][1];
97 GlobalIOList<Tuple2<scalar, vector>> table
112 times_[zoneI].setSize(table.size());
113 displacements_[zoneI].setSize(table.size());
117 times_[zoneI][j] = table[j].first();
118 displacements_[zoneI][j] = table[j].second();
131 SortableList<scalar> zoneCoordinates(2*faceZoneToTable.size());
133 forAll(faceZoneToTable, i)
135 const word& zoneName = faceZoneToTable[i][0];
136 const faceZone& fz = fZones[zoneName];
140 for (
const label pointi : fz().meshPoints())
142 const scalar coord =
points0()[pointi][dir];
148 zoneCoordinates[2*i] =
limits.min();
149 zoneCoordinates[2*i+1] =
limits.max();
153 Pout<<
"direction " << dir <<
" : "
154 <<
"zone " << zoneName
155 <<
" ranges from coordinate "
160 zoneCoordinates.sort();
163 zoneCoordinates[0] -= SMALL;
164 zoneCoordinates.last() += SMALL;
169 scalar minCoord, maxCoord;
179 Pout<<
"direction " << dir <<
" : "
180 <<
"mesh ranges from coordinate "
181 << minCoord <<
" to " << maxCoord
189 labelList& rangeZone = rangeToZone_[dir];
191 List<scalarField>& rangeWeights = rangeToWeights_[dir];
194 rangeZone.setSize(zoneCoordinates.size());
197 if (minCoord < zoneCoordinates[0])
199 label sz = rangeZone.size();
200 rangeToCoord.setSize(sz+1);
201 rangeZone.setSize(sz+1);
202 rangeToCoord[rangeI] = minCoord-SMALL;
203 rangeZone[rangeI] = -1;
207 Pout<<
"direction " << dir <<
" : "
208 <<
"range " << rangeI <<
" at coordinate "
209 << rangeToCoord[rangeI] <<
" from min of mesh "
210 << rangeZone[rangeI] <<
endl;
214 forAll(zoneCoordinates, i)
216 rangeToCoord[rangeI] = zoneCoordinates[i];
217 rangeZone[rangeI] = zoneCoordinates.indices()[i]/2;
221 Pout<<
"direction " << dir <<
" : "
222 <<
"range " << rangeI <<
" at coordinate "
223 << rangeToCoord[rangeI]
224 <<
" from zone " << rangeZone[rangeI] <<
endl;
228 if (maxCoord > zoneCoordinates.last())
230 label sz = rangeToCoord.size();
231 rangeToCoord.setSize(sz+1);
232 rangeZone.setSize(sz+1);
233 rangeToCoord[sz] = maxCoord+SMALL;
238 Pout<<
"direction " << dir <<
" : "
239 <<
"range " << rangeI <<
" at coordinate "
240 << rangeToCoord[sz] <<
" from max of mesh "
241 << rangeZone[sz] <<
endl;
252 forAll(meshCoords, pointi)
254 label rangeI =
findLower(rangeToCoord, meshCoords[pointi]);
256 if (rangeI == -1 || rangeI == rangeToCoord.size()-1)
259 <<
"Did not find point " <<
points0()[pointi]
260 <<
" coordinate " << meshCoords[pointi]
261 <<
" in ranges " << rangeToCoord
264 nRangePoints[rangeI]++;
269 for (label rangeI = 0; rangeI < rangeToCoord.size()-1; rangeI++)
272 Pout<<
"direction " << dir <<
" : "
273 <<
"range from " << rangeToCoord[rangeI]
274 <<
" to " << rangeToCoord[rangeI+1]
275 <<
" contains " << nRangePoints[rangeI]
276 <<
" points." <<
endl;
281 rangePoints.setSize(nRangePoints.size());
282 rangeWeights.setSize(nRangePoints.size());
283 forAll(rangePoints, rangeI)
285 rangePoints[rangeI].setSize(nRangePoints[rangeI]);
286 rangeWeights[rangeI].setSize(nRangePoints[rangeI]);
289 forAll(meshCoords, pointi)
291 label rangeI =
findLower(rangeToCoord, meshCoords[pointi]);
292 label&
nPoints = nRangePoints[rangeI];
293 rangePoints[rangeI][
nPoints] = pointi;
294 rangeWeights[rangeI][
nPoints] =
295 (meshCoords[pointi]-rangeToCoord[rangeI])
296 / (rangeToCoord[rangeI+1]-rangeToCoord[rangeI]);
305Foam::displacementInterpolationMotionSolver::
306displacementInterpolationMotionSolver
318Foam::displacementInterpolationMotionSolver::
319displacementInterpolationMotionSolver
341 <<
"The number of points in the mesh seems to have changed." <<
endl
342 <<
"In constant/polyMesh there are " <<
points0().
size()
343 <<
" points; in the current mesh there are " <<
mesh().
nPoints()
348 auto& curPoints = tcurPoints.ref();
354 if (times_[zoneI].size())
358 mesh().time().value(),
360 displacements_[zoneI]
366 Pout<<
"Zone displacements:" << zoneDisp <<
endl;
373 const labelList& rangeZone = rangeToZone_[dir];
377 for (label rangeI = 0; rangeI < rangeZone.size()-1; rangeI++)
379 const labelList& rPoints = rangePoints[rangeI];
380 const scalarField& rWeights = rangeWeights[rangeI];
383 label minZoneI = rangeZone[rangeI];
386 scalar minDisp = (minZoneI == -1 ? 0.0 : zoneDisp[minZoneI][dir]);
387 label maxZoneI = rangeZone[rangeI+1];
390 scalar maxDisp = (maxZoneI == -1 ? 0.0 : zoneDisp[maxZoneI][dir]);
394 label pointi = rPoints[i];
395 scalar w = rWeights[i];
397 curPoints[pointi][dir] += (1.0-w)*minDisp+w*maxDisp;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
label size() const noexcept
The number of elements in list.
IOList with global data (so optionally read from master).
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
const Time & time() const noexcept
Return Time associated with the objectRegistry.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void setSize(label n)
Alias for resize().
void size(const label n)
Older name for setAddressableSize.
static constexpr direction nComponents
Number of components in this vector space.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream.
Mesh motion solver for an fvMesh.
virtual tmp< pointField > curPoints() const
Return point location obtained from the current motion field.
Virtual base class for displacement motion solver.
pointVectorField & pointDisplacement() noexcept
Return reference to the point motion displacement field.
displacementMotionSolver(const displacementMotionSolver &)=delete
No copy construct.
Virtual base class for mesh motion solver.
const polyMesh & mesh() const
Return reference to mesh.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
pointField & points0() noexcept
Return reference to the reference ('0') pointField.
Mesh consisting of general polyhedral cells.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
label nPoints() const noexcept
Number of mesh points.
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Interpolates y values from one curve to another with a different x distribution.
Namespace for handling debugging switches.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
MinMax< scalar > scalarMinMax
A scalar min/max range.
vectorIOField pointIOField
pointIOField is a vectorIOField.
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
Field< Type > interpolateXY(const scalarField &xNew, const scalarField &xOld, const Field< Type > &yOld)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define forAll(list, i)
Loop across all elements in list.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))