interpolatingSolidBodyMotionSolver: Added mapMesh function

to support run-time mesh-to-mesh mapping.  The points0 are reset to the points
of the new mesh, i.e. the displacement is assumed 0 after mapping and the motion
functions need to take this into account.
This commit is contained in:
Henry Weller
2022-08-17 11:53:07 +01:00
parent c0950a2aa9
commit d6c62a9e7a
5 changed files with 75 additions and 40 deletions

View File

@ -201,7 +201,7 @@ void Foam::componentDisplacementMotionSolver::topoChange
void Foam::componentDisplacementMotionSolver::mapMesh(const polyMeshMap& map)
{
points0_ == mesh().points().component(cmpt_);
points0_ = mesh().points().component(cmpt_);
pointDisplacement_ == Zero;
}

View File

@ -160,7 +160,7 @@ void Foam::points0MotionSolver::topoChange(const polyTopoChangeMap& map)
void Foam::points0MotionSolver::mapMesh(const polyMeshMap& map)
{
points0_ == mesh().points();
points0_.primitiveFieldRef() = mesh().points();
}

View File

@ -110,6 +110,7 @@ public:
virtual void topoChange(const polyTopoChangeMap&);
//- Update from another mesh using the given map
// Resets points0 to the points of the new mesh
virtual void mapMesh(const polyMeshMap&);
//- Update corresponding to the given distribution map

View File

@ -43,6 +43,46 @@ namespace Foam
}
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
void Foam::interpolatingSolidBodyMotionSolver::calcScale()
{
const pointMesh& pMesh = pointMesh::New(mesh());
pointPatchDist pDist(pMesh, patchSet_, points0());
// Scaling: 1 up to di then linear down to 0 at do away from patches
scale_.primitiveFieldRef() =
min
(
max
(
(do_ - pDist.primitiveField())/(do_ - di_),
scalar(0)
),
scalar(1)
);
// Convert the scale function to a cosine
scale_.primitiveFieldRef() =
min
(
max
(
0.5
- 0.5
*cos(scale_.primitiveField()
*Foam::constant::mathematical::pi),
scalar(0)
),
scalar(1)
);
pointConstraints::New(pMesh).constrain(scale_);
scale_.write();
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::interpolatingSolidBodyMotionSolver::interpolatingSolidBodyMotionSolver
@ -67,50 +107,14 @@ Foam::interpolatingSolidBodyMotionSolver::interpolatingSolidBodyMotionSolver
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
IOobject::NO_WRITE
),
pointMesh::New(mesh),
dimensionedScalar(dimless, 0)
)
{
// Calculate scaling factor everywhere
{
const pointMesh& pMesh = pointMesh::New(mesh);
pointPatchDist pDist(pMesh, patchSet_, points0());
// Scaling: 1 up to di then linear down to 0 at do away from patches
scale_.primitiveFieldRef() =
min
(
max
(
(do_ - pDist.primitiveField())/(do_ - di_),
scalar(0)
),
scalar(1)
);
// Convert the scale function to a cosine
scale_.primitiveFieldRef() =
min
(
max
(
0.5
- 0.5
*cos(scale_.primitiveField()
*Foam::constant::mathematical::pi),
scalar(0)
),
scalar(1)
);
pointConstraints::New(pMesh).constrain(scale_);
scale_.write();
}
calcScale();
}
@ -131,6 +135,7 @@ Foam::interpolatingSolidBodyMotionSolver::curPoints() const
tmp<pointField> tpoints(new pointField(points0));
pointField& points = tpoints.ref();
Info << points.size() << " " << scale_.size() << endl;
forAll(points, pointi)
{
// Move non-stationary points
@ -157,4 +162,25 @@ Foam::interpolatingSolidBodyMotionSolver::curPoints() const
}
void Foam::interpolatingSolidBodyMotionSolver::topoChange
(
const polyTopoChangeMap& map
)
{
// Pending implementation of the inverse transformation of points0
NotImplemented;
}
void Foam::interpolatingSolidBodyMotionSolver::mapMesh(const polyMeshMap& map)
{
InfoInFunction << endl;
points0MotionSolver::mapMesh(map);
// scale is resized by the meshToMesh mapper
scale_ = Zero;
calcScale();
}
// ************************************************************************* //

View File

@ -75,6 +75,8 @@ class interpolatingSolidBodyMotionSolver
//- Current interpolation scale (1 at patches, 0 at distance_)
pointScalarField scale_;
void calcScale();
public:
@ -112,6 +114,12 @@ public:
virtual void solve()
{}
//- Update local data for topology changes
virtual void topoChange(const polyTopoChangeMap&);
//- Update from another mesh using the given map
virtual void mapMesh(const polyMeshMap&);
// Member Operators