Go to the documentation of this file.
60 const Foam::scalar Foam::perfectInterface::tol_ = 1
e-3;
76 ctrs[patchFacei] = pp[patchFacei].centre(
points);
85 Foam::perfectInterface::perfectInterface
90 const word& faceZoneName,
91 const word& masterPatchName,
92 const word& slavePatchName
96 faceZoneID_(faceZoneName, mme.
mesh().faceZones()),
102 Foam::perfectInterface::perfectInterface
114 mme.
mesh().faceZones()
138 Pout<<
"bool perfectInterface::changeTopology() const "
139 <<
"for object " <<
name() <<
" : "
140 <<
"Inactive" <<
endl;
160 const polyMesh&
mesh = topoChanger().mesh();
165 const edgeList& edges0 = pp0.edges();
168 const labelList& meshPts0 = pp0.meshPoints();
169 const labelList& meshPts1 = pp1.meshPoints();
174 scalar minLen = GREAT;
178 minLen =
min(minLen, edges0[edgeI].
mag(pts0));
180 scalar typDim = tol_*minLen;
184 Pout<<
"typDim:" << typDim <<
" edges0:" << edges0.size()
185 <<
" pts0:" << pts0.size() <<
" pts1:" << pts1.size()
186 <<
" pp0:" << pp0.size() <<
" pp1:" << pp1.size() <<
endl;
196 renumberPoints[i] = i;
213 <<
"Points on patch sides do not match to within tolerance "
219 renumberPoints[meshPts1[i]] = meshPts0[from1To0Points[i]];
231 calcFaceCentres(pp0),
232 calcFaceCentres(pp1),
241 <<
"Face centres of patch sides do not match to within tolerance "
256 label meshPointi = meshPts1[i];
258 if (meshPointi != renumberPoints[meshPointi])
262 affectedFaces.insert(
pFaces);
267 affectedFaces.erase(pp1.addressing()[i]);
275 label facei = pp0.addressing()[i];
277 if (affectedFaces.erase(facei))
280 <<
"Found face " << facei <<
" vertices "
281 <<
mesh.
faces()[facei] <<
" whose points are"
282 <<
" used both by master patch and slave patch" <<
endl;
288 for (
const label facei : affectedFaces)
292 face newFace(
f.size());
296 newFace[fp] = renumberPoints[
f[fp]];
314 bool zoneFlip =
false;
320 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
344 label meshPointi = meshPts1[i];
346 if (meshPointi != renumberPoints[meshPointi])
348 ref.setAction(polyRemovePoint(meshPointi));
356 label facei = pp1.addressing()[i];
357 ref.setAction(polyRemoveFace(facei));
369 label facei = pp0.addressing()[i];
373 face newFace(
f.size());
377 newFace[fp] = renumberPoints[
f[fp]];
382 label pp1Facei = pp1.addressing()[from0To1Faces[i]];
410 newFace.reverseFace(),
430 Pout<<
"bool perfectInterface::setRefinement(polyTopoChange&) const : "
431 <<
"for object " <<
name() <<
" : "
432 <<
"masterPatchID_:" << masterPatchID_
433 <<
" slavePatchID_:" << slavePatchID_
434 <<
" faceZoneID_:" << faceZoneID_ <<
endl;
439 masterPatchID_.active()
440 && slavePatchID_.active()
441 && faceZoneID_.active()
444 const polyMesh&
mesh = topoChanger().mesh();
447 const polyPatch& patch0 =
patches[masterPatchID_.index()];
448 const polyPatch& patch1 =
patches[slavePatchID_.index()];
462 setRefinement(pp0, pp1,
ref);
489 << faceZoneID_.name() <<
nl
490 << masterPatchID_.name() <<
nl
491 << slavePatchID_.name() <<
endl;
503 os.
writeEntry(
"masterPatchName", masterPatchID_.name());
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
List< label > labelList
A List of labels.
vectorField pointField
pointField is a vectorField.
virtual const pointField & points() const
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< edge > edgeList
A List of edges.
A class for handling words, derived from Foam::string.
const boolList & flipMap() const noexcept
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const edgeList & edges() const
Class containing data for point removal.
List of mesh modifiers defining the mesh dynamics.
virtual void modifyMotionPoints(pointField &motionPoints) const
const labelListList & pointFaces() const
Direct mesh changes based on v1.3 polyTopoChange syntax.
List< bool > boolList
A List of bools.
virtual void updateMesh(const mapPolyMesh &)
Class describing modification of a face.
const polyBoundaryMesh & boundaryMesh() const
Ostream & endl(Ostream &os)
virtual Ostream & beginBlock(const keyType &kw)
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
virtual void setRefinement(polyTopoChange &) const
A HashTable with keys but without contents that is similar to std::unordered_set.
label min(const labelHashSet &set, label minValue=labelMax)
virtual void writeDict(Ostream &) const
Mesh consisting of general polyhedral cells.
A class for handling keywords in dictionaries.
Determine correspondence between points. See below.
Generic templated field type.
A subset of mesh faces organised as a primitive patch.
A patch is a list of labels that address the faces in the global face list.
virtual const labelList & faceOwner() const
A List with indirect addressing.
virtual bool changeTopology() const
const faceZoneMesh & faceZones() const noexcept
label whichPatch(const label faceIndex) const
label whichZone(const label objectIndex) const
virtual Ostream & endBlock()
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
OBJstream os(runTime.globalPath()/outputName)
Macros for easy insertion into run-time selection tables.
const labelIOList & zoneID
label whichFace(const label globalCellID) const
errorManipArg< error, int > exit(error &err, const int errNo=1)
const Field< point_type > & localPoints() const
Virtual base class for mesh modifiers.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
bool isInternalFace(const label faceIndex) const noexcept
virtual const faceList & faces() const
#define FatalErrorInFunction
const polyMesh & mesh() const
Class containing data for face removal.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
const dimensionedScalar e
labelList identity(const label len, label start=0)
bool insert(const Key &key)
const polyBoundaryMesh & patches
Ostream & writeEntry(const keyType &key, const T &value)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
word name(const expressions::valueTypeCode typeCode)
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A face is a list of labels corresponding to mesh vertices.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
const labelList & meshPoints() const
defineTypeNameAndDebug(combustionModel, 0)
bool matchPoints(const UList< point > &pts0, const UList< point > &pts1, const UList< scalar > &matchDistance, const bool verbose, labelList &from0To1, const point &origin=point::zero)
#define WarningInFunction
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
virtual const labelList & faceNeighbour() const
A list of faces which address into the list of points.
virtual void write(Ostream &) const