fvMeshTopoChanger::meshToMesh: Added support for repeating and cycling the set of meshes

Description
    fvMeshTopoChanger which maps the fields to a new mesh or sequence of meshes

    which can optionally be mapped to repeatedly for example in multi-cycle
    engine cases or cycled through for symmetric forward and reverse motion.

Usage
    \table
        Property  | Description                   | Required | Default value
        libs      | Libraries to load             | no       |
        times     | List of times for the meshes  | yes      |
        repeat    | Repetition period             | no       |
        cycle     | Cycle period                  | no       |
        begin     | Begin time for the meshes     | no       | Time::beginTime()
        timeDelta | Time tolerance used for time -> index | yes      |
    \endtable

    Examples of the mesh-to-mesh mapping for the multi-cycle
    tutorials/incompressibleFluid/movingCone case:
    \verbatim
    topoChanger
    {
        type    meshToMesh;

        libs    ("libmeshToMeshTopoChanger.so");

        times   (0.0015 0.003);

        cycle   #calc "1.0/300.0";
        begin   0;

        timeDelta 1e-6;
    }
    \endverbatim
This commit is contained in:
Henry Weller
2023-06-09 17:30:51 +01:00
parent 5233335924
commit 2a39553e28
6 changed files with 166 additions and 40 deletions

View File

