44#include "MarchingCubes.h"
53int main(
int argc,
char *argv[])
57 "Re-sample surfaces used in foamyHexMesh operation"
70 Info<<
"Reading surfaces as specified in the foamyHexMeshDict and"
71 <<
" writing re-sampled " <<
n <<
" to " << exportName
93 "cvSearchableSurfaces",
100 foamyHexMeshDict.subDict(
"geometry"),
101 foamyHexMeshDict.getOrDefault(
"singleRegionName",
true)
104 Info<<
"Geometry read in = "
105 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
115 foamyHexMeshDict.subDict(
"surfaceConformation")
118 Info<<
"Set up geometry in = "
119 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
127 Info<<
"Meshing to bounding box " << bb <<
nl <<
endl;
137 MarchingCubes mc(span.x(), span.y(), span.z() ) ;
138 mc.set_resolution(
n.x(),
n.y(),
n.z());
147 for(
int k = 0 ;
k < mc.size_z() ;
k++ )
149 pt.
z() = bb.
min().
z() +
k*d.z();
150 for(
int j = 0 ; j < mc.size_y() ; j++ )
152 pt.
y() = bb.
min().
y() + j*d.y();
153 for(
int i = 0 ; i < mc.size_x() ; i++ )
155 pt.
x() = bb.
min().
x() + i*d.x();
161 Info<<
"Generated " <<
points.size() <<
" sampling points in = "
162 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
167 const labelList& surfaces = geometryToConformTo.surfaces();
177 searchableSurface::OUTSIDE,
185 for(
int k = 0 ;
k < mc.size_z() ;
k++ )
187 for(
int j = 0 ; j < mc.size_y() ; j++ )
189 for(
int i = 0 ; i < mc.size_x() ; i++ )
191 mc.set_data(
float(signedDist[pointi++]), i, j,
k);
196 Info<<
"Determined signed distance in = "
197 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
202 Info<<
"Constructed iso surface in = "
203 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
214 List<labelledTri> tris(mc.ntrigs());
234 bb.
min().
x() + v.x*span.x()/mc.size_x(),
235 bb.
min().
y() + v.y*span.y()/mc.size_y(),
236 bb.
min().
z() + v.z*span.z()/mc.size_z()
246 const wordList& regions = geometry[surfaces[i]].regions();
247 regionOffsets[i] = nRegions;
248 nRegions += regions.
size();
256 const wordList& regions = geometry[surfaces[i]].regions();
262 geometry[surfaces[i]].
name() +
"_" + regions[regionI],
272 Info<<
"Extracted triSurface in = "
273 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
278 List<pointIndexHit> hitInfo;
280 geometryToConformTo.findSurfaceNearest
299 if (hitSurfaces[triI] == surfI)
301 surfInfo.append(hitInfo[triI]);
302 surfIndices.append(triI);
308 geometry[surfaces[surfI]].getRegion(surfInfo, region);
312 label triI = surfIndices[i];
313 s[triI].region() = regionOffsets[surfI]+region[i];
318 Info<<
"Re-patched surface in = "
319 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
336 Info<<
"writing surfMesh:\n "
337 << smesh.searchableSurface::objectPath() <<
nl <<
endl;
338 smesh.searchableSurface::write();
340 Info<<
"Written surface in = "
341 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
CGAL::Triangle_3< K > Triangle
graph_traits< Graph >::vertex_descriptor Vertex
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
void size(const label n)
Older name for setAddressableSize.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
static void addNote(const string ¬e)
Add extra notes for the usage information.
void grow(const scalar delta)
Expand box by adjusting min/max by specified amount in each dimension.
const point & min() const noexcept
Minimum describing the bounding box.
vector span() const
The bounding box span (from minimum to maximum).
A class for handling file names.
Identifies a surface patch/zone by name and index, with geometric type.
A triFace with additional (region) index.
static void signedDistance(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest, const pointField &samples, const scalarField &nearestDistSqr, const volumeType illegalHandling, labelList &nearestSurfaces, scalarField &distance)
Find signed distance to nearest surface. Outside is positive.
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Implements a timeout mechanism via sigalarm.
Standard boundBox with extra functionality for use in octree.
IOoject and searching on triSurface.
Triangulated surface description with patch information.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const polyBoundaryMesh & patches
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))
List< word > wordList
List of word.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Vector< label > labelVector
Vector of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
pointField vertices(const blockVertexList &bvl)
vector point
Point is a vector.
List< geometricSurfacePatch > geometricSurfacePatchList
List of geometricSurfacePatch.
vectorField pointField
pointField is a vectorField.
cpuTimePosix cpuTime
Selection of preferred clock mechanism for the elapsed cpu time.
constexpr char nl
The newline '\n' character (0x0a).
std::vector< Triangle > triangles
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.