ENH: patchDist: helper class to determine distance on patch to neighbouring patches

This commit is contained in:
mattijs 2011-12-16 11:28:44 +00:00
parent df5feff482
commit af25e71e1c
4 changed files with 385 additions and 61 deletions

View File

@ -32,6 +32,7 @@ Description
#include "volFields.H"
#include "PatchEdgeFaceWave.H"
#include "patchEdgeFaceInfo.H"
#include "patchDist.H"
using namespace Foam;
@ -49,9 +50,10 @@ int main(int argc, char *argv[])
// Get name of patch
const word patchName = args[1];
const polyPatch& patch = patches[patchName];
// 1. Walk from a single edge
{
// Data on all edges and faces
List<patchEdgeFaceInfo> allEdgeInfo(patch.nEdges());
List<patchEdgeFaceInfo> allFaceInfo(patch.size());
@ -122,8 +124,42 @@ int main(int argc, char *argv[])
<< endl;
vsf.write();
}
// 2. Use a wrapper to walk from all boundary edges on selected patches
{
labelHashSet otherPatchIDs(identity(mesh.boundaryMesh().size()));
otherPatchIDs.erase(patch.index());
Info<< "Walking on patch " << patch.index()
<< " from edges shared with patches " << otherPatchIDs
<< endl;
patchDist pwd(patch, otherPatchIDs);
// Extract as patchField
volScalarField vsf
(
IOobject
(
"otherPatchDist",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("otherPatchDist", dimLength, 0.0)
);
vsf.boundaryField()[patch.index()] = pwd;
Info<< "Writing otherPatchDist volScalarField to " << runTime.value()
<< endl;
vsf.write();
}
Info<< "\nEnd\n" << endl;
return 0;
}

View File

@ -41,6 +41,7 @@ $(pWave)/pointEdgePoint.C
patchWave = $(algorithms)/PatchEdgeFaceWave
$(patchWave)/PatchEdgeFaceWaveName.C
$(patchWave)/patchEdgeFaceInfo.C
$(patchWave)/patchDist.C
meshWave = $(algorithms)/MeshWave
$(meshWave)/MeshWaveName.C

View File

@ -0,0 +1,173 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "patchDist.H"
#include "PatchEdgeFaceWave.H"
#include "syncTools.H"
#include "polyMesh.H"
#include "patchEdgeFaceInfo.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::patchDist::patchDist
(
const polyPatch& patch,
const labelHashSet& nbrPatchIDs
)
:
patch_(patch),
nbrPatchIDs_(nbrPatchIDs),
nUnset_(0)
{
patchDist::correct();
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::patchDist::~patchDist()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::patchDist::correct()
{
// Mark all edge connected to a nbrPatch.
label nBnd = 0;
forAllConstIter(labelHashSet, nbrPatchIDs_, iter)
{
label nbrPatchI = iter.key();
const polyPatch& nbrPatch = patch_.boundaryMesh()[nbrPatchI];
nBnd += nbrPatch.nEdges()-nbrPatch.nInternalEdges();
}
// Mark all edges. Note: should use HashSet but have no syncTools
// functionality for these.
EdgeMap<label> nbrEdges(2*nBnd);
forAllConstIter(labelHashSet, nbrPatchIDs_, iter)
{
label nbrPatchI = iter.key();
const polyPatch& nbrPatch = patch_.boundaryMesh()[nbrPatchI];
const labelList& nbrMp = nbrPatch.meshPoints();
for
(
label edgeI = nbrPatch.nInternalEdges();
edgeI < nbrPatch.nEdges();
edgeI++
)
{
const edge& e = nbrPatch.edges()[edgeI];
const edge meshE = edge(nbrMp[e[0]], nbrMp[e[1]]);
nbrEdges.insert(meshE, nbrPatchI);
}
}
// Make sure these boundary edges are marked everywhere.
syncTools::syncEdgeMap
(
patch_.boundaryMesh().mesh(),
nbrEdges,
maxEqOp<label>()
);
// Data on all edges and faces
List<patchEdgeFaceInfo> allEdgeInfo(patch_.nEdges());
List<patchEdgeFaceInfo> allFaceInfo(patch_.size());
// Initial seed
label nBndEdges = patch_.nEdges() - patch_.nInternalEdges();
DynamicList<label> initialEdges(2*nBndEdges);
DynamicList<patchEdgeFaceInfo> initialEdgesInfo(2*nBndEdges);
// Seed all my edges that are also nbrEdges
const labelList& mp = patch_.meshPoints();
for
(
label edgeI = patch_.nInternalEdges();
edgeI < patch_.nEdges();
edgeI++
)
{
const edge& e = patch_.edges()[edgeI];
const edge meshE = edge(mp[e[0]], mp[e[1]]);
EdgeMap<label>::const_iterator edgeFnd = nbrEdges.find(meshE);
if (edgeFnd != nbrEdges.end())
{
initialEdges.append(edgeI);
initialEdgesInfo.append
(
patchEdgeFaceInfo
(
e.centre(patch_.localPoints()),
0.0
)
);
}
}
// Walk
PatchEdgeFaceWave
<
primitivePatch,
patchEdgeFaceInfo
> calc
(
patch_.boundaryMesh().mesh(),
patch_,
initialEdges,
initialEdgesInfo,
allEdgeInfo,
allFaceInfo,
returnReduce(patch_.nEdges(), sumOp<label>())
);
// Extract into *this
setSize(patch_.size());
nUnset_ = 0;
forAll(allFaceInfo, faceI)
{
if (allFaceInfo[faceI].valid(calc.data()))
{
operator[](faceI) = Foam::sqrt(allFaceInfo[faceI].distSqr());
}
else
{
nUnset_++;
}
}
}
// ************************************************************************* //

View File

@ -0,0 +1,114 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::patchDist
Description
Like wallDist but calculates on patch distance to nearest neighbouring
patches. Uses PatchEdgeFaceWave to do actual calculation.
SourceFiles
patchDist.C
\*---------------------------------------------------------------------------*/
#ifndef patchDist_H
#define patchDist_H
#include "scalarField.H"
#include "HashSet.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
//class polyMesh;
class polyPatch;
/*---------------------------------------------------------------------------*\
Class patchDist Declaration
\*---------------------------------------------------------------------------*/
class patchDist
:
public scalarField
{
private:
// Private Member Data
//- Patch to operate on
const polyPatch& patch_;
//- Patches to determine the distance to
const labelHashSet nbrPatchIDs_;
//- Number of unset faces.
label nUnset_;
public:
// Constructors
//- Construct from patch and neighbour patches.
patchDist
(
const polyPatch& pp,
const labelHashSet& nbrPatchIDs
);
//- Destructor
virtual ~patchDist();
// Member Functions
const scalarField& y() const
{
return *this;
}
label nUnset() const
{
return nUnset_;
}
//- Correct for mesh geom/topo changes
virtual void correct();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //