voxelMeshSearch.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "voxelMeshSearch.H"
29#include "polyMesh.H"
30#include "IOobject.H"
31#include "fvMesh.H"
32#include "block.H"
33#include "emptyPolyPatch.H"
34#include "processorPolyPatch.H"
35#include "fvMeshTools.H"
36
37/* * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * */
38
39namespace Foam
40{
42}
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
48(
49 const labelVector& nDivs
50)
51{
52 return labelVector(1, nDivs.x(), nDivs.x()*nDivs.y());
53}
54
55
57(
58 const labelVector& nDivs,
59 const labelVector& voxel
60)
61{
62 return voxel.x()+voxel.y()*nDivs.x()+voxel.z()*nDivs.x()*nDivs.y();
63}
64
65
67(
68 const labelVector& nDivs,
69 label voxeli
70)
71{
72 const label nxy = nDivs.x()*nDivs.y();
73
74 labelVector voxel;
75 voxel.z() = voxeli/nxy;
76 voxeli = voxeli % nxy;
77 voxel.y() = voxeli/nDivs.x();
78 voxel.x() = voxeli%nDivs.x();
79
80 return voxel;
81}
82
83
85(
86 const boundBox& bb,
87 const labelVector& g,
88 const point& pt
89)
90{
91 const vector s(cmptDivide(bb.span(), vector(g.x(), g.y(), g.z())));
92
94 (
95 floor((pt.x()-bb.min().x())/s.x()),
96 floor((pt.y()-bb.min().y())/s.y()),
97 floor((pt.z()-bb.min().z())/s.z())
98 );
99
100 return v;
101}
102
103
105(
106 const boundBox& bb,
107 const labelVector& g,
108 const point& pt,
109 const bool clip
110)
111{
112 const vector s(cmptDivide(bb.span(), vector(g.x(), g.y(), g.z())));
113
115 (
116 floor((pt.x()-bb.min().x())/s.x()),
117 floor((pt.y()-bb.min().y())/s.y()),
118 floor((pt.z()-bb.min().z())/s.z())
119 );
120
121 if (clip)
122 {
123 v[0] = max(0, min(g[0]-1, v[0]));
124 v[1] = max(0, min(g[1]-1, v[1]));
125 v[2] = max(0, min(g[2]-1, v[2]));
126 }
127 else if
128 (
129 v[0] < 0
130 || v[1] < 0
131 || v[2] < 0
132 || v[0] >= g[0]
133 || v[1] >= g[1]
134 || v[2] >= g[2]
135 )
136 {
137 return -1;
138 }
139
140 return index(g, v);
141}
142
143
145(
146 const boundBox& bb,
147 const labelVector& g,
148 const labelVector& voxel
149)
150{
151 const vector s(cmptDivide(bb.span(), vector(g.x(), g.y(), g.z())));
152
153 return bb.min()+0.5*s+point(voxel[0]*s[0], voxel[1]*s[1], voxel[2]*s[2]);
154}
155
156
158(
159 OBJstream& os,
160 const boundBox& bb,
161 const labelVector& g
162)
163{
164 const vector s(cmptDivide(bb.span(), vector(g.x(), g.y(), g.z())));
165
166 for (label i = 1; i < g[0]; i++)
167 {
168 for (label j = 0; j < g[1]; j++)
169 {
170 for (label k = 0; k < g[2]; k++)
171 {
172 point p1(bb.min()+point((i-1)*s[0], j*s[1], k*s[2]));
173 point p2(bb.min()+point(i*s[0], j*s[1], k*s[2]));
174 os.write(linePointRef(p1, p2));
175 }
176 }
177 }
178 for (label i = 0; i < g[0]; i++)
179 {
180 for (label j = 1; j < g[1]; j++)
181 {
182 for (label k = 0; k < g[2]; k++)
183 {
184 point p1(bb.min()+point(i*s[0], (j-1)*s[1], k*s[2]));
185 point p2(bb.min()+point(i*s[0], j*s[1], k*s[2]));
186 os.write(linePointRef(p1, p2));
187 }
188 }
189 }
190 for (label i = 0; i < g[0]; i++)
191 {
192 for (label j = 0; j < g[1]; j++)
193 {
194 for (label k = 1; k < g[2]; k++)
195 {
196 point p1(bb.min()+point(i*s[0], j*s[1], (k-1)*s[2]));
197 point p2(bb.min()+point(i*s[0], j*s[1], k*s[2]));
198 os.write(linePointRef(p1, p2));
199 }
200 }
201 }
202}
203
204
205Foam::label Foam::voxelMeshSearch::searchProcPatch
206(
207 const label faceID,
208 const point& searchPoint
209) const
210{
211 const pointField& cellCentres = mesh_.cellCentres();
212 const polyBoundaryMesh& bMeshes = mesh_.boundaryMesh();
213
214 label patchi = bMeshes.patchID()[faceID-mesh_.nInternalFaces()];
215 const polyPatch& bMeshPatch = bMeshes[patchi];
216
217 if (!isA<processorPolyPatch>(bMeshPatch))
218 {
219 return -1;
220 }
221 else
222 {
223 // Find nearest cell. Linear search since cheaper than constructing
224 // tree?
225 const labelUList& faceCells = bMeshPatch.faceCells();
226 scalar minProximity = GREAT;
227
228 label nearestCellI = -1;
229 forAll(faceCells, i)
230 {
231 const point& cc = cellCentres[faceCells[i]];
232 scalar proximity = magSqr(cc-searchPoint);
233 if (proximity < minProximity)
234 {
235 minProximity = proximity;
236 nearestCellI = faceCells[i];
237 }
238 }
239 return nearestCellI;
240 }
241}
242
243
244Foam::label Foam::voxelMeshSearch::findIntersectedFace
245(
246 const label cellI,
247 const point& p
248) const
249{
250 // Return -1 or the label of the face intersected when tracking from
251 // p to centre of cellI
252
253 const faceList& faces = mesh_.faces();
254 const pointField& faceCentres = mesh_.faceCentres();
255 const pointField& points = mesh_.points();
256
257 const point& cc = mesh_.cellCentres()[cellI];
258 const labelList& cFaces = mesh_.cells()[cellI];
259
260 const vector q(cc-p);
261
262 forAll(cFaces, cFacei)
263 {
264 label facei = cFaces[cFacei];
265
266 pointHit hitInfo = faces[facei].intersection
267 (
268 p,
269 q,
270 faceCentres[facei],
271 points,
273 );
274
275 if (hitInfo.hit() && (hitInfo.distance() < 1))
276 {
277 return facei;
278 }
279 }
280 return -1;
281}
282
283
284// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
285
287(
288 const polyMesh& mesh,
289 const bool doUpdate
290)
291:
292 mesh_(mesh)
293{
294 // Determine number of voxels from number of cells in mesh
295 const labelVector& dim = mesh_.geometricD();
296
297 // Guarantee at least one voxel
298 label nCells = max(1, mesh_.nCells());
299
300 label nDivs = -1;
301 if (mesh_.nGeometricD() == 1)
302 {
303 nDivs = nCells;
304 }
305 else if (mesh_.nGeometricD() == 2)
306 {
307 nDivs = label(Foam::sqrt(scalar(nCells)));
308 }
309 else
310 {
311 nDivs = label(Foam::cbrt(scalar(nCells)));
312 }
313
314 nDivs_ = labelVector(nDivs, nDivs, nDivs);
315 forAll(dim, i)
316 {
317 if (dim[i] == -1)
318 {
319 nDivs_[i] = 1;
320 }
321 }
322
323 // Redo the local bounding box
324 localBb_ = boundBox(mesh_.points(), false);
325
326 const point eps(1e-10, 1e-10, 1e-10);
327
328 localBb_.min() = localBb_.min()-eps;
329 localBb_.max() = localBb_.max()+eps;
330
331 if (debug)
332 {
333 Pout<< "voxelMeshSearch : mesh:" << mesh_.name()
334 << " nDivs:" << nDivs_ << endl;
335 }
336
337 if (doUpdate)
338 {
339 update();
340 }
341}
342
343
345(
346 const polyMesh& mesh,
347 const boundBox& localBb,
348 const labelVector& nDivs,
349 const bool doUpdate
350)
351:
352 mesh_(mesh),
353 localBb_(localBb),
354 nDivs_(nDivs)
355{
356 if (doUpdate)
357 {
358 update();
359 }
360}
361
362
363// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
364
366{
367 // Initialise seed cell array
368
369 seedCell_.setSize(cmptProduct(nDivs_));
370 seedCell_ = -1;
371
372
373 // Find seed cells
374 const pointField& points = mesh_.points();
375 const labelListList& cellPoints = mesh_.cellPoints();
376
377 forAll(cellPoints, celli)
378 {
379 const labelList& cPoints = cellPoints[celli];
380
381 // Get cell bounding box
382 boundBox bb(points, cPoints, false);
383
384 fill(seedCell_, localBb_, nDivs_, bb, celli);
385 }
386
387
388 if (debug)
389 {
390 Pout<< "voxelMeshSearch : mesh:" << mesh_.name()
391 << " nDivs:" << nDivs_
392 << " localBb:" << localBb_ << endl;
393 }
394
395
398 //const pointField& cellCentres = mesh_.cellCentres();
399 //forAll(cellCentres, celli)
400 //{
401 // label voxeli = index(cellCentres[celli]);
402 // seedCell_[voxeli] = celli;
403 //}
404
405 return true;
406}
407
408
409Foam::label Foam::voxelMeshSearch::findCell(const point& p) const
410{
411 // First check if the point is contained in the bounding box, else exit
412 if (!localBb_.contains(p))
413 {
414 return -1;
415 }
416
417
418 // Locate the voxel index for this point. Do not clip.
419 label voxeli = index(localBb_, nDivs_, p, false);
420
421 // The point may still be inside the bb but outside the actual domain.
422 if (voxeli < 0)
423 {
424 return -1;
425 }
426 else
427 {
428 // Inverse map to compute the seed cell.
429 label celli = seedCell_[voxeli];
430
431 if (celli < 0)
432 {
433 return -1;
434 }
435 else
436 {
437 // Simplified, non-parallel tracking from cell centre of
438 // celli to wanted location p. Note that the cell thus
439 // found does not have to be the absolute 'correct' one as
440 // long as at least one of the processors finds a cell.
441
442 track_.clear();
443 while (true)
444 {
445 if (track_.size() < 5)
446 {
447 track_.append(celli);
448 }
449
450 // I am in celli now. How many faces do I have ?
451 label facei = findIntersectedFace(celli, p);
452
453 if (facei == -1)
454 {
455 return celli;
456 }
457
458 const label startOfTrack(max(0, track_.size()-5));
459
460 label nextCell;
461 if (mesh_.isInternalFace(facei))
462 {
463 label own = mesh_.faceOwner()[facei];
464 label nei = mesh_.faceNeighbour()[facei];
465 nextCell = (own == celli ? nei : own);
466
467 if (track_.found(nextCell, startOfTrack))
468 {
469 return celli;
470 }
471 }
472 else
473 {
474 nextCell = searchProcPatch(facei, p);
475
476 if (nextCell == -1 || nextCell == celli)
477 {
478 return nextCell;
479 }
480 else if (track_.found(nextCell, startOfTrack))
481 {
482 return -1; // point is really out
483 }
484 }
485
486 celli = nextCell;
487 }
488 return -1;
489 }
490 }
491}
492
493
495(
496 const IOobject& io
497) const
498{
500
504 {
505 //Info<< "Creating block" << endl;
506
507 block b
508 (
510 localBb_.points(),
513 nDivs_,
515 );
516
517 cellShapes = b.shapes();
518
519 //Info<< "Creating boundary faces" << endl;
520
521 boundary.setSize(b.boundaryPatches().size());
522 forAll(boundary, patchi)
523 {
524 faceList faces(b.boundaryPatches()[patchi].size());
525 forAll(faces, facei)
526 {
527 faces[facei] = face(b.boundaryPatches()[patchi][facei]);
528 }
529 boundary[patchi].transfer(faces);
530 }
531
532 points.transfer(const_cast<pointField&>(b.points()));
533 }
534
535 //Info<< "Creating patch dictionaries" << endl;
537 forAll(patchNames, patchi)
538 {
539 patchNames[patchi] = polyPatch::defaultName(patchi);
540 }
541
543 forAll(boundaryDicts, patchi)
544 {
545 boundaryDicts.set(patchi, new dictionary());
546 dictionary& patchDict = boundaryDicts[patchi];
547 patchDict.add("type", emptyPolyPatch::typeName);
548 }
549
550 //Info<< "Creating polyMesh" << endl;
551 IOobject polyIO(io);
554 (
555 //IOobject
556 //(
557 // polyMesh::defaultRegion,
558 // runTime.constant(),
559 // runTime,
560 // IOobject::NO_READ
561 //),
562 polyIO,
563 std::move(points),
565 boundary,
568 "defaultFaces",
570 false
571 );
572
573 //Info<< "Writing polyMesh" << endl;
574 mesh.write();
575
576 //Info<< "Reading fvMesh" << endl;
577
579 (
580 io.db(),
581 io.name()
582 );
583 IOobject fvIO(io);
585
586 return autoPtr<fvMesh>::New(fvIO);
587}
588
589
590// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
label k
const uniformDimensionedVectorField & g
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
readOption readOpt() const noexcept
The read option.
Definition: IOobjectI.H:164
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:500
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
OFstream that keeps track of vertices.
Definition: OBJstream.H:61
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
scalar centre() const
Mid-point location, zero for an empty list.
Definition: PDRblockI.H:67
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:61
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:127
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:73
An analytical geometric cellShape.
Definition: cellShape.H:72
pointField points(const UList< point > &meshPoints) const
The points corresponding to this shape.
Definition: cellShapeI.H:151
static word defaultName
The default cloud name: defaultCloud.
Definition: cloud.H:90
reference ref() const
A reference to the entry (Error if not found)
Definition: dictionary.H:225
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:640
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static void createDummyFvMeshFiles(const objectRegistry &parent, const word &regionName, const bool verbose=false)
Create additional fvSchemes/fvSolution files.
Definition: fvMeshTools.C:1064
virtual bool write(const bool valid=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1079
const vector & offset() const noexcept
Offset vector (from patch faces to destination mesh objects)
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const labelList & patchID() const
Per boundary face label the patch index.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:883
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:872
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:371
label nCells() const noexcept
Number of mesh cells.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Fast, non-parallel searching in mesh without use of octree.
static labelVector index3(const labelVector &nDivs, const label voxeli)
Combined voxel index to individual indices.
label findCell(const point &) const
Find a cell.
bool update()
Update lookup tables for geometry changes.
const labelVector & nDivs() const
Number of voxels for local mesh.
static void writeGrid(OBJstream &, const boundBox &, const labelVector &)
Debug: write all edges.
autoPtr< fvMesh > makeMesh(const IOobject &) const
Debug: construct fvMesh. Note: writes a dummy mesh to.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const cellModel & hex
Info<< "Creating cells"<< endl;cellShapes=b.shapes();Info<< "Creating boundary faces"<< endl;boundary.setSize(b.boundaryPatches().size());forAll(boundary, patchi) { faceList faces(b.boundaryPatches()[patchi].size());forAll(faces, facei) { faces[facei]=face(b.boundaryPatches()[patchi][facei]);} boundary[patchi].transfer(faces);} points.transfer(const_cast< pointField & >(b.points()));}Info<< "Creating patch dictionaries"<< endl;wordList patchNames(boundary.size());forAll(patchNames, patchi){ patchNames[patchi]=polyPatch::defaultName(patchi);}PtrList< dictionary > boundaryDicts(boundary.size())
cellShapeList cellShapes
faceListList boundary
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionSet clip(const dimensionSet &ds1, const dimensionSet &ds2)
Definition: dimensionSet.C:278
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
PtrList< blockFace > blockFaceList
A PtrList of blockFaces.
Definition: blockFaceList.H:47
List< label > labelList
A List of labels.
Definition: List.H:66
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
Vector< label > labelVector
Vector of labels.
Definition: labelVector.H:51
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Cmpt cmptProduct(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:598
vector point
Point is a vector.
Definition: point.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
PointHit< point > pointHit
A PointIndexHit for 3D points.
Definition: pointHit.H:44
dimensionedScalar cbrt(const dimensionedScalar &ds)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
wordList patchNames(nPatches)
volScalarField & b
Definition: createFields.H:27
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.