Go to the documentation of this file.
54 void Foam::sampledCuttingPlane::checkBoundsIntersection
66 if (!clipBb.overlaps(
meshBb))
71 <<
"Bounds " << clipBb
72 <<
" do not overlap the mesh bounding box " <<
meshBb
77 if (!clipBb.intersects(pln))
82 <<
"Plane "<< pln <<
" does not intersect the bounds "
89 if (!
meshBb.intersects(pln))
94 <<
"Plane "<< pln <<
" does not intersect the mesh bounds "
101 void Foam::sampledCuttingPlane::setDistanceFields(
const plane& pln)
108 const fvMesh&
mesh = cellDistance.mesh();
120 fld[i] = pln.signedDistance(cc[i]);
126 volScalarField::Boundary& cellDistanceBf =
127 cellDistance.boundaryFieldRef();
129 forAll(cellDistanceBf, patchi)
133 isA<emptyFvPatchScalarField>
135 cellDistanceBf[patchi]
142 new calculatedFvPatchScalarField
153 fld.setSize(pp.size());
156 fld[i] = pln.signedDistance(cc[i]);
167 fld[i] = pln.signedDistance(cc[i]);
184 pointDistance_[i] = pln.signedDistance(pts[i]);
190 void Foam::sampledCuttingPlane::combineSurfaces
192 PtrList<isoSurfaceBase>& isoSurfPtrs
195 isoSurfacePtr_.reset(
nullptr);
203 && isoSurfPtrs.size() == 1
207 isoSurfacePtr_.reset(isoSurfPtrs.release(0));
209 else if (isoSurfPtrs.size() == 1)
211 autoPtr<isoSurfaceBase> ptr(isoSurfPtrs.release(0));
215 meshCells_.transfer(surf.meshCells());
226 for (
const auto& surf : isoSurfPtrs)
228 nFaces += surf.size();
229 nPoints += surf.points().size();
234 meshCells_.resize(nFaces);
240 forAll(isoSurfPtrs, surfi)
242 autoPtr<isoSurfaceBase> ptr(isoSurfPtrs.release(surfi));
245 SubList<face> subFaces(newFaces, surf.size(), nFaces);
246 SubList<point> subPoints(newPoints, surf.points().size(),
nPoints);
247 SubList<label> subCells(meshCells_, surf.size(), nFaces);
249 newZones[surfi] = surfZone
257 subFaces = surf.surfFaces();
258 subPoints = surf.points();
259 subCells = surf.meshCells();
263 for (face&
f : subFaces)
265 for (label& pointi :
f)
272 nFaces += subFaces.size();
278 std::move(newPoints),
283 surface_.transfer(combined);
287 if (subMeshPtr_ && meshCells_.size())
290 UIndirectList<label>(subMeshPtr_->cellMap(), meshCells_);
295 void Foam::sampledCuttingPlane::createGeometry()
299 Pout<<
"sampledCuttingPlane::createGeometry :updating geometry."
306 isoSurfacePtr_.reset(
nullptr);
312 pointDistance_.clear();
313 cellDistancePtr_.clear();
315 const bool hasCellZones =
318 const fvMesh& fvm =
static_cast<const fvMesh&
>(this->
mesh());
327 subMeshPtr_.reset(
nullptr);
330 if (!ignoreCellsPtr_)
332 ignoreCellsPtr_.reset(
new bitSet);
336 bitSet select(
mesh().cellZones().selection(zoneNames_));
338 if (select.any() && !select.all())
343 *ignoreCellsPtr_ = std::move(select);
354 ignoreCellsPtr_->clearStorage();
358 ignoreCellsPtr_.reset(
new bitSet);
362 if (!subMeshPtr_ && hasCellZones)
364 const label exposedPatchi =
367 bitSet cellsToSelect(
mesh().cellZones().selection(zoneNames_));
370 <<
"Allocating subset of size "
371 << cellsToSelect.count() <<
" with exposed faces into patch "
372 << exposedPatchi <<
endl;
378 const boundBox& clipBb = isoParams_.getClipBounds();
379 if (clipBb.valid() && cellsToSelect.any())
381 const auto& cellCentres = fvm.C();
383 for (
const label celli : cellsToSelect)
385 const point& cc = cellCentres[celli];
387 if (!clipBb.contains(cc))
389 cellsToSelect.unset(celli);
394 <<
"Bounded subset of size "
395 << cellsToSelect.count() <<
endl;
400 new fvMeshSubset(fvm, cellsToSelect, exposedPatchi)
410 ? subMeshPtr_->subMesh()
414 checkBoundsIntersection(plane_,
mesh.
bounds());
420 cellDistancePtr_.reset
439 setDistanceFields(plane_);
443 Pout<<
"Writing cell distance:" << cellDistance.objectPath() <<
endl;
444 cellDistance.
write();
459 pointDist.primitiveFieldRef() = pointDistance_;
461 Pout<<
"Writing point distance:" << pointDist.objectPath() <<
endl;
468 PtrList<isoSurfaceBase> isoSurfPtrs(offsets_.size());
487 combineSurfaces(isoSurfPtrs);
493 cellDistancePtr_.reset(
nullptr);
494 pointDistance_.clear();
510 const polyMesh&
mesh,
511 const dictionary&
dict
520 isoSurfaceParams::ALGO_TOPO,
521 isoSurfaceParams::filterType::DIAGCELL
523 average_(
dict.getOrDefault(
"average", false)),
524 simpleSubMesh_(
dict.getOrDefault(
"simpleSubMesh", false)),
531 isoSurfacePtr_(nullptr),
533 subMeshPtr_(nullptr),
534 ignoreCellsPtr_(nullptr),
535 cellDistancePtr_(nullptr),
538 dict.readIfPresent(
"offsets", offsets_);
540 if (offsets_.empty())
543 offsets_.first() =
Zero;
546 if (offsets_.size() > 1)
548 const label nOrig = offsets_.size();
552 if (nOrig != offsets_.size())
555 <<
"Removed non-unique offsets" <<
nl;
562 simpleSubMesh_ =
false;
565 if (offsets_.size() > 1)
568 <<
"Multiple offsets with iso-surface (point) not supported"
569 <<
" since needs original interpolators." <<
nl
577 if (!
dict.readIfPresent(
"zones", zoneNames_) &&
dict.found(
"zone"))
580 dict.readEntry(
"zone", zoneNames_.first());
585 dict.readIfPresent(
"exposedPatchName", exposedPatchName_);
588 <<
"Restricting to cellZone(s) " <<
flatOutput(zoneNames_)
589 <<
" with exposed internal faces into patch "
607 Pout<<
"sampledCuttingPlane::expire :"
608 <<
" needsUpdate:" << needsUpdate_ <<
endl;
613 isoSurfacePtr_.reset(
nullptr);
633 Pout<<
"sampledCuttingPlane::update :"
634 <<
" needsUpdate:" << needsUpdate_ <<
endl;
644 needsUpdate_ =
false;
655 return sampleOnFaces(sampler);
665 return sampleOnFaces(sampler);
675 return sampleOnFaces(sampler);
685 return sampleOnFaces(sampler);
695 return sampleOnFaces(sampler);
705 return sampleOnPoints(interpolator);
715 return sampleOnPoints(interpolator);
725 return sampleOnPoints(interpolator);
735 return sampleOnPoints(interpolator);
745 return sampleOnPoints(interpolator);
751 os <<
"sampledCuttingPlane: " <<
name() <<
" :"
752 <<
" plane:" << plane_
757 os <<
" faces:" << faces().size()
758 <<
" points:" <<
points().size();
fvPatchField< scalar > fvPatchScalarField
vectorField pointField
pointField is a vectorField.
virtual const pointField & points() const
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
@ DIAGCELL
Remove pyramid edge points, face-diagonals.
static autoPtr< isoSurfaceBase > New(const isoSurfaceParams ¶ms, const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const bitSet &ignoreCells=bitSet())
A class for handling words, derived from Foam::string.
label findIndex(const wordRe &key) const
void resize(const label len)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
A class for managing temporary objects.
static constexpr const zero Zero
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
algorithmType algorithm() const noexcept
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
static word timeName(const scalar t, const int precision=precision_)
void inplaceUniqueSort(ListType &input)
const polyBoundaryMesh & boundaryMesh() const
Ostream & endl(Ostream &os)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
Mesh consisting of general polyhedral cells.
label nPoints() const noexcept
virtual void print(Ostream &os, int level=0) const
virtual bool needsUpdate() const
bool interpolate() const noexcept
const boundBox & getClipBounds() const noexcept
const cellZoneMesh & cellZones() const noexcept
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
An abstract class for surfaces with sampling.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
label findPatchID(const word &patchName, const bool allowNotFound=true) const
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const volVectorField & C() const
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 word & name() const noexcept
virtual bool write(const token &tok)=0
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
sampledCuttingPlane(const word &name, const polyMesh &mesh, const dictionary &dict)
errorManipArg< error, int > exit(error &err, const int errNo=1)
const boundBox & bounds() const
const fvBoundaryMesh & boundary() const
static word defaultName(const label n=-1)
const vectorField & cellCentres() const
MeshedSurface< face > meshedSurface
List< face > faceList
A List of faces.
const polyMesh & mesh() const noexcept
const vectorField & faceCentres() const
Ostream & print(Ostream &os, UIntType value, char off='0', char on='1')
virtual void clearGeom() const
List< surfZone > surfZoneList
word name(const expressions::valueTypeCode typeCode)
const Time & time() const
SubField< Type > subField
#define FatalIOErrorInFunction(ios)
virtual tmp< scalarField > sample(const interpolation< scalar > &sampler) const
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
#define IOWarningInFunction(ios)
vector point
Point is a vector.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
const Boundary & boundaryField() const
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField