mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: mapFields: expose AMI mapping methods
This commit is contained in:
@ -41,7 +41,8 @@ void mapConsistentMesh
|
|||||||
(
|
(
|
||||||
const fvMesh& meshSource,
|
const fvMesh& meshSource,
|
||||||
const fvMesh& meshTarget,
|
const fvMesh& meshTarget,
|
||||||
const meshToMesh::interpolationMethod& mapMethod,
|
const word& mapMethod,
|
||||||
|
const word& AMIMapMethod,
|
||||||
const bool subtract,
|
const bool subtract,
|
||||||
const HashSet<word>& selectedFields,
|
const HashSet<word>& selectedFields,
|
||||||
const bool noLagrangian
|
const bool noLagrangian
|
||||||
@ -50,7 +51,7 @@ void mapConsistentMesh
|
|||||||
Info<< nl << "Consistently creating and mapping fields for time "
|
Info<< nl << "Consistently creating and mapping fields for time "
|
||||||
<< meshSource.time().timeName() << nl << endl;
|
<< meshSource.time().timeName() << nl << endl;
|
||||||
|
|
||||||
meshToMesh interp(meshSource, meshTarget, mapMethod);
|
meshToMesh interp(meshSource, meshTarget, mapMethod, AMIMapMethod);
|
||||||
|
|
||||||
if (subtract)
|
if (subtract)
|
||||||
{
|
{
|
||||||
@ -79,7 +80,8 @@ void mapSubMesh
|
|||||||
const fvMesh& meshTarget,
|
const fvMesh& meshTarget,
|
||||||
const HashTable<word>& patchMap,
|
const HashTable<word>& patchMap,
|
||||||
const wordList& cuttingPatches,
|
const wordList& cuttingPatches,
|
||||||
const meshToMesh::interpolationMethod& mapMethod,
|
const word& mapMethod,
|
||||||
|
const word& AMIMapMethod,
|
||||||
const bool subtract,
|
const bool subtract,
|
||||||
const HashSet<word>& selectedFields,
|
const HashSet<word>& selectedFields,
|
||||||
const bool noLagrangian
|
const bool noLagrangian
|
||||||
@ -93,6 +95,7 @@ void mapSubMesh
|
|||||||
meshSource,
|
meshSource,
|
||||||
meshTarget,
|
meshTarget,
|
||||||
mapMethod,
|
mapMethod,
|
||||||
|
AMIMapMethod,
|
||||||
patchMap,
|
patchMap,
|
||||||
cuttingPatches
|
cuttingPatches
|
||||||
);
|
);
|
||||||
@ -186,6 +189,12 @@ int main(int argc, char *argv[])
|
|||||||
"word",
|
"word",
|
||||||
"specify the mapping method (direct|mapNearest|cellVolumeWeight)"
|
"specify the mapping method (direct|mapNearest|cellVolumeWeight)"
|
||||||
);
|
);
|
||||||
|
argList::addOption
|
||||||
|
(
|
||||||
|
"patchMapMethod",
|
||||||
|
"word",
|
||||||
|
"specify the patch mapping method (direct|mapNearest|faceAreaWeight)"
|
||||||
|
);
|
||||||
argList::addBoolOption
|
argList::addBoolOption
|
||||||
(
|
(
|
||||||
"subtract",
|
"subtract",
|
||||||
@ -231,17 +240,50 @@ int main(int argc, char *argv[])
|
|||||||
|
|
||||||
const bool consistent = args.optionFound("consistent");
|
const bool consistent = args.optionFound("consistent");
|
||||||
|
|
||||||
meshToMesh::interpolationMethod mapMethod =
|
|
||||||
meshToMesh::imCellVolumeWeight;
|
|
||||||
|
|
||||||
if (args.optionFound("mapMethod"))
|
word mapMethod = meshToMesh::interpolationMethodNames_
|
||||||
|
[
|
||||||
|
meshToMesh::imCellVolumeWeight
|
||||||
|
];
|
||||||
|
|
||||||
|
if (args.optionReadIfPresent("mapMethod", mapMethod))
|
||||||
{
|
{
|
||||||
mapMethod = meshToMesh::interpolationMethodNames_[args["mapMethod"]];
|
Info<< "Mapping method: " << mapMethod << endl;
|
||||||
|
|
||||||
Info<< "Mapping method: "
|
|
||||||
<< meshToMesh::interpolationMethodNames_[mapMethod] << endl;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
word patchMapMethod;
|
||||||
|
if (meshToMesh::interpolationMethodNames_.found(mapMethod))
|
||||||
|
{
|
||||||
|
// Lookup corresponding AMI method
|
||||||
|
meshToMesh::interpolationMethod method =
|
||||||
|
meshToMesh::interpolationMethodNames_[mapMethod];
|
||||||
|
|
||||||
|
patchMapMethod = AMIPatchToPatchInterpolation::interpolationMethodToWord
|
||||||
|
(
|
||||||
|
meshToMesh::interpolationMethodAMI(method)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
// Optionally override
|
||||||
|
if (args.optionFound("patchMapMethod"))
|
||||||
|
{
|
||||||
|
patchMapMethod = args["patchMapMethod"];
|
||||||
|
|
||||||
|
Info<< "Patch mapping method: " << patchMapMethod << endl;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
if (patchMapMethod.empty())
|
||||||
|
{
|
||||||
|
FatalErrorIn(args.executable())
|
||||||
|
<< "No valid patchMapMethod for method " << mapMethod
|
||||||
|
<< ". Please supply one through the 'patchMapMethod' option"
|
||||||
|
<< exit(FatalError);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
const bool subtract = args.optionFound("subtract");
|
const bool subtract = args.optionFound("subtract");
|
||||||
if (subtract)
|
if (subtract)
|
||||||
{
|
{
|
||||||
@ -314,6 +356,7 @@ int main(int argc, char *argv[])
|
|||||||
meshSource,
|
meshSource,
|
||||||
meshTarget,
|
meshTarget,
|
||||||
mapMethod,
|
mapMethod,
|
||||||
|
patchMapMethod,
|
||||||
subtract,
|
subtract,
|
||||||
selectedFields,
|
selectedFields,
|
||||||
noLagrangian
|
noLagrangian
|
||||||
@ -328,6 +371,7 @@ int main(int argc, char *argv[])
|
|||||||
patchMap,
|
patchMap,
|
||||||
addProcessorPatches(meshTarget, cuttingPatches),
|
addProcessorPatches(meshTarget, cuttingPatches),
|
||||||
mapMethod,
|
mapMethod,
|
||||||
|
patchMapMethod,
|
||||||
subtract,
|
subtract,
|
||||||
selectedFields,
|
selectedFields,
|
||||||
noLagrangian
|
noLagrangian
|
||||||
|
|||||||
@ -536,66 +536,13 @@ void Foam::AMIInterpolation<SourcePatch, TargetPatch>::agglomerate
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
||||||
|
|
||||||
template<class SourcePatch, class TargetPatch>
|
template<class SourcePatch, class TargetPatch>
|
||||||
Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
void Foam::AMIInterpolation<SourcePatch, TargetPatch>::constructFromSurface
|
||||||
(
|
(
|
||||||
const SourcePatch& srcPatch,
|
const SourcePatch& srcPatch,
|
||||||
const TargetPatch& tgtPatch,
|
const TargetPatch& tgtPatch,
|
||||||
const faceAreaIntersect::triangulationMode& triMode,
|
const autoPtr<searchableSurface>& surfPtr
|
||||||
const bool requireMatch,
|
|
||||||
const interpolationMethod& method,
|
|
||||||
const scalar lowWeightCorrection,
|
|
||||||
const bool reverseTarget
|
|
||||||
)
|
)
|
||||||
:
|
|
||||||
method_(method),
|
|
||||||
reverseTarget_(reverseTarget),
|
|
||||||
requireMatch_(requireMatch),
|
|
||||||
singlePatchProc_(-999),
|
|
||||||
lowWeightCorrection_(lowWeightCorrection),
|
|
||||||
srcAddress_(),
|
|
||||||
srcWeights_(),
|
|
||||||
srcWeightsSum_(),
|
|
||||||
tgtAddress_(),
|
|
||||||
tgtWeights_(),
|
|
||||||
tgtWeightsSum_(),
|
|
||||||
triMode_(triMode),
|
|
||||||
srcMapPtr_(NULL),
|
|
||||||
tgtMapPtr_(NULL)
|
|
||||||
{
|
|
||||||
update(srcPatch, tgtPatch);
|
|
||||||
}
|
|
||||||
|
|
||||||
|
|
||||||
template<class SourcePatch, class TargetPatch>
|
|
||||||
Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
|
||||||
(
|
|
||||||
const SourcePatch& srcPatch,
|
|
||||||
const TargetPatch& tgtPatch,
|
|
||||||
const autoPtr<searchableSurface>& surfPtr,
|
|
||||||
const faceAreaIntersect::triangulationMode& triMode,
|
|
||||||
const bool requireMatch,
|
|
||||||
const interpolationMethod& method,
|
|
||||||
const scalar lowWeightCorrection,
|
|
||||||
const bool reverseTarget
|
|
||||||
)
|
|
||||||
:
|
|
||||||
method_(method),
|
|
||||||
reverseTarget_(reverseTarget),
|
|
||||||
requireMatch_(requireMatch),
|
|
||||||
singlePatchProc_(-999),
|
|
||||||
lowWeightCorrection_(lowWeightCorrection),
|
|
||||||
srcAddress_(),
|
|
||||||
srcWeights_(),
|
|
||||||
srcWeightsSum_(),
|
|
||||||
tgtAddress_(),
|
|
||||||
tgtWeights_(),
|
|
||||||
tgtWeightsSum_(),
|
|
||||||
triMode_(triMode),
|
|
||||||
srcMapPtr_(NULL),
|
|
||||||
tgtMapPtr_(NULL)
|
|
||||||
{
|
{
|
||||||
if (surfPtr.valid())
|
if (surfPtr.valid())
|
||||||
{
|
{
|
||||||
@ -658,6 +605,134 @@ Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
template<class SourcePatch, class TargetPatch>
|
||||||
|
Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
||||||
|
(
|
||||||
|
const SourcePatch& srcPatch,
|
||||||
|
const TargetPatch& tgtPatch,
|
||||||
|
const faceAreaIntersect::triangulationMode& triMode,
|
||||||
|
const bool requireMatch,
|
||||||
|
const interpolationMethod& method,
|
||||||
|
const scalar lowWeightCorrection,
|
||||||
|
const bool reverseTarget
|
||||||
|
)
|
||||||
|
:
|
||||||
|
methodName_(interpolationMethodToWord(method)),
|
||||||
|
reverseTarget_(reverseTarget),
|
||||||
|
requireMatch_(requireMatch),
|
||||||
|
singlePatchProc_(-999),
|
||||||
|
lowWeightCorrection_(lowWeightCorrection),
|
||||||
|
srcAddress_(),
|
||||||
|
srcWeights_(),
|
||||||
|
srcWeightsSum_(),
|
||||||
|
tgtAddress_(),
|
||||||
|
tgtWeights_(),
|
||||||
|
tgtWeightsSum_(),
|
||||||
|
triMode_(triMode),
|
||||||
|
srcMapPtr_(NULL),
|
||||||
|
tgtMapPtr_(NULL)
|
||||||
|
{
|
||||||
|
update(srcPatch, tgtPatch);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class SourcePatch, class TargetPatch>
|
||||||
|
Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
||||||
|
(
|
||||||
|
const SourcePatch& srcPatch,
|
||||||
|
const TargetPatch& tgtPatch,
|
||||||
|
const faceAreaIntersect::triangulationMode& triMode,
|
||||||
|
const bool requireMatch,
|
||||||
|
const word& methodName,
|
||||||
|
const scalar lowWeightCorrection,
|
||||||
|
const bool reverseTarget
|
||||||
|
)
|
||||||
|
:
|
||||||
|
methodName_(methodName),
|
||||||
|
reverseTarget_(reverseTarget),
|
||||||
|
requireMatch_(requireMatch),
|
||||||
|
singlePatchProc_(-999),
|
||||||
|
lowWeightCorrection_(lowWeightCorrection),
|
||||||
|
srcAddress_(),
|
||||||
|
srcWeights_(),
|
||||||
|
srcWeightsSum_(),
|
||||||
|
tgtAddress_(),
|
||||||
|
tgtWeights_(),
|
||||||
|
tgtWeightsSum_(),
|
||||||
|
triMode_(triMode),
|
||||||
|
srcMapPtr_(NULL),
|
||||||
|
tgtMapPtr_(NULL)
|
||||||
|
{
|
||||||
|
update(srcPatch, tgtPatch);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class SourcePatch, class TargetPatch>
|
||||||
|
Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
||||||
|
(
|
||||||
|
const SourcePatch& srcPatch,
|
||||||
|
const TargetPatch& tgtPatch,
|
||||||
|
const autoPtr<searchableSurface>& surfPtr,
|
||||||
|
const faceAreaIntersect::triangulationMode& triMode,
|
||||||
|
const bool requireMatch,
|
||||||
|
const interpolationMethod& method,
|
||||||
|
const scalar lowWeightCorrection,
|
||||||
|
const bool reverseTarget
|
||||||
|
)
|
||||||
|
:
|
||||||
|
methodName_(interpolationMethodToWord(method)),
|
||||||
|
reverseTarget_(reverseTarget),
|
||||||
|
requireMatch_(requireMatch),
|
||||||
|
singlePatchProc_(-999),
|
||||||
|
lowWeightCorrection_(lowWeightCorrection),
|
||||||
|
srcAddress_(),
|
||||||
|
srcWeights_(),
|
||||||
|
srcWeightsSum_(),
|
||||||
|
tgtAddress_(),
|
||||||
|
tgtWeights_(),
|
||||||
|
tgtWeightsSum_(),
|
||||||
|
triMode_(triMode),
|
||||||
|
srcMapPtr_(NULL),
|
||||||
|
tgtMapPtr_(NULL)
|
||||||
|
{
|
||||||
|
constructFromSurface(srcPatch, tgtPatch, surfPtr);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
template<class SourcePatch, class TargetPatch>
|
||||||
|
Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
||||||
|
(
|
||||||
|
const SourcePatch& srcPatch,
|
||||||
|
const TargetPatch& tgtPatch,
|
||||||
|
const autoPtr<searchableSurface>& surfPtr,
|
||||||
|
const faceAreaIntersect::triangulationMode& triMode,
|
||||||
|
const bool requireMatch,
|
||||||
|
const word& methodName,
|
||||||
|
const scalar lowWeightCorrection,
|
||||||
|
const bool reverseTarget
|
||||||
|
)
|
||||||
|
:
|
||||||
|
methodName_(methodName),
|
||||||
|
reverseTarget_(reverseTarget),
|
||||||
|
requireMatch_(requireMatch),
|
||||||
|
singlePatchProc_(-999),
|
||||||
|
lowWeightCorrection_(lowWeightCorrection),
|
||||||
|
srcAddress_(),
|
||||||
|
srcWeights_(),
|
||||||
|
srcWeightsSum_(),
|
||||||
|
tgtAddress_(),
|
||||||
|
tgtWeights_(),
|
||||||
|
tgtWeightsSum_(),
|
||||||
|
triMode_(triMode),
|
||||||
|
srcMapPtr_(NULL),
|
||||||
|
tgtMapPtr_(NULL)
|
||||||
|
{
|
||||||
|
constructFromSurface(srcPatch, tgtPatch, surfPtr);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
template<class SourcePatch, class TargetPatch>
|
template<class SourcePatch, class TargetPatch>
|
||||||
Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
||||||
(
|
(
|
||||||
@ -666,7 +741,7 @@ Foam::AMIInterpolation<SourcePatch, TargetPatch>::AMIInterpolation
|
|||||||
const labelList& targetRestrictAddressing
|
const labelList& targetRestrictAddressing
|
||||||
)
|
)
|
||||||
:
|
:
|
||||||
method_(fineAMI.method_),
|
methodName_(fineAMI.methodName_),
|
||||||
reverseTarget_(fineAMI.reverseTarget_),
|
reverseTarget_(fineAMI.reverseTarget_),
|
||||||
requireMatch_(fineAMI.requireMatch_),
|
requireMatch_(fineAMI.requireMatch_),
|
||||||
singlePatchProc_(fineAMI.singlePatchProc_),
|
singlePatchProc_(fineAMI.singlePatchProc_),
|
||||||
@ -883,7 +958,7 @@ void Foam::AMIInterpolation<SourcePatch, TargetPatch>::update
|
|||||||
(
|
(
|
||||||
AMIMethod<SourcePatch, TargetPatch>::New
|
AMIMethod<SourcePatch, TargetPatch>::New
|
||||||
(
|
(
|
||||||
interpolationMethodToWord(method_),
|
methodName_,
|
||||||
srcPatch,
|
srcPatch,
|
||||||
newTgtPatch,
|
newTgtPatch,
|
||||||
srcMagSf_,
|
srcMagSf_,
|
||||||
@ -1000,7 +1075,7 @@ void Foam::AMIInterpolation<SourcePatch, TargetPatch>::update
|
|||||||
(
|
(
|
||||||
AMIMethod<SourcePatch, TargetPatch>::New
|
AMIMethod<SourcePatch, TargetPatch>::New
|
||||||
(
|
(
|
||||||
interpolationMethodToWord(method_),
|
methodName_,
|
||||||
srcPatch,
|
srcPatch,
|
||||||
tgtPatch,
|
tgtPatch,
|
||||||
srcMagSf_,
|
srcMagSf_,
|
||||||
|
|||||||
@ -110,7 +110,7 @@ private:
|
|||||||
// Private data
|
// Private data
|
||||||
|
|
||||||
//- Interpolation method
|
//- Interpolation method
|
||||||
interpolationMethod method_;
|
const word methodName_;
|
||||||
|
|
||||||
//- Flag to indicate that the two patches are co-directional and
|
//- Flag to indicate that the two patches are co-directional and
|
||||||
// that the orientation of the target patch should be reversed
|
// that the orientation of the target patch should be reversed
|
||||||
@ -250,7 +250,7 @@ private:
|
|||||||
);
|
);
|
||||||
|
|
||||||
|
|
||||||
// Constructor helper
|
// Constructor helpers
|
||||||
|
|
||||||
static void agglomerate
|
static void agglomerate
|
||||||
(
|
(
|
||||||
@ -269,6 +269,12 @@ private:
|
|||||||
autoPtr<mapDistribute>& tgtMap
|
autoPtr<mapDistribute>& tgtMap
|
||||||
);
|
);
|
||||||
|
|
||||||
|
void constructFromSurface
|
||||||
|
(
|
||||||
|
const SourcePatch& srcPatch,
|
||||||
|
const TargetPatch& tgtPatch,
|
||||||
|
const autoPtr<searchableSurface>& surfPtr
|
||||||
|
);
|
||||||
|
|
||||||
public:
|
public:
|
||||||
|
|
||||||
@ -286,6 +292,19 @@ public:
|
|||||||
const bool reverseTarget = false
|
const bool reverseTarget = false
|
||||||
);
|
);
|
||||||
|
|
||||||
|
//- Construct from components
|
||||||
|
AMIInterpolation
|
||||||
|
(
|
||||||
|
const SourcePatch& srcPatch,
|
||||||
|
const TargetPatch& tgtPatch,
|
||||||
|
const faceAreaIntersect::triangulationMode& triMode,
|
||||||
|
const bool requireMatch = true,
|
||||||
|
const word& methodName =
|
||||||
|
interpolationMethodToWord(imFaceAreaWeight),
|
||||||
|
const scalar lowWeightCorrection = -1,
|
||||||
|
const bool reverseTarget = false
|
||||||
|
);
|
||||||
|
|
||||||
//- Construct from components, with projection surface
|
//- Construct from components, with projection surface
|
||||||
AMIInterpolation
|
AMIInterpolation
|
||||||
(
|
(
|
||||||
@ -299,6 +318,20 @@ public:
|
|||||||
const bool reverseTarget = false
|
const bool reverseTarget = false
|
||||||
);
|
);
|
||||||
|
|
||||||
|
//- Construct from components, with projection surface
|
||||||
|
AMIInterpolation
|
||||||
|
(
|
||||||
|
const SourcePatch& srcPatch,
|
||||||
|
const TargetPatch& tgtPatch,
|
||||||
|
const autoPtr<searchableSurface>& surf,
|
||||||
|
const faceAreaIntersect::triangulationMode& triMode,
|
||||||
|
const bool requireMatch = true,
|
||||||
|
const word& methodName =
|
||||||
|
interpolationMethodToWord(imFaceAreaWeight),
|
||||||
|
const scalar lowWeightCorrection = -1,
|
||||||
|
const bool reverseTarget = false
|
||||||
|
);
|
||||||
|
|
||||||
//- Construct from agglomeration of AMIInterpolation. Agglomeration
|
//- Construct from agglomeration of AMIInterpolation. Agglomeration
|
||||||
// passed in as new coarse size and addressing from fine from coarse
|
// passed in as new coarse size and addressing from fine from coarse
|
||||||
AMIInterpolation
|
AMIInterpolation
|
||||||
|
|||||||
@ -119,6 +119,7 @@ void Foam::meshToMesh::normaliseWeights
|
|||||||
|
|
||||||
void Foam::meshToMesh::calcAddressing
|
void Foam::meshToMesh::calcAddressing
|
||||||
(
|
(
|
||||||
|
const word& methodName,
|
||||||
const polyMesh& src,
|
const polyMesh& src,
|
||||||
const polyMesh& tgt
|
const polyMesh& tgt
|
||||||
)
|
)
|
||||||
@ -127,7 +128,7 @@ void Foam::meshToMesh::calcAddressing
|
|||||||
(
|
(
|
||||||
meshToMeshMethod::New
|
meshToMeshMethod::New
|
||||||
(
|
(
|
||||||
interpolationMethodNames_[method_],
|
methodName,
|
||||||
src,
|
src,
|
||||||
tgt
|
tgt
|
||||||
)
|
)
|
||||||
@ -150,11 +151,11 @@ void Foam::meshToMesh::calcAddressing
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
void Foam::meshToMesh::calculate()
|
void Foam::meshToMesh::calculate(const word& methodName)
|
||||||
{
|
{
|
||||||
Info<< "Creating mesh-to-mesh addressing for " << srcRegion_.name()
|
Info<< "Creating mesh-to-mesh addressing for " << srcRegion_.name()
|
||||||
<< " and " << tgtRegion_.name() << " regions using "
|
<< " and " << tgtRegion_.name() << " regions using "
|
||||||
<< interpolationMethodNames_[method_] << endl;
|
<< methodName << endl;
|
||||||
|
|
||||||
singleMeshProc_ = calcDistribution(srcRegion_, tgtRegion_);
|
singleMeshProc_ = calcDistribution(srcRegion_, tgtRegion_);
|
||||||
|
|
||||||
@ -241,7 +242,7 @@ void Foam::meshToMesh::calculate()
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
calcAddressing(srcRegion_, newTgt);
|
calcAddressing(methodName, srcRegion_, newTgt);
|
||||||
|
|
||||||
// per source cell the target cell address in newTgt mesh
|
// per source cell the target cell address in newTgt mesh
|
||||||
forAll(srcToTgtCellAddr_, i)
|
forAll(srcToTgtCellAddr_, i)
|
||||||
@ -320,7 +321,7 @@ void Foam::meshToMesh::calculate()
|
|||||||
}
|
}
|
||||||
else
|
else
|
||||||
{
|
{
|
||||||
calcAddressing(srcRegion_, tgtRegion_);
|
calcAddressing(methodName, srcRegion_, tgtRegion_);
|
||||||
|
|
||||||
normaliseWeights
|
normaliseWeights
|
||||||
(
|
(
|
||||||
@ -342,12 +343,9 @@ void Foam::meshToMesh::calculate()
|
|||||||
|
|
||||||
|
|
||||||
Foam::AMIPatchToPatchInterpolation::interpolationMethod
|
Foam::AMIPatchToPatchInterpolation::interpolationMethod
|
||||||
Foam::meshToMesh::interpolationMethodAMI
|
Foam::meshToMesh::interpolationMethodAMI(const interpolationMethod method)
|
||||||
(
|
|
||||||
const interpolationMethod method
|
|
||||||
) const
|
|
||||||
{
|
{
|
||||||
switch (method_)
|
switch (method)
|
||||||
{
|
{
|
||||||
case imDirect:
|
case imDirect:
|
||||||
{
|
{
|
||||||
@ -374,7 +372,7 @@ Foam::meshToMesh::interpolationMethodAMI
|
|||||||
"const interpolationMethod method"
|
"const interpolationMethod method"
|
||||||
") const"
|
") const"
|
||||||
)
|
)
|
||||||
<< "Unhandled enumeration " << method_
|
<< "Unhandled enumeration " << method
|
||||||
<< abort(FatalError);
|
<< abort(FatalError);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -383,90 +381,66 @@ Foam::meshToMesh::interpolationMethodAMI
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
const Foam::PtrList<Foam::AMIPatchToPatchInterpolation>&
|
void Foam::meshToMesh::calculatePatchAMIs(const word& AMIMethodName)
|
||||||
Foam::meshToMesh::patchAMIs() const
|
|
||||||
{
|
{
|
||||||
if (patchAMIs_.empty())
|
if (!patchAMIs_.empty())
|
||||||
{
|
{
|
||||||
const word amiMethod =
|
FatalErrorIn("meshToMesh::calculatePatchAMIs()")
|
||||||
AMIPatchToPatchInterpolation::interpolationMethodToWord
|
<< "patch AMI already calculated"
|
||||||
(
|
<< exit(FatalError);
|
||||||
interpolationMethodAMI(method_)
|
|
||||||
);
|
|
||||||
|
|
||||||
patchAMIs_.setSize(srcPatchID_.size());
|
|
||||||
|
|
||||||
forAll(srcPatchID_, i)
|
|
||||||
{
|
|
||||||
label srcPatchI = srcPatchID_[i];
|
|
||||||
label tgtPatchI = tgtPatchID_[i];
|
|
||||||
|
|
||||||
const polyPatch& srcPP = srcRegion_.boundaryMesh()[srcPatchI];
|
|
||||||
const polyPatch& tgtPP = tgtRegion_.boundaryMesh()[tgtPatchI];
|
|
||||||
|
|
||||||
Info<< "Creating AMI between source patch " << srcPP.name()
|
|
||||||
<< " and target patch " << tgtPP.name()
|
|
||||||
<< " using " << amiMethod
|
|
||||||
<< endl;
|
|
||||||
|
|
||||||
Info<< incrIndent;
|
|
||||||
|
|
||||||
patchAMIs_.set
|
|
||||||
(
|
|
||||||
i,
|
|
||||||
new AMIPatchToPatchInterpolation
|
|
||||||
(
|
|
||||||
srcPP,
|
|
||||||
tgtPP,
|
|
||||||
faceAreaIntersect::tmMesh,
|
|
||||||
false,
|
|
||||||
interpolationMethodAMI(method_),
|
|
||||||
-1,
|
|
||||||
true // flip target patch since patch normals are aligned
|
|
||||||
)
|
|
||||||
);
|
|
||||||
|
|
||||||
Info<< decrIndent;
|
|
||||||
}
|
|
||||||
}
|
}
|
||||||
|
|
||||||
return patchAMIs_;
|
patchAMIs_.setSize(srcPatchID_.size());
|
||||||
|
|
||||||
|
forAll(srcPatchID_, i)
|
||||||
|
{
|
||||||
|
label srcPatchI = srcPatchID_[i];
|
||||||
|
label tgtPatchI = tgtPatchID_[i];
|
||||||
|
|
||||||
|
const polyPatch& srcPP = srcRegion_.boundaryMesh()[srcPatchI];
|
||||||
|
const polyPatch& tgtPP = tgtRegion_.boundaryMesh()[tgtPatchI];
|
||||||
|
|
||||||
|
Info<< "Creating AMI between source patch " << srcPP.name()
|
||||||
|
<< " and target patch " << tgtPP.name()
|
||||||
|
<< " using " << AMIMethodName
|
||||||
|
<< endl;
|
||||||
|
|
||||||
|
Info<< incrIndent;
|
||||||
|
|
||||||
|
patchAMIs_.set
|
||||||
|
(
|
||||||
|
i,
|
||||||
|
new AMIPatchToPatchInterpolation
|
||||||
|
(
|
||||||
|
srcPP,
|
||||||
|
tgtPP,
|
||||||
|
faceAreaIntersect::tmMesh,
|
||||||
|
false,
|
||||||
|
AMIMethodName,
|
||||||
|
-1,
|
||||||
|
true // flip target patch since patch normals are aligned
|
||||||
|
)
|
||||||
|
);
|
||||||
|
|
||||||
|
Info<< decrIndent;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
void Foam::meshToMesh::constructNoCuttingPatches
|
||||||
|
|
||||||
Foam::meshToMesh::meshToMesh
|
|
||||||
(
|
(
|
||||||
const polyMesh& src,
|
const word& methodName,
|
||||||
const polyMesh& tgt,
|
const word& AMIMethodName,
|
||||||
const interpolationMethod& method,
|
const bool interpAllPatches
|
||||||
bool interpAllPatches
|
|
||||||
)
|
)
|
||||||
:
|
|
||||||
srcRegion_(src),
|
|
||||||
tgtRegion_(tgt),
|
|
||||||
srcPatchID_(),
|
|
||||||
tgtPatchID_(),
|
|
||||||
patchAMIs_(),
|
|
||||||
cuttingPatches_(),
|
|
||||||
srcToTgtCellAddr_(),
|
|
||||||
tgtToSrcCellAddr_(),
|
|
||||||
srcToTgtCellWght_(),
|
|
||||||
tgtToSrcCellWght_(),
|
|
||||||
method_(method),
|
|
||||||
V_(0.0),
|
|
||||||
singleMeshProc_(-1),
|
|
||||||
srcMapPtr_(NULL),
|
|
||||||
tgtMapPtr_(NULL)
|
|
||||||
{
|
{
|
||||||
if (interpAllPatches)
|
if (interpAllPatches)
|
||||||
{
|
{
|
||||||
const polyBoundaryMesh& srcBM = src.boundaryMesh();
|
const polyBoundaryMesh& srcBM = srcRegion_.boundaryMesh();
|
||||||
const polyBoundaryMesh& tgtBM = tgt.boundaryMesh();
|
const polyBoundaryMesh& tgtBM = tgtRegion_.boundaryMesh();
|
||||||
|
|
||||||
DynamicList<label> srcPatchID(src.boundaryMesh().size());
|
DynamicList<label> srcPatchID(srcBM.size());
|
||||||
DynamicList<label> tgtPatchID(tgt.boundaryMesh().size());
|
DynamicList<label> tgtPatchID(tgtBM.size());
|
||||||
forAll(srcBM, patchI)
|
forAll(srcBM, patchI)
|
||||||
{
|
{
|
||||||
const polyPatch& pp = srcBM[patchI];
|
const polyPatch& pp = srcBM[patchI];
|
||||||
@ -474,7 +448,7 @@ Foam::meshToMesh::meshToMesh
|
|||||||
{
|
{
|
||||||
srcPatchID.append(pp.index());
|
srcPatchID.append(pp.index());
|
||||||
|
|
||||||
label tgtPatchI = tgt.boundaryMesh().findPatchID(pp.name());
|
label tgtPatchI = tgtBM.findPatchID(pp.name());
|
||||||
|
|
||||||
if (tgtPatchI != -1)
|
if (tgtPatchI != -1)
|
||||||
{
|
{
|
||||||
@ -504,10 +478,116 @@ Foam::meshToMesh::meshToMesh
|
|||||||
}
|
}
|
||||||
|
|
||||||
// calculate volume addressing and weights
|
// calculate volume addressing and weights
|
||||||
calculate();
|
calculate(methodName);
|
||||||
|
|
||||||
// calculate patch addressing and weights
|
// calculate patch addressing and weights
|
||||||
(void)patchAMIs();
|
calculatePatchAMIs(AMIMethodName);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
void Foam::meshToMesh::constructFromCuttingPatches
|
||||||
|
(
|
||||||
|
const word& methodName,
|
||||||
|
const word& AMIMethodName,
|
||||||
|
const HashTable<word>& patchMap,
|
||||||
|
const wordList& cuttingPatches
|
||||||
|
)
|
||||||
|
{
|
||||||
|
srcPatchID_.setSize(patchMap.size());
|
||||||
|
tgtPatchID_.setSize(patchMap.size());
|
||||||
|
|
||||||
|
label i = 0;
|
||||||
|
forAllConstIter(HashTable<word>, patchMap, iter)
|
||||||
|
{
|
||||||
|
const word& tgtPatchName = iter.key();
|
||||||
|
const word& srcPatchName = iter();
|
||||||
|
|
||||||
|
const polyPatch& srcPatch = srcRegion_.boundaryMesh()[srcPatchName];
|
||||||
|
const polyPatch& tgtPatch = tgtRegion_.boundaryMesh()[tgtPatchName];
|
||||||
|
|
||||||
|
srcPatchID_[i] = srcPatch.index();
|
||||||
|
tgtPatchID_[i] = tgtPatch.index();
|
||||||
|
i++;
|
||||||
|
}
|
||||||
|
|
||||||
|
// calculate volume addressing and weights
|
||||||
|
calculate(methodName);
|
||||||
|
|
||||||
|
// calculate patch addressing and weights
|
||||||
|
calculatePatchAMIs(AMIMethodName);
|
||||||
|
|
||||||
|
// set IDs of cutting patches on target mesh
|
||||||
|
cuttingPatches_.setSize(cuttingPatches.size());
|
||||||
|
forAll(cuttingPatches_, i)
|
||||||
|
{
|
||||||
|
const word& patchName = cuttingPatches[i];
|
||||||
|
cuttingPatches_[i] = tgtRegion_.boundaryMesh().findPatchID(patchName);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
Foam::meshToMesh::meshToMesh
|
||||||
|
(
|
||||||
|
const polyMesh& src,
|
||||||
|
const polyMesh& tgt,
|
||||||
|
const interpolationMethod& method,
|
||||||
|
bool interpAllPatches
|
||||||
|
)
|
||||||
|
:
|
||||||
|
srcRegion_(src),
|
||||||
|
tgtRegion_(tgt),
|
||||||
|
srcPatchID_(),
|
||||||
|
tgtPatchID_(),
|
||||||
|
patchAMIs_(),
|
||||||
|
cuttingPatches_(),
|
||||||
|
srcToTgtCellAddr_(),
|
||||||
|
tgtToSrcCellAddr_(),
|
||||||
|
srcToTgtCellWght_(),
|
||||||
|
tgtToSrcCellWght_(),
|
||||||
|
V_(0.0),
|
||||||
|
singleMeshProc_(-1),
|
||||||
|
srcMapPtr_(NULL),
|
||||||
|
tgtMapPtr_(NULL)
|
||||||
|
{
|
||||||
|
constructNoCuttingPatches
|
||||||
|
(
|
||||||
|
interpolationMethodNames_[method],
|
||||||
|
AMIPatchToPatchInterpolation::interpolationMethodToWord
|
||||||
|
(
|
||||||
|
interpolationMethodAMI(method)
|
||||||
|
),
|
||||||
|
interpAllPatches
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::meshToMesh::meshToMesh
|
||||||
|
(
|
||||||
|
const polyMesh& src,
|
||||||
|
const polyMesh& tgt,
|
||||||
|
const word& methodName,
|
||||||
|
const word& AMIMethodName,
|
||||||
|
bool interpAllPatches
|
||||||
|
)
|
||||||
|
:
|
||||||
|
srcRegion_(src),
|
||||||
|
tgtRegion_(tgt),
|
||||||
|
srcPatchID_(),
|
||||||
|
tgtPatchID_(),
|
||||||
|
patchAMIs_(),
|
||||||
|
cuttingPatches_(),
|
||||||
|
srcToTgtCellAddr_(),
|
||||||
|
tgtToSrcCellAddr_(),
|
||||||
|
srcToTgtCellWght_(),
|
||||||
|
tgtToSrcCellWght_(),
|
||||||
|
V_(0.0),
|
||||||
|
singleMeshProc_(-1),
|
||||||
|
srcMapPtr_(NULL),
|
||||||
|
tgtMapPtr_(NULL)
|
||||||
|
{
|
||||||
|
constructNoCuttingPatches(methodName, AMIMethodName, interpAllPatches);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -530,42 +610,56 @@ Foam::meshToMesh::meshToMesh
|
|||||||
tgtToSrcCellAddr_(),
|
tgtToSrcCellAddr_(),
|
||||||
srcToTgtCellWght_(),
|
srcToTgtCellWght_(),
|
||||||
tgtToSrcCellWght_(),
|
tgtToSrcCellWght_(),
|
||||||
method_(method),
|
|
||||||
V_(0.0),
|
V_(0.0),
|
||||||
singleMeshProc_(-1),
|
singleMeshProc_(-1),
|
||||||
srcMapPtr_(NULL),
|
srcMapPtr_(NULL),
|
||||||
tgtMapPtr_(NULL)
|
tgtMapPtr_(NULL)
|
||||||
{
|
{
|
||||||
srcPatchID_.setSize(patchMap.size());
|
constructFromCuttingPatches
|
||||||
tgtPatchID_.setSize(patchMap.size());
|
(
|
||||||
|
interpolationMethodNames_[method],
|
||||||
|
AMIPatchToPatchInterpolation::interpolationMethodToWord
|
||||||
|
(
|
||||||
|
interpolationMethodAMI(method)
|
||||||
|
),
|
||||||
|
patchMap,
|
||||||
|
cuttingPatches
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
label i = 0;
|
|
||||||
forAllConstIter(HashTable<word>, patchMap, iter)
|
|
||||||
{
|
|
||||||
const word& tgtPatchName = iter.key();
|
|
||||||
const word& srcPatchName = iter();
|
|
||||||
|
|
||||||
const polyPatch& srcPatch = srcRegion_.boundaryMesh()[srcPatchName];
|
Foam::meshToMesh::meshToMesh
|
||||||
const polyPatch& tgtPatch = tgtRegion_.boundaryMesh()[tgtPatchName];
|
(
|
||||||
|
const polyMesh& src,
|
||||||
srcPatchID_[i] = srcPatch.index();
|
const polyMesh& tgt,
|
||||||
tgtPatchID_[i] = tgtPatch.index();
|
const word& methodName, // internal mapping
|
||||||
i++;
|
const word& AMIMethodName, // boundary mapping
|
||||||
}
|
const HashTable<word>& patchMap,
|
||||||
|
const wordList& cuttingPatches
|
||||||
// calculate volume addressing and weights
|
)
|
||||||
calculate();
|
:
|
||||||
|
srcRegion_(src),
|
||||||
// calculate patch addressing and weights
|
tgtRegion_(tgt),
|
||||||
(void)patchAMIs();
|
srcPatchID_(),
|
||||||
|
tgtPatchID_(),
|
||||||
// set IDs of cutting patches on target mesh
|
patchAMIs_(),
|
||||||
cuttingPatches_.setSize(cuttingPatches.size());
|
cuttingPatches_(),
|
||||||
forAll(cuttingPatches_, i)
|
srcToTgtCellAddr_(),
|
||||||
{
|
tgtToSrcCellAddr_(),
|
||||||
const word& patchName = cuttingPatches[i];
|
srcToTgtCellWght_(),
|
||||||
cuttingPatches_[i] = tgt.boundaryMesh().findPatchID(patchName);
|
tgtToSrcCellWght_(),
|
||||||
}
|
V_(0.0),
|
||||||
|
singleMeshProc_(-1),
|
||||||
|
srcMapPtr_(NULL),
|
||||||
|
tgtMapPtr_(NULL)
|
||||||
|
{
|
||||||
|
constructFromCuttingPatches
|
||||||
|
(
|
||||||
|
methodName,
|
||||||
|
AMIMethodName,
|
||||||
|
patchMap,
|
||||||
|
cuttingPatches
|
||||||
|
);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|||||||
@ -92,7 +92,7 @@ private:
|
|||||||
List<label> tgtPatchID_;
|
List<label> tgtPatchID_;
|
||||||
|
|
||||||
//- List of AMIs between source and target patches
|
//- List of AMIs between source and target patches
|
||||||
mutable PtrList<AMIPatchToPatchInterpolation> patchAMIs_;
|
PtrList<AMIPatchToPatchInterpolation> patchAMIs_;
|
||||||
|
|
||||||
//- Cutting patches whose values are set using a zero-gradient condition
|
//- Cutting patches whose values are set using a zero-gradient condition
|
||||||
List<label> cuttingPatches_;
|
List<label> cuttingPatches_;
|
||||||
@ -109,9 +109,6 @@ private:
|
|||||||
//- Target to source cell interpolation weights
|
//- Target to source cell interpolation weights
|
||||||
scalarListList tgtToSrcCellWght_;
|
scalarListList tgtToSrcCellWght_;
|
||||||
|
|
||||||
//- Interpolation method
|
|
||||||
interpolationMethod method_;
|
|
||||||
|
|
||||||
//- Cell total volume in overlap region [m3]
|
//- Cell total volume in overlap region [m3]
|
||||||
scalar V_;
|
scalar V_;
|
||||||
|
|
||||||
@ -143,22 +140,41 @@ private:
|
|||||||
scalarListList& wght
|
scalarListList& wght
|
||||||
) const;
|
) const;
|
||||||
|
|
||||||
//- Calculate the addressing between overalping regions of src and tgt
|
//- Calculate the addressing between overlapping regions of src and tgt
|
||||||
// meshes
|
// meshes
|
||||||
void calcAddressing(const polyMesh& src, const polyMesh& tgt);
|
void calcAddressing
|
||||||
|
(
|
||||||
|
const word& methodName,
|
||||||
|
const polyMesh& src,
|
||||||
|
const polyMesh& tgt
|
||||||
|
);
|
||||||
|
|
||||||
//- Calculate - main driver function
|
//- Calculate - main driver function
|
||||||
void calculate();
|
void calculate(const word& methodName);
|
||||||
|
|
||||||
//- Conversion between mesh and patch interpolation methods
|
//- Calculate patch overlap
|
||||||
AMIPatchToPatchInterpolation::interpolationMethod
|
void calculatePatchAMIs(const word& amiMethodName);
|
||||||
interpolationMethodAMI
|
|
||||||
|
//- Constructor helper
|
||||||
|
void constructNoCuttingPatches
|
||||||
(
|
(
|
||||||
const interpolationMethod method
|
const word& methodName,
|
||||||
) const;
|
const word& AMIMethodName,
|
||||||
|
const bool interpAllPatches
|
||||||
|
);
|
||||||
|
|
||||||
|
//- Constructor helper
|
||||||
|
void constructFromCuttingPatches
|
||||||
|
(
|
||||||
|
const word& methodName,
|
||||||
|
const word& AMIMethodName,
|
||||||
|
const HashTable<word>& patchMap,
|
||||||
|
const wordList& cuttingPatches
|
||||||
|
);
|
||||||
|
|
||||||
//- Return the list of AMIs between source and target patches
|
//- Return the list of AMIs between source and target patches
|
||||||
const PtrList<AMIPatchToPatchInterpolation>& patchAMIs() const;
|
inline const PtrList<AMIPatchToPatchInterpolation>&
|
||||||
|
patchAMIs() const;
|
||||||
|
|
||||||
|
|
||||||
// Parallel operations
|
// Parallel operations
|
||||||
@ -237,6 +253,15 @@ public:
|
|||||||
const bool interpAllPatches = true
|
const bool interpAllPatches = true
|
||||||
);
|
);
|
||||||
|
|
||||||
|
//- Construct from source and target meshes, generic mapping methods
|
||||||
|
meshToMesh
|
||||||
|
(
|
||||||
|
const polyMesh& src,
|
||||||
|
const polyMesh& tgt,
|
||||||
|
const word& methodName, // internal mapping
|
||||||
|
const word& AMIMethodName, // boundary mapping
|
||||||
|
const bool interpAllPatches = true
|
||||||
|
);
|
||||||
|
|
||||||
//- Construct from source and target meshes
|
//- Construct from source and target meshes
|
||||||
meshToMesh
|
meshToMesh
|
||||||
@ -249,6 +274,18 @@ public:
|
|||||||
);
|
);
|
||||||
|
|
||||||
|
|
||||||
|
//- Construct from source and target meshes, generic mapping methods
|
||||||
|
meshToMesh
|
||||||
|
(
|
||||||
|
const polyMesh& src,
|
||||||
|
const polyMesh& tgt,
|
||||||
|
const word& methodName, // internal mapping
|
||||||
|
const word& AMIMethodName, // boundary mapping
|
||||||
|
const HashTable<word>& patchMap,
|
||||||
|
const wordList& cuttingPatches
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
//- Destructor
|
//- Destructor
|
||||||
virtual ~meshToMesh();
|
virtual ~meshToMesh();
|
||||||
|
|
||||||
@ -278,6 +315,13 @@ public:
|
|||||||
//- Return const access to the overlap volume
|
//- Return const access to the overlap volume
|
||||||
inline scalar V() const;
|
inline scalar V() const;
|
||||||
|
|
||||||
|
//- Conversion between mesh and patch interpolation methods
|
||||||
|
static AMIPatchToPatchInterpolation::interpolationMethod
|
||||||
|
interpolationMethodAMI
|
||||||
|
(
|
||||||
|
const interpolationMethod method
|
||||||
|
);
|
||||||
|
|
||||||
|
|
||||||
// Evaluation
|
// Evaluation
|
||||||
|
|
||||||
|
|||||||
@ -73,4 +73,11 @@ inline Foam::scalar Foam::meshToMesh::V() const
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
inline const Foam::PtrList<Foam::AMIPatchToPatchInterpolation>&
|
||||||
|
Foam::meshToMesh::patchAMIs() const
|
||||||
|
{
|
||||||
|
return patchAMIs_;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
// ************************************************************************* //
|
// ************************************************************************* //
|
||||||
|
|||||||
Reference in New Issue
Block a user