Go to the documentation of this file.
44 #include "MarchingCubes.h"
53 int 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 bb.
min() -= smallVec;
128 bb.
max() += smallVec;
131 Info<<
"Meshing to bounding box " << bb <<
nl <<
endl;
141 MarchingCubes mc(span.x(), span.y(), span.z() ) ;
142 mc.set_resolution(
n.x(),
n.y(),
n.z());
151 for(
int k = 0 ;
k < mc.size_z() ;
k++ )
153 pt.
z() = bb.
min().
z() +
k*d.z();
154 for(
int j = 0 ; j < mc.size_y() ; j++ )
156 pt.
y() = bb.
min().
y() + j*d.y();
157 for(
int i = 0 ; i < mc.size_x() ; i++ )
159 pt.
x() = bb.
min().
x() + i*d.x();
165 Info<<
"Generated " <<
points.size() <<
" sampling points in = "
166 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
171 const labelList& surfaces = geometryToConformTo.surfaces();
181 searchableSurface::OUTSIDE,
189 for(
int k = 0 ;
k < mc.size_z() ;
k++ )
191 for(
int j = 0 ; j < mc.size_y() ; j++ )
193 for(
int i = 0 ; i < mc.size_x() ; i++ )
195 mc.set_data(
float(signedDist[pointi++]), i, j,
k);
200 Info<<
"Determined signed distance in = "
201 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
206 Info<<
"Constructed iso surface in = "
207 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
217 Triangle* triangles = mc.triangles();
238 bb.
min().
x() + v.x*span.x()/mc.size_x(),
239 bb.
min().
y() + v.y*span.y()/mc.size_y(),
240 bb.
min().
z() + v.z*span.z()/mc.size_z()
246 labelList regionOffsets(surfaces.size());
250 const wordList& regions = geometry[surfaces[i]].regions();
251 regionOffsets[i] = nRegions;
252 nRegions += regions.size();
260 const wordList& regions = geometry[surfaces[i]].regions();
266 geometry[surfaces[i]].
name() +
"_" + regions[regionI],
276 Info<<
"Extracted triSurface in = "
277 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
284 geometryToConformTo.findSurfaceNearest
303 if (hitSurfaces[triI] == surfI)
305 surfInfo.append(hitInfo[triI]);
306 surfIndices.append(triI);
312 geometry[surfaces[surfI]].getRegion(surfInfo, region);
316 label triI = surfIndices[i];
317 s[triI].region() = regionOffsets[surfI]+region[i];
322 Info<<
"Re-patched surface in = "
323 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
340 Info<<
"writing surfMesh:\n "
341 << smesh.searchableSurface::objectPath() <<
nl <<
endl;
342 smesh.searchableSurface::write();
344 Info<<
"Written surface in = "
345 <<
timer.cpuTimeIncrement() <<
" s." <<
nl <<
endl;
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
CGAL::Triangle_3< K > Triangle
A class for handling file names.
const point & max() 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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Identifies a surface patch/zone by name and index, with geometric type.
Standard boundBox with extra functionality for use in octree.
graph_traits< Graph >::vertex_descriptor Vertex
static void addNote(const string ¬e)
IOoject and searching on triSurface.
Ostream & endl(Ostream &os)
Starts timing CPU usage and return elapsed time from start.
T get(const label index) const
static void addArgument(const string &argName, const string &usage="")
static void noFunctionObjects(bool addWithOption=false)
Generic templated field type.
Triangulated surface description with patch information.
static void signedDistance(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest, const pointField &samples, const scalarField &nearestDistSqr, const volumeType illegalHandling, labelList &nearestSurfaces, scalarField &distance)
const point & min() const
pointField vertices(const blockVertexList &bvl)
const word & system() const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static int myProcNo(const label communicator=worldComm)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Implements a timeout mechanism via sigalarm.
const polyBoundaryMesh & patches
A triFace with additional (region) index.
word name(const expressions::valueTypeCode typeCode)
vector point
Point is a vector.
const word & constant() const
Foam::argList args(argc, argv)