Commit ddfe6c39 authored by mattijs's avatar mattijs
Browse files

user overrideable tolerances

parent 0d343e4d
......@@ -105,6 +105,8 @@ int main(int argc, char *argv[])
Foam::argList::validOptions.insert("overwrite", "");
Foam::argList::validOptions.insert("toleranceDict", "file with tolerances");
# include "setRootCase.H"
# include "createTime.H"
runTime.functionObjects().off();
......@@ -168,6 +170,22 @@ int main(int argc, char *argv[])
<< "If this is not the case use the -partial option" << nl << endl;
}
// set up the tolerances for the sliding mesh
dictionary slidingTolerances;
if (args.options().found("toleranceDict"))
{
IOdictionary toleranceFile(
IOobject(
args.options()["toleranceDict"],
runTime.constant(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
)
);
slidingTolerances += toleranceFile;
}
// Check for non-empty master and slave patches
checkPatch(mesh.boundaryMesh(), masterPatchName);
checkPatch(mesh.boundaryMesh(), slavePatchName);
......@@ -320,6 +338,11 @@ int main(int argc, char *argv[])
true // couple/decouple mode
)
);
static_cast<slidingInterface&>(stitcher[0]).setTolerances
(
slidingTolerances,
true
);
}
......
......@@ -43,7 +43,7 @@ License
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
const Foam::scalar Foam::slidingInterface::edgeCoPlanarTol_ = 0.8;
const Foam::scalar Foam::slidingInterface::edgeCoPlanarTolDefault_ = 0.8;
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
......
......@@ -27,10 +27,8 @@ License
#include "slidingInterface.H"
#include "polyTopoChanger.H"
#include "polyMesh.H"
#include "primitiveMesh.H"
#include "polyTopoChange.H"
#include "addToRunTimeSelectionTable.H"
#include "triPointRef.H"
#include "plane.H"
// Index of debug signs:
......@@ -173,6 +171,14 @@ Foam::slidingInterface::slidingInterface
attached_(false),
projectionAlgo_(algo),
trigger_(false),
pointMergeTol_(pointMergeTolDefault_),
edgeMergeTol_(edgeMergeTolDefault_),
nFacesPerSlaveEdge_(nFacesPerSlaveEdgeDefault_),
edgeFaceEscapeLimit_(edgeFaceEscapeLimitDefault_),
integralAdjTol_(integralAdjTolDefault_),
edgeMasterCatchFraction_(edgeMasterCatchFractionDefault_),
edgeCoPlanarTol_(edgeCoPlanarTolDefault_),
edgeEndCutoffTol_(edgeEndCutoffTolDefault_),
cutFaceMasterPtr_(NULL),
cutFaceSlavePtr_(NULL),
masterFaceCellsPtr_(NULL),
......@@ -280,6 +286,9 @@ Foam::slidingInterface::slidingInterface
masterPointEdgeHitsPtr_(NULL),
projectedSlavePointsPtr_(NULL)
{
// Optionally default tolerances from dictionary
setTolerances(dict);
checkDefinition();
// If the interface is attached, the master and slave face zone addressing
......@@ -686,6 +695,63 @@ const Foam::pointField& Foam::slidingInterface::pointProjection() const
return *projectedSlavePointsPtr_;
}
void Foam::slidingInterface::setTolerances(const dictionary&dict, bool report)
{
pointMergeTol_ = dict.lookupOrDefault<scalar>
(
"pointMergeTol",
pointMergeTol_
);
edgeMergeTol_ = dict.lookupOrDefault<scalar>
(
"edgeMergeTol",
edgeMergeTol_
);
nFacesPerSlaveEdge_ = dict.lookupOrDefault<scalar>
(
"nFacesPerSlaveEdge",
nFacesPerSlaveEdge_
);
edgeFaceEscapeLimit_ = dict.lookupOrDefault<scalar>
(
"edgeFaceEscapeLimit",
edgeFaceEscapeLimit_
);
integralAdjTol_ = dict.lookupOrDefault<scalar>
(
"integralAdjTol",
integralAdjTol_
);
edgeMasterCatchFraction_ = dict.lookupOrDefault<scalar>
(
"edgeMasterCatchFraction",
edgeMasterCatchFraction_
);
edgeCoPlanarTol_ = dict.lookupOrDefault<scalar>
(
"edgeCoPlanarTol",
edgeCoPlanarTol_
);
edgeEndCutoffTol_ = dict.lookupOrDefault<scalar>
(
"edgeEndCutoffTol",
edgeEndCutoffTol_
);
if (report)
{
Info<< "Sliding interface parameters:" << nl
<< "pointMergeTol : " << pointMergeTol_ << nl
<< "edgeMergeTol : " << edgeMergeTol_ << nl
<< "nFacesPerSlaveEdge : " << nFacesPerSlaveEdge_ << nl
<< "edgeFaceEscapeLimit : " << edgeFaceEscapeLimit_ << nl
<< "integralAdjTol : " << integralAdjTol_ << nl
<< "edgeMasterCatchFraction : " << edgeMasterCatchFraction_ << nl
<< "edgeCoPlanarTol : " << edgeCoPlanarTol_ << nl
<< "edgeEndCutoffTol : " << edgeEndCutoffTol_ << endl;
}
}
void Foam::slidingInterface::write(Ostream& os) const
{
......@@ -703,6 +769,14 @@ void Foam::slidingInterface::write(Ostream& os) const
}
// To write out all those tolerances
#define WRITE_NON_DEFAULT(name) \
if( name ## _ != name ## Default_ )\
{ \
os << " " #name " " << name ## _ << token::END_STATEMENT << nl; \
}
void Foam::slidingInterface::writeDict(Ostream& os) const
{
os << nl << name() << nl << token::BEGIN_BLOCK << nl
......@@ -743,6 +817,15 @@ void Foam::slidingInterface::writeDict(Ostream& os) const
<< token::END_STATEMENT << nl;
}
WRITE_NON_DEFAULT(pointMergeTol)
WRITE_NON_DEFAULT(edgeMergeTol)
WRITE_NON_DEFAULT(nFacesPerSlaveEdge)
WRITE_NON_DEFAULT(edgeFaceEscapeLimit)
WRITE_NON_DEFAULT(integralAdjTol)
WRITE_NON_DEFAULT(edgeMasterCatchFraction)
WRITE_NON_DEFAULT(edgeCoPlanarTol)
WRITE_NON_DEFAULT(edgeEndCutoffTol)
os << token::END_BLOCK << endl;
}
......
......@@ -129,6 +129,32 @@ private:
//- Trigger topological change
mutable bool trigger_;
// Tolerances. Initialised to static ones below.
//- Point merge tolerance
scalar pointMergeTol_;
//- Edge merge tolerance
scalar edgeMergeTol_;
//- Estimated number of faces an edge goes through
label nFacesPerSlaveEdge_;
//- Edge-face interaction escape limit
label edgeFaceEscapeLimit_;
//- Integral match point adjustment tolerance
scalar integralAdjTol_;
//- Edge intersection master catch fraction
scalar edgeMasterCatchFraction_;
//- Edge intersection co-planar tolerance
scalar edgeCoPlanarTol_;
//- Edge end cut-off tolerance
scalar edgeEndCutoffTol_;
// Private addressing data.
......@@ -256,28 +282,28 @@ private:
// Static data members
//- Point merge tolerance
static const scalar pointMergeTol_;
static const scalar pointMergeTolDefault_;
//- Edge merge tolerance
static const scalar edgeMergeTol_;
static const scalar edgeMergeTolDefault_;
//- Estimated number of faces an edge goes through
static const label nFacesPerSlaveEdge_;
static const label nFacesPerSlaveEdgeDefault_;
//- Edge-face interaction escape limit
static const label edgeFaceEscapeLimit_;
static const label edgeFaceEscapeLimitDefault_;
//- Integral match point adjustment tolerance
static const scalar integralAdjTol_;
static const scalar integralAdjTolDefault_;
//- Edge intersection master catch fraction
static const scalar edgeMasterCatchFraction_;
static const scalar edgeMasterCatchFractionDefault_;
//- Edge intersection co-planar tolerance
static const scalar edgeCoPlanarTol_;
static const scalar edgeCoPlanarTolDefault_;
//- Edge end cut-off tolerance
static const scalar edgeEndCutoffTol_;
static const scalar edgeEndCutoffTolDefault_;
public:
......@@ -350,6 +376,8 @@ public:
//- Return projected points for a slave patch
const pointField& pointProjection() const;
//- Set the tolerances from the values in a dictionary
void setTolerances(const dictionary&, bool report=false);
//- Write
virtual void write(Ostream&) const;
......
......@@ -26,22 +26,19 @@ License
#include "slidingInterface.H"
#include "polyMesh.H"
#include "primitiveMesh.H"
#include "line.H"
#include "triPointRef.H"
#include "plane.H"
#include "polyTopoChanger.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
const Foam::scalar Foam::slidingInterface::pointMergeTol_ = 0.05;
const Foam::scalar Foam::slidingInterface::edgeMergeTol_ = 0.01;
const Foam::label Foam::slidingInterface::nFacesPerSlaveEdge_ = 5;
const Foam::label Foam::slidingInterface::edgeFaceEscapeLimit_ = 10;
const Foam::scalar Foam::slidingInterface::pointMergeTolDefault_ = 0.05;
const Foam::scalar Foam::slidingInterface::edgeMergeTolDefault_ = 0.01;
const Foam::label Foam::slidingInterface::nFacesPerSlaveEdgeDefault_ = 5;
const Foam::label Foam::slidingInterface::edgeFaceEscapeLimitDefault_ = 10;
const Foam::scalar Foam::slidingInterface::integralAdjTol_ = 0.05;
const Foam::scalar Foam::slidingInterface::edgeMasterCatchFraction_ = 0.4;
const Foam::scalar Foam::slidingInterface::edgeEndCutoffTol_ = 0.0001;
const Foam::scalar Foam::slidingInterface::integralAdjTolDefault_ = 0.05;
const Foam::scalar Foam::slidingInterface::edgeMasterCatchFractionDefault_ = 0.4;
const Foam::scalar Foam::slidingInterface::edgeEndCutoffTolDefault_ = 0.0001;
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
......
Supports Markdown
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment