mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
Initial autoMesh merge
This commit is contained in:
667
src/autoMesh/autoHexMesh/meshRefinement/meshRefinement.H
Normal file
667
src/autoMesh/autoHexMesh/meshRefinement/meshRefinement.H
Normal file
@ -0,0 +1,667 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 1991-2007 OpenCFD Ltd.
|
||||
\\/ 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 2 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, write to the Free Software Foundation,
|
||||
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
||||
|
||||
Class
|
||||
meshRefinement
|
||||
|
||||
Description
|
||||
Helper class which maintains intersections of (changing) mesh with
|
||||
(static) surfaces.
|
||||
Maintains
|
||||
- per face any intersections of this edge with any of the surfaces
|
||||
|
||||
SourceFiles
|
||||
meshRefinement.C
|
||||
meshRefinementBaffles.C
|
||||
meshRefinementMerge.C
|
||||
meshRefinementRefine.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef meshRefinement_H
|
||||
#define meshRefinement_H
|
||||
|
||||
#include "hexRef8.H"
|
||||
#include "mapPolyMesh.H"
|
||||
#include "autoPtr.H"
|
||||
#include "pointIOField.H"
|
||||
#include "labelPair.H"
|
||||
#include "indirectPrimitivePatch.H"
|
||||
#include "pointFieldsFwd.H"
|
||||
#include "Tuple2.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
// Class forward declarations
|
||||
class fvMesh;
|
||||
class mapDistributePolyMesh;
|
||||
class decompositionMethod;
|
||||
class refinementSurfaces;
|
||||
class removeCells;
|
||||
class featureEdgeMesh;
|
||||
class fvMeshDistribute;
|
||||
class searchableSurface;
|
||||
class regionSplit;
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class meshRefinement Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class meshRefinement
|
||||
{
|
||||
|
||||
public:
|
||||
|
||||
// Public data types
|
||||
|
||||
//- Enumeration for debug dumping
|
||||
enum writeFlag
|
||||
{
|
||||
MESH = 1,
|
||||
SCALARLEVELS = 2,
|
||||
OBJINTERSECTIONS = 4
|
||||
};
|
||||
|
||||
|
||||
//- Enumeration for how the userdata is to be mapped upon refinement.
|
||||
enum mapType
|
||||
{
|
||||
MASTERONLY = 1, // maintain master only
|
||||
KEEPALL = 2, // have slaves (upon refinement) from master
|
||||
REMOVE = 4 // set value to -1 any face that has been refined
|
||||
};
|
||||
|
||||
|
||||
private:
|
||||
|
||||
// Private data
|
||||
|
||||
//- Reference to mesh
|
||||
fvMesh& mesh_;
|
||||
|
||||
//- tolerance used for sorting coordinates (used in 'less' routine)
|
||||
const scalar tol_;
|
||||
|
||||
const refinementSurfaces& surfaces_;
|
||||
|
||||
//- refinement engine
|
||||
hexRef8 meshCutter_;
|
||||
|
||||
//- per cc-cc vector the index of the surface hit
|
||||
labelIOList surfaceIndex_;
|
||||
|
||||
//- user supplied face based data.
|
||||
List<Tuple2<mapType, labelList> > userFaceData_;
|
||||
|
||||
|
||||
// Private Member Functions
|
||||
|
||||
//- Reorder list according to map.
|
||||
template<class T>
|
||||
static void updateList
|
||||
(
|
||||
const labelList& newToOld,
|
||||
const T& nullValue,
|
||||
List<T>& elems
|
||||
);
|
||||
|
||||
//- Add patchfield of given type to all fields on mesh
|
||||
template<class GeoField>
|
||||
static void addPatchFields(fvMesh&, const word& patchFieldType);
|
||||
|
||||
//- Reorder patchfields of all fields on mesh
|
||||
template<class GeoField>
|
||||
static void reorderPatchFields(fvMesh&, const labelList& oldToNew);
|
||||
|
||||
//- Find out which faces have changed given cells (old mesh labels)
|
||||
// that were marked for refinement.
|
||||
static labelList getChangedFaces
|
||||
(
|
||||
const mapPolyMesh&,
|
||||
const labelList& oldCellsToRefine
|
||||
);
|
||||
|
||||
//- Calculate coupled boundary end vector and refinement level
|
||||
void calcNeighbourData
|
||||
(
|
||||
labelList& neiLevel,
|
||||
pointField& neiCc
|
||||
) const;
|
||||
////- Calculate coupled boundary end vector and refinement level. Sort
|
||||
//// so both sides of coupled face have data in same order.
|
||||
//void calcCanonicalBoundaryData
|
||||
//(
|
||||
// labelList& leftLevel,
|
||||
// pointField& leftCc,
|
||||
// labelList& rightLevel,
|
||||
// pointField& rightCc
|
||||
//) const;
|
||||
|
||||
//- Find any intersection of surface. Store in surfaceIndex_.
|
||||
void updateIntersections(const labelList& changedFaces);
|
||||
|
||||
//- Set instance of all local IOobjects
|
||||
void setInstance(const fileName&);
|
||||
|
||||
//- Remove cells. Put exposedFaces into exposedPatchIDs.
|
||||
autoPtr<mapPolyMesh> doRemoveCells
|
||||
(
|
||||
const labelList& cellsToRemove,
|
||||
const labelList& exposedFaces,
|
||||
const labelList& exposedPatchIDs,
|
||||
removeCells& cellRemover
|
||||
);
|
||||
|
||||
// Get cells which are inside any closed surface. Note that
|
||||
// all closed surfaces
|
||||
// will have already been oriented to have keepPoint outside.
|
||||
labelList getInsideCells(const word&) const;
|
||||
|
||||
// Do all to remove inside cells
|
||||
autoPtr<mapPolyMesh> removeInsideCells
|
||||
(
|
||||
const string& msg,
|
||||
const label exposedPatchI
|
||||
);
|
||||
|
||||
//- Used in decomposeCombineRegions. Given global region per cell
|
||||
// determines master processor/cell for regions straddling
|
||||
// procboundaries.
|
||||
void getRegionMaster
|
||||
(
|
||||
const boolList& blockedFace,
|
||||
const regionSplit& globalRegion,
|
||||
Map<label>& regionToMaster
|
||||
) const;
|
||||
|
||||
|
||||
// Refinement candidate selection
|
||||
|
||||
//- Mark cell for refinement (if not already marked). Return false
|
||||
// if refinelimit hit. Keeps running count (in nRefine) of cells
|
||||
// marked for refinement
|
||||
static bool markForRefine
|
||||
(
|
||||
const label markValue,
|
||||
const label nAllowRefine,
|
||||
label& cellValue,
|
||||
label& nRefine
|
||||
);
|
||||
|
||||
//- Calculate list of cells to refine based on intersection of
|
||||
// features.
|
||||
label markFeatureRefinement
|
||||
(
|
||||
const point& keepPoint,
|
||||
const PtrList<featureEdgeMesh>& featureMeshes,
|
||||
const labelList& featureLevels,
|
||||
const label nAllowRefine,
|
||||
|
||||
labelList& refineCell,
|
||||
label& nRefine
|
||||
) const;
|
||||
|
||||
//- Mark cells for refinement-shells based refinement.
|
||||
label markInternalRefinement
|
||||
(
|
||||
const PtrList<searchableSurface>&,
|
||||
const labelList& shellLevels,
|
||||
const boolList& shellRefineInside,
|
||||
const label nAllowRefine,
|
||||
labelList& refineCell,
|
||||
label& nRefine
|
||||
) const;
|
||||
|
||||
//- Mark cells for surface intersection based refinement.
|
||||
label markSurfaceRefinement
|
||||
(
|
||||
const label nAllowRefine,
|
||||
const labelList& neiLevel,
|
||||
const pointField& neiCc,
|
||||
labelList& refineCell,
|
||||
label& nRefine
|
||||
) const;
|
||||
|
||||
//- Mark cell if local curvature > curvature or
|
||||
// markDifferingRegions = true and intersections with different
|
||||
// regions.
|
||||
bool checkCurvature
|
||||
(
|
||||
const labelList& globalToPatch,
|
||||
const scalar curvature,
|
||||
const bool markDifferingRegions,
|
||||
const label nAllowRefine,
|
||||
const label newSurfI,
|
||||
const label newTriI,
|
||||
const label cellI,
|
||||
|
||||
label& maxCellSurfI,
|
||||
label& maxCellTriI,
|
||||
|
||||
labelList& refineCell,
|
||||
label& nRefine
|
||||
) const;
|
||||
|
||||
|
||||
//- Mark cells for surface curvature based refinement. Marks if
|
||||
// local curvature > curvature or if on different regions
|
||||
// (markDifferingRegions)
|
||||
label markSurfaceCurvatureRefinement
|
||||
(
|
||||
const labelList& globalToPatch,
|
||||
const scalar curvature,
|
||||
const bool markDifferingRegions,
|
||||
const label nAllowRefine,
|
||||
const labelList& neiLevel,
|
||||
const pointField& neiCc,
|
||||
labelList& refineCell,
|
||||
label& nRefine
|
||||
) const;
|
||||
|
||||
// Baffle handling
|
||||
|
||||
//- Determine patches for baffles
|
||||
void getBafflePatches
|
||||
(
|
||||
const labelList& globalToPatch,
|
||||
const labelList& neiLevel,
|
||||
const pointField& neiCc,
|
||||
labelList& ownPatch,
|
||||
labelList& neiPatch
|
||||
) const;
|
||||
|
||||
//- Determine patch for baffle using some heuristic (and not
|
||||
// surface)
|
||||
label getBafflePatch
|
||||
(
|
||||
const labelList& facePatch,
|
||||
const label faceI
|
||||
) const;
|
||||
|
||||
//- Repatches external face or creates baffle for internal face
|
||||
// with user specified patches (might be different for both sides).
|
||||
// Returns label of added face.
|
||||
label createBaffle
|
||||
(
|
||||
const label faceI,
|
||||
const label ownPatch,
|
||||
const label neiPatch,
|
||||
polyTopoChange& meshMod
|
||||
) const;
|
||||
|
||||
//- Helper function to mark face as being on 'boundary'. Used by
|
||||
// markFacesOnProblemCells
|
||||
void markBoundaryFace
|
||||
(
|
||||
const label faceI,
|
||||
boolList& isBoundaryFace,
|
||||
boolList& isBoundaryEdge,
|
||||
boolList& isBoundaryPoint
|
||||
) const;
|
||||
|
||||
//- Returns list with for every internal face -1 or the patch
|
||||
// they should be baffled into.
|
||||
labelList markFacesOnProblemCells
|
||||
(
|
||||
const labelList& globalToPatch
|
||||
) const;
|
||||
|
||||
|
||||
// Baffle merging
|
||||
|
||||
//- Extract those baffles (duplicate) faces that are on the edge
|
||||
// of a baffle region. These are candidates for merging.
|
||||
List<labelPair> filterDuplicateFaces(const List<labelPair>&) const;
|
||||
|
||||
|
||||
// Zone handling
|
||||
|
||||
//- Finds zone per cell for cells inside closed named surfaces.
|
||||
// (uses geometric test for insideness)
|
||||
void findCellZoneGeometric
|
||||
(
|
||||
const labelList& closedNamedSurfaces,
|
||||
const labelList& namedSurfaceIndex,
|
||||
const labelList& surfaceToCellZone,
|
||||
labelList& cellToZone
|
||||
) const;
|
||||
|
||||
//- Finds zone per cell. Uses topological walk with all faces
|
||||
// marked in namedSurfaceIndex regarded as blocked.
|
||||
void findCellZoneTopo
|
||||
(
|
||||
const point& keepPoint,
|
||||
const labelList& namedSurfaceIndex,
|
||||
const labelList& surfaceToCellZone,
|
||||
labelList& cellToZone
|
||||
) const;
|
||||
|
||||
|
||||
//- Disallow default bitwise copy construct
|
||||
meshRefinement(const meshRefinement&);
|
||||
|
||||
//- Disallow default bitwise assignment
|
||||
void operator=(const meshRefinement&);
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
ClassName("meshRefinement");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
meshRefinement
|
||||
(
|
||||
fvMesh& mesh,
|
||||
const scalar tol,
|
||||
const refinementSurfaces&
|
||||
);
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
// Access
|
||||
|
||||
//- reference to mesh
|
||||
const fvMesh& mesh() const
|
||||
{
|
||||
return mesh_;
|
||||
}
|
||||
fvMesh& mesh()
|
||||
{
|
||||
return mesh_;
|
||||
}
|
||||
|
||||
//- reference to surface search engines
|
||||
const refinementSurfaces& surfaces() const
|
||||
{
|
||||
return surfaces_;
|
||||
}
|
||||
|
||||
//- reference to meshcutting engine
|
||||
const hexRef8& meshCutter() const
|
||||
{
|
||||
return meshCutter_;
|
||||
}
|
||||
|
||||
//- per start-end edge the index of the surface hit
|
||||
const labelList& surfaceIndex() const
|
||||
{
|
||||
return surfaceIndex_;
|
||||
}
|
||||
|
||||
labelList& surfaceIndex()
|
||||
{
|
||||
return surfaceIndex_;
|
||||
}
|
||||
|
||||
//- Additional face data that is maintained across
|
||||
// topo changes. Every entry is a list over all faces.
|
||||
// Bit of a hack. Additional flag to say whether to maintain master
|
||||
// only (false) or increase set to account for face-from-face.
|
||||
const List<Tuple2<mapType, labelList> >& userFaceData() const
|
||||
{
|
||||
return userFaceData_;
|
||||
}
|
||||
|
||||
List<Tuple2<mapType, labelList> >& userFaceData()
|
||||
{
|
||||
return userFaceData_;
|
||||
}
|
||||
|
||||
|
||||
// Other
|
||||
|
||||
//- Count number of intersections (local)
|
||||
label countHits() const;
|
||||
|
||||
//- Helper function to get decomposition such that all connected
|
||||
// regions get moved onto one processor. Used to prevent baffles
|
||||
// straddling processor boundaries. explicitConnections is to
|
||||
// keep pairs of non-coupled boundary faces together
|
||||
// (e.g. to keep baffles together)
|
||||
labelList decomposeCombineRegions
|
||||
(
|
||||
const boolList& blockedFace,
|
||||
const List<labelPair>& explicitConnections,
|
||||
decompositionMethod&
|
||||
) const;
|
||||
|
||||
//- Get faces with intersection.
|
||||
labelList intersectedFaces() const;
|
||||
|
||||
//- Get points on surfaces with intersection and boundary faces.
|
||||
labelList intersectedPoints() const;
|
||||
|
||||
//- Get added patches (inverse of globalToPatch)
|
||||
static labelList addedPatches(const labelList& globalToPatch);
|
||||
|
||||
//- Create patch from set of patches
|
||||
static autoPtr<indirectPrimitivePatch> makePatch
|
||||
(
|
||||
const polyMesh&,
|
||||
const labelList&
|
||||
);
|
||||
|
||||
//- Helper function to make a pointVectorField with correct
|
||||
// bcs for mesh movement:
|
||||
// - adaptPatchIDs : fixedValue
|
||||
// - processor : calculated (so free to move)
|
||||
// - cyclic/wedge/symmetry : slip
|
||||
// - other : slip
|
||||
static tmp<pointVectorField> makeDisplacementField
|
||||
(
|
||||
const pointMesh& pMesh,
|
||||
const labelList& adaptPatchIDs
|
||||
);
|
||||
|
||||
// Refinement
|
||||
|
||||
//- Calculate list of cells to refine.
|
||||
labelList refineCandidates
|
||||
(
|
||||
const point& keepPoint,
|
||||
const labelList& globalToPatch,
|
||||
const scalar curvature,
|
||||
|
||||
const PtrList<featureEdgeMesh>& featureMeshes,
|
||||
const labelList& featureLevels,
|
||||
|
||||
const PtrList<searchableSurface>&,
|
||||
const labelList& shellLevels,
|
||||
const boolList& shellRefineInside,
|
||||
|
||||
const bool featureRefinement,
|
||||
const bool internalRefinement,
|
||||
const bool surfaceRefinement,
|
||||
const bool curvatureRefinement,
|
||||
const label maxGlobalCells,
|
||||
const label maxLocalCells
|
||||
) const;
|
||||
|
||||
//- Refine some cells
|
||||
autoPtr<mapPolyMesh> refine(const labelList& cellsToRefine);
|
||||
|
||||
//- Refine some cells and rebalance
|
||||
autoPtr<mapDistributePolyMesh> refineAndBalance
|
||||
(
|
||||
const string& msg,
|
||||
decompositionMethod& decomposer,
|
||||
fvMeshDistribute& distributor,
|
||||
const labelList& cellsToRefine
|
||||
);
|
||||
|
||||
|
||||
// Baffle handling
|
||||
|
||||
//- Split off unreachable areas of mesh.
|
||||
void baffleAndSplitMesh
|
||||
(
|
||||
const bool handleSnapProblems,
|
||||
const bool mergeFreeStanding,
|
||||
Time& runTime,
|
||||
const labelList& globalToPatch,
|
||||
const point& keepPoint
|
||||
);
|
||||
|
||||
//- Split off (with optional buffer layers) unreachable areas
|
||||
// of mesh. Does not introduce baffles.
|
||||
autoPtr<mapPolyMesh> splitMesh
|
||||
(
|
||||
const label nBufferLayers,
|
||||
const labelList& globalToPatch,
|
||||
const point& keepPoint
|
||||
);
|
||||
|
||||
//- Find boundary points that connect to more than one cell
|
||||
// region and split them.
|
||||
autoPtr<mapPolyMesh> dupNonManifoldPoints();
|
||||
|
||||
//- Create baffle for every internal face where ownPatch != -1
|
||||
autoPtr<mapPolyMesh> createBaffles
|
||||
(
|
||||
const labelList& ownPatch,
|
||||
const labelList& neiPatch
|
||||
);
|
||||
|
||||
//- Return a list of coupled face pairs, i.e. faces that
|
||||
// use the same vertices.
|
||||
List<labelPair> getDuplicateFaces(const labelList& testFaces) const;
|
||||
|
||||
//- Merge baffles. Gets pairs of faces.
|
||||
autoPtr<mapPolyMesh> mergeBaffles(const List<labelPair>&);
|
||||
|
||||
//- Put faces/cells into zones according to surface specification.
|
||||
// Returns null if no zone surfaces present. Region containing
|
||||
// the keepPoint will not be put into a cellZone.
|
||||
autoPtr<mapPolyMesh> zonify(const point& keepPoint);
|
||||
|
||||
|
||||
// Other topo changes
|
||||
|
||||
//- Helper function to add patch to mesh
|
||||
static label addPatch(fvMesh&, const word& name, const word& type);
|
||||
|
||||
//- Split mesh. Keep part containing point.
|
||||
autoPtr<mapPolyMesh> splitMeshRegions(const point& keepPoint);
|
||||
|
||||
//- Update local numbering for mesh redistribution
|
||||
void distribute(const mapDistributePolyMesh&);
|
||||
|
||||
//- Update for external change to mesh. changedFaces are in new mesh
|
||||
// face labels.
|
||||
void updateMesh(const mapPolyMesh&, const labelList& changedFaces);
|
||||
|
||||
|
||||
// Restoring : is where other processes delete and reinsert data.
|
||||
|
||||
//- Signal points/face/cells for which to store data
|
||||
void storeData
|
||||
(
|
||||
const labelList& pointsToStore,
|
||||
const labelList& facesToStore,
|
||||
const labelList& cellsToStore
|
||||
);
|
||||
|
||||
//- Update local numbering + undo
|
||||
// Data to restore given as new pointlabel + stored pointlabel
|
||||
// (i.e. what was in pointsToStore)
|
||||
void updateMesh
|
||||
(
|
||||
const mapPolyMesh&,
|
||||
const labelList& changedFaces,
|
||||
const Map<label>& pointsToRestore,
|
||||
const Map<label>& facesToRestore,
|
||||
const Map<label>& cellsToRestore
|
||||
);
|
||||
|
||||
|
||||
//- Merge faces on the same patch (usually from exposing refinement)
|
||||
// Returns global number of faces merged.
|
||||
label mergePatchFaces
|
||||
(
|
||||
const scalar minCos,
|
||||
const scalar concaveCos
|
||||
//const dictionary& motionDict,
|
||||
//const labelList& globalToPatch
|
||||
);
|
||||
|
||||
//- Remove points not used by any face or points used
|
||||
// by only two faces where the edges are in line
|
||||
autoPtr<mapPolyMesh> mergeEdges(const scalar minCos);
|
||||
|
||||
|
||||
// Debug/IO
|
||||
|
||||
//- Debugging: check that all faces still obey start()>end()
|
||||
void checkData();
|
||||
|
||||
//- Compare two lists over all boundary faces
|
||||
template<class T>
|
||||
void testSyncBoundaryFaceList
|
||||
(
|
||||
const scalar tol,
|
||||
const string&,
|
||||
const UList<T>&,
|
||||
const UList<T>&
|
||||
) const;
|
||||
|
||||
//- Print some mesh stats.
|
||||
void printMeshInfo(const bool, const string&) const;
|
||||
|
||||
//- Write mesh and all data
|
||||
bool write() const;
|
||||
|
||||
//- Write refinement level as volScalarFields for postprocessing
|
||||
void dumpRefinementLevel() const;
|
||||
|
||||
//- Debug: Write intersection information to OBJ format
|
||||
void dumpIntersections(const fileName& prefix) const;
|
||||
|
||||
//- Do any one of above IO functions. flag is combination of
|
||||
// writeFlag values.
|
||||
void write(const label flag, const fileName&) const;
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#ifdef NoRepository
|
||||
# include "meshRefinementTemplates.C"
|
||||
#endif
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
||||
Reference in New Issue
Block a user