@ -53,6 +53,44 @@ namespace fvMeshTopoChangers
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
bool Foam::fvMeshTopoChangers::meshToMesh::forward() const
{
return
!cycle_
|| int((mesh().time().userTimeValue() - begin_)/cycle_) % 2 == 0;
}
Foam::scalar Foam::fvMeshTopoChangers::meshToMesh::meshTime() const
{
const Time& time = mesh().time();
if (repeat_ > 0)
{
return begin_ + fmod(time.userTimeValue() - begin_, repeat_);
}
else if (cycle_ > 0)
{
if (forward())
{
return
begin_
+ fmod(time.userTimeValue() - begin_, cycle_);
}
else
{
return
begin_ + cycle_
- fmod(time.userTimeValue() - begin_, cycle_);
}
}
else
{
return time.userTimeValue();
}
}
void Foam::fvMeshTopoChangers::meshToMesh::interpolateUfs()
{
// Interpolate U to Uf
@ -88,6 +126,9 @@ Foam::fvMeshTopoChangers::meshToMesh::meshToMesh
dict_(dict),
times_(dict.lookup("times")),
timeDelta_(dict.lookup<scalar>("timeDelta")),
begin_(dict.lookupOrDefault("begin", mesh().time().beginTime().value())),
repeat_(dict.lookupOrDefault("repeat", 0.0)),
cycle_(dict.lookupOrDefault("cycle", 0.0)),
timeIndex_(-1)
{
forAll(times_, i)
@ -105,43 +146,82 @@ Foam::fvMeshTopoChangers::meshToMesh::~meshToMesh()
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::scalar Foam::fvMeshTopoChangers::meshToMesh::timeToNextMesh() const
{
const Time& time = mesh().time();
if (repeat_ || cycle_ || time.userTimeValue() + timeDelta_ < times_.last())
{
const scalar meshTime = this->meshTime();
if (!cycle_ || int((time.userTimeValue() - begin_)/cycle_) % 2 == 0)
{
forAll(times_, i)
{
if (times_[i] > meshTime + timeDelta_)
{
return time.userTimeToTime(times_[i] - meshTime);
}
}
}
else
{
forAllReverse(times_, i)
{
if (times_[i] < meshTime - timeDelta_)
{
return time.userTimeToTime(meshTime - times_[i]);
}
}
}
}
return vGreat;
}
bool Foam::fvMeshTopoChangers::meshToMesh::update()
{
const Time& time = mesh().time();
// Add the meshToMeshAdjustTimeStepFunctionObject functionObject
// if not already present
if
(
mesh().time().functionObjects().findObjectID("meshToMeshAdjustTimeStep")
time.functionObjects().findObjectID("meshToMeshAdjustTimeStep")
== -1
)
{
const_cast<Time&>(mesh().time()).functionObjects().append
const_cast<Time&>(time).functionObjects().append
(
new functionObjects::meshToMeshAdjustTimeStepFunctionObject
(
"meshToMeshAdjustTimeStep",
mesh().time(),
time,
dict_
)
);
}
// Only refine on the first call in a time-step
if (timeIndex_ != mesh().time().timeIndex())
if (timeIndex_ != time.timeIndex())
{
timeIndex_ = mesh().time().timeIndex();
timeIndex_ = time.timeIndex();
}
else
{
return false;
}
const scalar userTime = mesh().time().userTimeValue();
// Obtain the mesh time from the user time
// repeating or cycling the mesh sequence if required
if (timeIndices_.found((userTime + timeDelta_/2)/timeDelta_))
const scalar meshTime = this->meshTime();
if (timeIndices_.found((meshTime + timeDelta_/2)/timeDelta_))
{
const word otherMeshDir =
"meshToMesh_" + mesh().time().timeName(userTime);
"meshToMesh_" + time.timeName(meshTime);
Info << "Mapping to mesh " << otherMeshDir << endl;
@ -150,8 +230,8 @@ bool Foam::fvMeshTopoChangers::meshToMesh::update()
IOobject
(
otherMeshDir,
mesh().time().constant(),
mesh().time(),
time.constant(),
time,
IOobject::MUST_READ
),
false,

View File

@ -25,7 +25,39 @@ Class
Foam::fvMeshTopoChangers::meshToMesh
Description
fvMeshTopoChanger which maps to new mesh
fvMeshTopoChanger which maps the fields to a new mesh or sequence of meshes
which can optionally be mapped to repeatedly for example in multi-cycle
engine cases or cycled through for symmetric forward and reverse motion.
Usage
\table
Property | Description | Required | Default value
libs | Libraries to load | no |
times | List of times for the meshes | yes |
repeat | Repetition period | no |
cycle | Cycle period | no |
begin | Begin time for the meshes | no | Time::beginTime()
timeDelta | Time tolerance used for time -> index | yes |
\endtable
Examples of the mesh-to-mesh mapping for the multi-cycle
tutorials/incompressibleFluid/movingCone case:
\verbatim
topoChanger
{
type meshToMesh;
libs ("libmeshToMeshTopoChanger.so");
times (0.0015 0.003);
cycle #calc "1.0/300.0";
begin 0;
timeDelta 1e-6;
}
\endverbatim
SourceFiles
fvMeshTopoChangersMeshToMesh.C
@ -66,12 +98,30 @@ class meshToMesh
//- Hash set of mesh mapping time indices
HashSet<int64_t, Hash<int64_t>> timeIndices_;
//- Optional begin time for the meshes
scalar begin_;
//- Optional repetition period
scalar repeat_;
//- Optional cycle period
scalar cycle_;
//- The time index used for updating
label timeIndex_;
// Private Member Functions
//- Return true if the set of meshes are being traversed in the forward
// sequence or false if cycling is currently traversing the meshes in
// revers order
bool forward() const;
//- Return the user time mapped to the mesh sequence
// handling the repeat or cycle option
scalar meshTime() const;
//- Interpolate U's to Uf's
void interpolateUfs();
@ -97,15 +147,7 @@ public:
// Member Functions
const scalarList& times() const
{
return times_;
}
scalar timeDelta() const
{
return timeDelta_;
}
scalar timeToNextMesh() const;
//- Update the mesh for both mesh motion and topology change
virtual bool update();

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2022-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -85,20 +85,7 @@ bool Foam::functionObjects::meshToMeshAdjustTimeStepFunctionObject::read
Foam::scalar
Foam::functionObjects::meshToMeshAdjustTimeStepFunctionObject::timeToNextWrite()
{
const scalarList& times = meshToMesh_.times();
if (time_.userTimeValue() + meshToMesh_.timeDelta() < times.last())
{
forAll(times, i)
{
if (times[i] > time_.userTimeValue() + meshToMesh_.timeDelta())
{
return time_.userTimeToTime(times[i] - time_.userTimeValue());
}
}
}
return vGreat;
return meshToMesh_.timeToNextMesh();
}

View File

@ -17,20 +17,34 @@ dimensions [0 1 -1 0 0 0 0];
internalField uniform 0;
movement
{
type uniformFixedValue;
uniformValue
{
type square;
amplitude 1;
frequency 150;
start 0;
level 0;
}
}
boundaryField
{
#includeEtc "caseDicts/setConstraintTypes"
movingWall
{
type uniformFixedValue;
uniformValue constant 1;
$movement;
}
farFieldMoving
{
type uniformFixedValue;
uniformValue constant 1;
$movement;
}
fixedWall

View File

@ -35,6 +35,9 @@ topoChanger
times (0.0015 0.003);
cycle #calc "1.0/300.0";
begin 0;
timeDelta 1e-6;
}

View File

@ -24,7 +24,7 @@ startTime 0;
stopAt endTime;
endTime 0.0034;
endTime 0.01;
deltaT 5e-06;