62 const bitSet& isMasterPoint,
63 const indirectPrimitivePatch&
pp,
64 const pointField& localPoints
71 Field<weightedPosition> wps
80 weightedPosition& wp = wps[verti];
81 for (
const label edgei : pointEdges[verti])
83 const label otherVerti = edges[edgei].otherVertex(verti);
85 if (isMasterPoint[meshPoints[otherVerti]])
88 wp.second() += localPoints[otherVerti];
95 tmp<pointField> tavg(
new pointField(wps.size()));
100 const weightedPosition& wp = wps[verti];
102 if (
mag(wp.first()) < VSMALL)
109 avg[verti] = wp.second()/wp.first();
117Foam::tmp<Foam::pointField>
118Foam::snappySnapDriver::smoothLambdaMuPatchDisplacement
129 const label
iters = 90;
130 const scalar
lambda = 0.33;
131 const scalar
mu = 0.34;
133 for (label iter = 0; iter <
iters; iter++)
137 (1 -
lambda)*newLocalPoints
142 (1 +
mu)*newLocalPoints
143 -
mu*avg(
mesh, isMasterPoint,
pp, newLocalPoints);
149Foam::tmp<Foam::scalarField>
150Foam::snappySnapDriver::wantedThickness
153 const scalar cellSizeFraction
158 const labelList& cellLevel = meshRefiner_.meshCutter().cellLevel();
162 const scalar edge0Len =
163 meshRefiner_.meshCutter().level0EdgeLength();
179 maxPointLevel[
f[fp]] =
max(maxPointLevel[
f[fp]], ownLevel);
194 const scalar edgeLen = edge0Len/(1<<maxPointLevel[pointi]);
195 thickness[pointi] = cellSizeFraction*edgeLen;
403Foam::autoPtr<Foam::displacementMotionSolver>
404Foam::snappySnapDriver::makeMotionSolver
425 tallDisp.ref().checkOut();
431 snapDict.get<
word>(
"solver"),
438 pMesh.thisDb().time().constant(),
452 pMesh.thisDb().time().constant(),
467void Foam::snappySnapDriver::setDisplacement
485 || (meshPoints.size() != patchDisp.size())
492 <<
" fld.size():" <<
fld.size()
493 <<
" patchDisp.size():" << patchDisp.size()
494 <<
" meshPoints.size():" << meshPoints.size()
495 <<
" mesh.nPoints():" << pMesh.size()
505 forAll(meshPoints, patchPointi)
507 const label meshPointi = meshPoints[patchPointi];
509 meshDisp[meshPointi] =
510 patchDisp[patchPointi]+
points[meshPointi]-
points0[meshPointi];
515 for (
const label patchi : adaptPatchIDs)
517 bfld[patchi] ==
pointField(meshDisp, bfld[patchi].
patch().meshPoints());
524 pcs.constrainDisplacement(
fld,
true);
528Foam::autoPtr<Foam::mapPolyMesh> Foam::snappySnapDriver::addBufferLayers
567 forAll(edgeGlobalFaces, edgei)
569 if (edgeGlobalFaces[edgei].size() > 2)
574 nPointLayers[
e[0]] = 0;
575 nPointLayers[
e[1]] = 0;
585 else if (edgeGlobalFaces[edgei].size() == 0)
608 if (!pLayers.found(label(1)))
610 nFaceLayers[facei] = 0;
645 addLayer.setRefinement
676 if (map.hasMotionPoints())
700 for (
const labelList& faceToCells : addedCells)
702 for (
const label celli : faceToCells)
707 meshRefiner_.updateMesh(map, isChangedFace.toc());
712 Info<<
"Writing mesh-with-layer to time "
713 << meshRefiner_.timeName() <<
endl;
722 meshRefiner_.timeName()
729void Foam::snappySnapDriver::doSnapBufferLayers
734 const scalar featureCos,
735 const scalar planarAngle,
743 <<
"Morphing phase" <<
nl
744 <<
"--------------" <<
nl
790 meshRefiner_.getZones
806 meshRefiner_.createZoneBaffles
816 meshRefiner_.dupNonManifoldBoundaryPoints();
822 meshRefiner_.subsetBaffles
831 bool doFeatures =
true;
833 if (snapParams.nFeatureSnap() > 0)
839 nFeatIter = snapParams.nFeatureSnap();
842 Info<<
"Snapping to features in " << nFeatIter
843 <<
" iterations ..." <<
endl;
847 const labelList adaptPatchIDs(meshRefiner_.meshedPatches());
866 /
"pp_initial_" + meshRefiner_.timeName() +
".obj"
870 ppPtr().localFaces(),
871 ppPtr().localPoints(),
885 smoothLambdaMuPatchDisplacement
892 pointField ppLocalPoints(ppPtr().localPoints()+patchDisp);
894 const bool smoothInternal =
true;
926 motionPtr().pointDisplacement()
937 patchDisp -= (ppLocalPoints-ppPtr().localPoints());
942 Info<<
"Writing smoothed mesh to time "
943 << meshRefiner_.timeName() <<
endl;
952 meshRefiner_.timeName()
959 Info<<
"Constructing mesh displacer ..." <<
endl;
960 Info<<
"Using mesh parameters " << motionDict <<
nl <<
endl;
979 motionPtr().pointDisplacement(),
987 Info<<
"Checking initial mesh ..." <<
endl;
996 Info<<
"Detected " << nInitErrors <<
" illegal faces"
997 <<
" (concave, zero area or negative cell pyramid volume)"
1001 Info<<
"Checked initial mesh in = "
1012 for (label iter = 0; iter < nFeatIter; iter++)
1015 <<
"Morph iteration " << iter <<
nl
1016 <<
"-----------------" <<
endl;
1026 const bool strictRegionSnap
1029 ? snapParams.strictRegionSnap()
1037 globalToMasterPatch_,
1038 globalToSlavePatch_,
1048 if (snapParams.detectNearSurfacesSnap())
1076 disp = calcNearestSurfaceFeature
1084 scalar(iter+1)/nFeatIter,
1112 patchFeaturePoint = ppLocalPoints+patchAttraction;
1119 /
"calcNearestSurfaceFeature"
1120 + meshRefiner_.timeName()
1123 forAll(ppLocalPoints, pointi)
1125 const point& pt = ppLocalPoints[pointi];
1141 meshRefiner_.doSplitFaces
1150 const Map<label> oldMeshPointMap(ppPtr->meshPointMap());
1151 Pout<<
"old pp points:" << ppPtr->nPoints()
1152 <<
" oldMeshPointMap:" << oldMeshPointMap.size()
1156 meshMoverPtr.clear();
1170 if (map.hasMotionPoints())
1189 splitFacesMap.insert(splitFaces[i], i);
1191 forAll(map.faceMap(), facei)
1193 if (splitFacesMap.find(map.faceMap()[facei]))
1195 changedFaces.append(facei);
1199 meshRefiner_.updateMesh
1202 meshRefiner_.growFaceCellFace(changedFaces)
1209 Info<<
"Writing split-faces mesh to time "
1210 << meshRefiner_.timeName() <<
endl;
1228 Info<<
"Updating for face-splitting" <<
endl;
1231 forAll(internalBaffles, i)
1234 baffle.first() = map.reverseFaceMap()[baffle.first()];
1235 baffle.second() = map.reverseFaceMap()[baffle.second()];
1240 motionPtr = makeMotionSolver
1254 motionPtr().pointDisplacement(),
1260 const auto&
mp = ppPtr->meshPoints();
1267 ppMap[i] = oldMeshPointMap[
mp[i]];
1296 Info<<
"DONE Updating for face-splitting" <<
endl;
1335 motionPtr().pointDisplacement()
1356 const vectorField newDisp(tnewPoints()-meshMoverPtr().oldPoints());
1358 meshMoverPtr().displacement().vectorField::operator=(newDisp);
1359 meshMoverPtr().setDisplacement(disp);
1362 const bool meshOk = scaleMesh
1373 <<
"Did not successfully snap mesh."
1374 <<
" Continuing to snap to resolve easy" <<
nl
1375 <<
" surfaces but the"
1376 <<
" resulting mesh will not satisfy your quality"
1377 <<
" constraints" <<
nl <<
endl;
1381 meshMoverPtr().correct();
1388 patchDisp -= (ppLocalPoints-ppPtr().localPoints());
1393 Info<<
"Writing partially moved mesh to time "
1394 << meshRefiner_.timeName() <<
endl;
1403 meshRefiner_.timeName()
1410 Info<<
nl <<
"Checking moved mesh ..." <<
endl;
1423 Info<<
"Detected " << nWrong
1425 <<
" (concave, zero area or negative cell pyramid volume)"
1432 for (
const label facei : wrongFaces)
1437 decomposeCell.set(own);
1444 decomposeCell.set(nei);
1448 Pout<<
"spliyyinG :" << decomposeCell.count() <<
" cells"
1453 tetDecomp.setRefinement
1461 const Map<label> oldMeshPointMap(ppPtr->meshPointMap());
1463 Pout<<
"old pp points:" << ppPtr->nPoints()
1464 <<
" oldMeshPointMap:" << oldMeshPointMap.size()
1468 meshMoverPtr.clear();
1482 if (map.hasMotionPoints())
1496 tetDecomp.updateMesh(map);
1499 forAll(map.cellMap(), celli)
1501 if (decomposeCell[map.cellMap()[celli]])
1513 Pout<<
"isChangedFace :" << decomposeCell.count() <<
" faces"
1518 meshRefiner_.updateMesh(map, isChangedFace.toc());
1524 Info<<
"Writing tet-decomp mesh to time "
1525 << meshRefiner_.timeName() <<
endl;
1542 Info<<
"Updating for tet-decomp" <<
endl;
1545 forAll(internalBaffles, i)
1548 baffle.first() = map.reverseFaceMap()[baffle.first()];
1549 baffle.second() = map.reverseFaceMap()[baffle.second()];
1554 Pout<<
"new pp points:" << ppPtr->nPoints()
1555 <<
" new meshPointMap:" << ppPtr->meshPointMap().size()
1557 const auto&
mp = ppPtr->meshPoints();
1564 const label oldMeshPointi = map.pointMap()[
mp[i]];
1565 const auto mpFnd = oldMeshPointMap.find(oldMeshPointi);
1598 snapDist = calcSnapDistance(
mesh, snapParams, ppPtr());
1601 Info<<
"DONE Updating for tet-decomp" <<
endl;
1834 const Map<label> oldMeshPointMap(ppPtr->meshPointMap());
1836 meshMoverPtr.clear();
1841 wantedThickness(ppPtr(), 1
e-3)
1859 forAll(addedPoints, oldPatchPointi)
1861 const label
newPointi = addedPoints[oldPatchPointi].last();
1862 const label newPatchPointi = ppPtr().meshPointMap()[
newPointi];
1863 ppMap[newPatchPointi] = oldPatchPointi;
1889 const bool snapToGeometry =
true;
1921 patchFeaturePoint-ppLocalPoints,
1924 motionPtr().pointDisplacement()
1933 /
"buffer_layer_new_projection"
1934 + meshRefiner_.timeName()
1937 forAll(ppLocalPoints, pointi)
1939 const point& pt = ppLocalPoints[pointi];
1940 str.write(
linePointRef(pt, patchFeaturePoint[pointi]));
1942 Pout<<
"** writing mapped attraction to " << str.name() <<
endl;
1955 patchDisp -= (ppLocalPoints-ppPtr().localPoints());
1960 Info<<
"Writing smoothed LAYER mesh to time "
1961 << meshRefiner_.timeName() <<
endl;
1970 meshRefiner_.timeName()
1989 Info<<
"Writing baffle-merged mesh to time "
1990 << meshRefiner_.timeName() <<
endl;
1999 meshRefiner_.timeName()
2007 labelList duplicateFace(getInternalOrBaffleDuplicateFace());
2009 repatchToSurface(snapParams, adaptPatchIDs, duplicateFace);
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
GeometricBoundaryField< vector, pointPatchField, pointMesh > Boundary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const Addr & addressing() const noexcept
The addressing used for the list.
label size() const noexcept
The number of elements in the list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A HashTable to objects of type <T> with a label key.
static FOAM_NO_DANGLING_REFERENCE const pointConstraints & New(const pointMesh &mesh, Args &&... args)
const Mesh & mesh() const noexcept
Reference to the mesh.
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
label nPoints() const
Number of points supporting patch faces.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
fileName path() const
The path for the case = rootPath/caseName.
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
label size() const noexcept
The number of entries in the list.
Adds layers of cells to outside (or inside) of polyMesh. Can optionally create stand-alone extruded m...
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp, const bitSet &orientation)
Per patch edge the pp faces (in global indices) using it.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
double cpuTimeIncrement() const
Return CPU time [seconds] since last call to cpuTimeIncrement(), resetCpuTimeIncrement().
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
static autoPtr< displacementMotionSolver > New(const word &solverTypeName, const polyMesh &, const IOdictionary &, const pointVectorField &pointDisplacement, const pointIOField &points0)
Select constructed from polyMesh, dictionary and components.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
A face is a list of labels corresponding to mesh vertices.
A class for handling file names.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
void clearOut(const bool isMeshUpdate=false)
Clear all geometry and addressing.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
static labelPairList findDuplicateFacePairs(const polyMesh &)
Helper routine to find all baffles (two boundary faces.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
writeType
Enumeration for what to write. Used as a bit-pattern.
debugType
Enumeration for what to debug. Used as a bit-pattern.
static tmp< pointVectorField > makeDisplacementField(const pointMesh &pMesh, const labelList &adaptPatchIDs)
Helper function to make a pointVectorField with correct.
FaceMergeType
Enumeration for what to do with co-planar patch faces on a single.
static autoPtr< indirectPrimitivePatch > makePatch(const polyMesh &, const labelList &)
Create patch from set of patches.
static void updateList(const labelList &newToOld, const T &nullValue, List< T > &elems)
Helper: reorder list according to map.
static writeType writeLevel()
Get/set write level.
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces, const bool dryRun=false)
Check mesh with mesh settings in dict. Collects incorrect faces.
Given a displacement moves the mesh by scaling the displacement back until there are no more mesh err...
static const weightedPosition zero
Accumulates point constraints through successive applications of the applyConstraint function.
Application of (multi-)patch point constraints.
Mesh representing a set of points created from polyMesh.
Mesh consisting of general polyhedral cells.
void setInstance(const fileName &instance, const IOobjectOption::writeOption wOpt=IOobject::AUTO_WRITE)
Set the instance for mesh files.
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
const fileName & facesInstance() const
Return the current instance directory for faces.
bool moving() const noexcept
Is mesh moving.
virtual const labelList & faceOwner() const
Return face owner.
virtual const labelList & faceNeighbour() const
Return face neighbour.
virtual const pointField & points() const
Return raw points.
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh").
Direct mesh changes based on v1.3 polyTopoChange syntax.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const cellList & cells() const
Simple container to keep together snap specific information.
static void determineSidePatches(meshRefinement &meshRefiner, const globalIndex &globalFaces, const labelListList &edgeGlobalFaces, const indirectPrimitivePatch &pp, labelList &edgePatchID, labelList &edgeZoneID, boolList &edgeFlip, labelList &inflateFaceID)
Helper: see what zones and patches edges should be extruded into.
Decomposes polyMesh into tets (or pyramids).
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for TET. (4 tri).
A class for managing temporary objects.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
static void syncPoints(const polyMesh &mesh, List< weightedPosition > &)
Synchronisation for mesh point positions.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
const std::string patch
OpenFOAM patch number as a std::string.
List< edge > edgeList
List of edge.
Pair< label > labelPair
A pair of labels.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
List< labelPair > labelPairList
List of labelPair.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
vectorIOField pointIOField
pointIOField is a vectorIOField.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
line< point, const point & > linePointRef
A line using referred points.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
List< bool > boolList
A List of bools.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
#define addProfiling(Name,...)
Define profiling trigger with specified name and description string. The description is generated by ...
loopControl iters(runTime, aMesh.solutionDict(), "solution")
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Unit conversion functions.
pointField points0(pointIOField(IOobject("points", mesh.time().constant(), polyMesh::meshSubDir, mesh, IOobject::MUST_READ, IOobject::NO_WRITE, IOobject::NO_REGISTER)))