simple project that make use of the MC and EMC algorithms
This commit is contained in:
parent
1509a9b434
commit
29e46cb2af
|
@ -0,0 +1,22 @@
|
|||
#ifndef ___DEFINITIONS
|
||||
#define ___DEFINITIONS
|
||||
|
||||
#include <vcg/simplex/vertex/with/afvmvn.h>
|
||||
#include <vcg/simplex/face/with/afavfn.h>
|
||||
#include <vcg/simplex/edge/with/ae.h>
|
||||
#include <vcg/complex/trimesh/base.h>
|
||||
#include <vcg/complex/trimesh/allocate.h>
|
||||
|
||||
typedef float ScalarType;
|
||||
|
||||
class Edge;
|
||||
class Face;
|
||||
class Vertex : public vcg::VertexAFVMVN< ScalarType, Edge, Face > {};
|
||||
class Face : public vcg::FaceAFAVFN< Vertex, Edge, Face> {};
|
||||
class Mesh : public vcg::tri::TriMesh< std::vector< Vertex>, std::vector< Face > > {};
|
||||
|
||||
typedef vcg::tri::Allocator< Mesh > Allocator;
|
||||
typedef vcg::Box3< int > BoundingBox;
|
||||
typedef Vertex* VertexPointer;
|
||||
|
||||
#endif //___DEFINITIONS
|
|
@ -0,0 +1,18 @@
|
|||
#ifndef __EXTRS_IMPLICIT
|
||||
#define __EXTRS_IMPLICIT
|
||||
|
||||
#include <vcg/space/point3.h>
|
||||
|
||||
class Implicit
|
||||
{
|
||||
public:
|
||||
Implicit() {};
|
||||
|
||||
virtual ~Implicit() {};
|
||||
|
||||
virtual float V(int x, int y, int z) const = 0;
|
||||
|
||||
virtual vcg::Point3f N(float x, float y, float z) const = 0;
|
||||
};
|
||||
|
||||
#endif // __EXTRS_IMPLICIT
|
|
@ -0,0 +1,95 @@
|
|||
#ifndef __VCGTEST_IMPLICITSPHERE
|
||||
#define __VCGTEST_IMPLICITSPHERE
|
||||
|
||||
class ImplicitSphere
|
||||
{
|
||||
public:
|
||||
ImplicitSphere()
|
||||
{
|
||||
_center.Zero();
|
||||
_radius = _sqr_radius = 0.0;
|
||||
};
|
||||
|
||||
ImplicitSphere(vcg::Point3f ¢er, float radius)
|
||||
{
|
||||
_center = center;
|
||||
_radius = radius;
|
||||
_sqr_radius = radius*radius;
|
||||
};
|
||||
|
||||
ImplicitSphere(const ImplicitSphere &sphere)
|
||||
{
|
||||
_center = sphere._center;
|
||||
_radius = sphere._radius;
|
||||
_sqr_radius = sphere._sqr_radius;
|
||||
};
|
||||
|
||||
ImplicitSphere& operator=(const ImplicitSphere &sphere)
|
||||
{
|
||||
if (this != &sphere)
|
||||
{
|
||||
_center = sphere._center;
|
||||
_radius = sphere._radius;
|
||||
_sqr_radius = sphere._sqr_radius;
|
||||
}
|
||||
return *this;
|
||||
};
|
||||
|
||||
~ImplicitSphere()
|
||||
{};
|
||||
|
||||
bool operator!=(const ImplicitSphere &sphere)
|
||||
{
|
||||
return (sphere._center!=_center && sphere._radius!=_radius);
|
||||
};
|
||||
|
||||
|
||||
float V(int x, int y, int z) const
|
||||
{
|
||||
vcg::Point3f point((float) x, (float) y, (float) z);
|
||||
return (_center-point).Norm() - _radius;
|
||||
};
|
||||
|
||||
bool DirectedDistance(const vcg::Point3i &p1, const vcg::Point3i &p2, vcg::Point3f &v, vcg::Point3f &n, float &dist)
|
||||
{
|
||||
vcg::Point3f orig, dir;
|
||||
orig.X() = (float) p1.X();
|
||||
orig.Y() = (float) p1.Y();
|
||||
orig.Z() = (float) p1.Z();
|
||||
dir.X() = (float) p2.X()-p1.X();
|
||||
dir.Y() = (float) p2.Y()-p1.Y();
|
||||
dir.Z() = (float) p2.Z()-p1.Z();
|
||||
|
||||
double a = dir.SquaredNorm();
|
||||
double b = 2.0*(dir*(orig - _center));
|
||||
double c = (orig - _center).SquaredNorm() - _radius*_radius;
|
||||
double d = b*b - 4.0*a*c;
|
||||
|
||||
if (d >= 0)
|
||||
{
|
||||
d = sqrt(d);
|
||||
|
||||
double t1 = (-b-d) / (2.0*a);
|
||||
double t2 = (-b+d) / (2.0*a);
|
||||
double t = 1.00001;
|
||||
if (t1 >= 0.0 && t1 < t) t = t1;
|
||||
if (t2 >= 0.0 && t2 < t) t = t2;
|
||||
|
||||
if (t != 1.00001)
|
||||
{
|
||||
v = (vcg::Point3f) (orig + dir*((float)t));
|
||||
n = (vcg::Point3f) ((v - _center) / _radius);
|
||||
dist = (float) ((dir*n) < 0.0) ? dir.Norm()*t : -dir.Norm()*t;
|
||||
return true;
|
||||
}
|
||||
}
|
||||
return false;
|
||||
};
|
||||
|
||||
private:
|
||||
vcg::Point3f _center;
|
||||
float _radius;
|
||||
float _sqr_radius;
|
||||
};
|
||||
|
||||
#endif // __VCGTEST_IMPLICITSPHERE
|
|
@ -0,0 +1,69 @@
|
|||
#ifndef __VCGTEST_SPHEREDIFFERENCE
|
||||
#define __VCGTEST_SPHEREDIFFERENCE
|
||||
|
||||
class SphereDifference
|
||||
{
|
||||
public:
|
||||
SphereDifference()
|
||||
{}
|
||||
|
||||
SphereDifference( const SphereDifference &sphere_difference)
|
||||
{
|
||||
_union = sphere_difference._union;
|
||||
_sphere = sphere_difference._sphere;
|
||||
}
|
||||
|
||||
SphereDifference( const SphereUnion &sphere_union, const ImplicitSphere &sphere)
|
||||
{
|
||||
_union = sphere_union;
|
||||
_sphere = sphere;
|
||||
}
|
||||
|
||||
float V(int x, int y, int z)
|
||||
{
|
||||
return vcg::math::Max<float>(_union.V(x, y, z), -_sphere.V(x, y, z));
|
||||
}
|
||||
|
||||
inline bool DirectedDistance(const vcg::Point3i p1, const vcg::Point3i p2, vcg::Point3f &p, vcg::Point3f &n, float &d)
|
||||
{
|
||||
vcg::Point3f v1, n1;
|
||||
vcg::Point3f v2, n2;
|
||||
float d1, d2;
|
||||
|
||||
bool ok1 = _union.DirectedDistance(p1, p2, v1, n1, d1);
|
||||
bool ok2 = _sphere.DirectedDistance(p1, p2, v2, n2, d2);
|
||||
d2 = -d2;
|
||||
|
||||
if (ok1 && ok2)
|
||||
{
|
||||
if (d1 > d2)
|
||||
ok2 = false;
|
||||
else
|
||||
ok1 = false;
|
||||
}
|
||||
|
||||
if (ok1)
|
||||
{
|
||||
p = v1;
|
||||
n = n1;
|
||||
d = d1;
|
||||
return true;
|
||||
}
|
||||
else if (ok2)
|
||||
{
|
||||
p = v2;
|
||||
n = n2;
|
||||
d = d2;
|
||||
return true;
|
||||
}
|
||||
|
||||
return false;
|
||||
}
|
||||
|
||||
private:
|
||||
SphereUnion _union;
|
||||
ImplicitSphere _sphere;
|
||||
|
||||
};
|
||||
|
||||
#endif // __VCGTEST_SPHEREDIFFERENCE
|
|
@ -0,0 +1,84 @@
|
|||
#ifndef __VCGTEST_SPHEREUNION
|
||||
#define __VCGTEST_SPHEREUNION
|
||||
|
||||
class SphereUnion
|
||||
{
|
||||
public:
|
||||
SphereUnion()
|
||||
{};
|
||||
|
||||
SphereUnion(const ImplicitSphere &sphere1, const ImplicitSphere &sphere2)
|
||||
{
|
||||
_sphere1 = sphere1;
|
||||
_sphere2 = sphere2;
|
||||
};
|
||||
|
||||
SphereUnion(const SphereUnion &sphere_union)
|
||||
{
|
||||
_sphere1 = sphere_union._sphere1;
|
||||
_sphere2 = sphere_union._sphere2;
|
||||
}
|
||||
|
||||
SphereUnion& operator=(const SphereUnion &sphere_union)
|
||||
{
|
||||
if (this != &sphere_union)
|
||||
{
|
||||
_sphere1 = sphere_union._sphere1;
|
||||
_sphere2 = sphere_union._sphere2;
|
||||
}
|
||||
return *this;
|
||||
}
|
||||
|
||||
bool operator!=(const SphereUnion &sphere_union)
|
||||
{
|
||||
bool comp1 = _sphere1 != sphere_union._sphere1;
|
||||
bool comp2 = _sphere2 != sphere_union._sphere2;
|
||||
return (comp1 && comp2);
|
||||
}
|
||||
|
||||
float V(int x, int y, int z)
|
||||
{
|
||||
return vcg::math::Min<float>(_sphere1.V(x, y, z), _sphere2.V(x, y, z));
|
||||
};
|
||||
|
||||
bool DirectedDistance(const vcg::Point3i &p1, const vcg::Point3i &p2, vcg::Point3f &v, vcg::Point3f &n, float &d)
|
||||
{
|
||||
vcg::Point3f v1, n1;
|
||||
vcg::Point3f v2, n2;
|
||||
float d1, d2;
|
||||
|
||||
bool ok1 = _sphere1.DirectedDistance(p1, p2, v1, n1, d1);
|
||||
bool ok2 = _sphere2.DirectedDistance(p1, p2, v2, n2, d2);
|
||||
|
||||
if (ok1 && ok2)
|
||||
{
|
||||
if (d1 < d2)
|
||||
ok2 = false;
|
||||
else
|
||||
ok1 = false;
|
||||
}
|
||||
|
||||
if (ok1)
|
||||
{
|
||||
v = v1;
|
||||
n = n1;
|
||||
d = d1;
|
||||
return true;
|
||||
}
|
||||
else if (ok2)
|
||||
{
|
||||
v = v2;
|
||||
n = n2;
|
||||
d = d2;
|
||||
return true;
|
||||
}
|
||||
else
|
||||
return false;
|
||||
};
|
||||
|
||||
private:
|
||||
ImplicitSphere _sphere1;
|
||||
ImplicitSphere _sphere2;
|
||||
};
|
||||
|
||||
#endif // __VCGTEST_SPHEREUNION
|
|
@ -0,0 +1,91 @@
|
|||
#ifndef __VCGTEST_VOLUME
|
||||
#define __VCGTEST_VOLUME
|
||||
|
||||
#include "ImplicitSphere.h"
|
||||
#include "SphereUnion.h"
|
||||
#include "SphereDifference.h"
|
||||
|
||||
class Volume
|
||||
{
|
||||
public:
|
||||
Volume()
|
||||
{
|
||||
ImplicitSphere s1(vcg::Point3f(-5.0, 0.0, 0.0), 10.0);
|
||||
ImplicitSphere s2(vcg::Point3f( 5.0, 5.0, 3.0), 7.0);
|
||||
ImplicitSphere s3(vcg::Point3f( 1.0, 0.0, 10.0), 6.0);
|
||||
SphereUnion sphere_union(s1, s2);
|
||||
SphereDifference sphere_difference(sphere_union, s3);
|
||||
_sphere_diff = sphere_difference;
|
||||
}
|
||||
|
||||
float V(const int pi, const int pj, const int pk)
|
||||
{
|
||||
return _sphere_diff.V(pi, pj, pk);
|
||||
}
|
||||
|
||||
void GetXIntercept(const vcg::Point3i &p1, const vcg::Point3i &p2, VertexPointer &v)
|
||||
{
|
||||
vcg::Point3f p, n;
|
||||
float d;
|
||||
if (_sphere_diff.DirectedDistance(p1, p2, p, n, d))
|
||||
{
|
||||
v->P() = p;
|
||||
v->N() = n;
|
||||
}
|
||||
else
|
||||
{
|
||||
float f1 = V(p1.X(), p1.Y(), p1.Z());
|
||||
float f2 = V(p2.X(), p2.Y(), p2.Z());
|
||||
float u = (float) f1/(f1-f2);
|
||||
v->P().X() = (float) p1.X()*(1-u) + u*p2.X();
|
||||
v->P().Y() = (float) p1.Y();
|
||||
v->P().Z() = (float) p1.Z();
|
||||
|
||||
}
|
||||
}
|
||||
void GetYIntercept(const vcg::Point3i &p1, const vcg::Point3i &p2, VertexPointer &v)
|
||||
{
|
||||
vcg::Point3f p, n;
|
||||
float d;
|
||||
if (_sphere_diff.DirectedDistance(p1, p2, p, n, d))
|
||||
{
|
||||
v->P() = p;
|
||||
v->N() = n;
|
||||
}
|
||||
else
|
||||
{
|
||||
float f1 = V(p1.X(), p1.Y(), p1.Z());
|
||||
float f2 = V(p2.X(), p2.Y(), p2.Z());
|
||||
float u = (float) f1/(f1-f2);
|
||||
v->P().X() = (float) p1.X();
|
||||
v->P().Y() = (float) p1.Y()*(1-u) + u*p2.Y();
|
||||
v->P().Z() = (float) p1.Z();
|
||||
|
||||
}
|
||||
}
|
||||
void GetZIntercept(const vcg::Point3i &p1, const vcg::Point3i &p2, VertexPointer &v)
|
||||
{
|
||||
vcg::Point3f p, n;
|
||||
float d;
|
||||
if (_sphere_diff.DirectedDistance(p1, p2, p, n, d))
|
||||
{
|
||||
v->P() = p;
|
||||
v->N() = n;
|
||||
}
|
||||
else
|
||||
{
|
||||
float f1 = V(p1.X(), p1.Y(), p1.Z());
|
||||
float f2 = V(p2.X(), p2.Y(), p2.Z());
|
||||
float u = (float) f1/(f1-f2);
|
||||
v->P().X() = (float) p1.X();
|
||||
v->P().Y() = (float) p1.Y();
|
||||
v->P().Z() = (float) p1.Z()*(1-u) + u*p2.Z();
|
||||
|
||||
}
|
||||
}
|
||||
|
||||
private:
|
||||
SphereDifference _sphere_diff;
|
||||
};
|
||||
|
||||
#endif // __VCGTEST_VOLUME
|
|
@ -0,0 +1,242 @@
|
|||
#ifndef __VCGTEST_WALKER
|
||||
#define __VCGTEST_WALKER
|
||||
|
||||
#include "Definitions.h"
|
||||
#include "Volume.h"
|
||||
|
||||
// La classe Walker implementa la politica di visita del volume; conoscendo l'ordine di visita del volume
|
||||
// è conveniente che il Walker stesso si faccia carico del caching dei dati utilizzati durante l'esecuzione
|
||||
// degli algoritmi MarchingCubes ed ExtendedMarchingCubes, in particolare il calcolo del volume ai vertici
|
||||
// delle celle e delle intersezioni della superficie con le celle. In questo esempio il volume da processare
|
||||
// viene suddiviso in fette; in questo modo se il volume ha dimensione h*l*w (rispettivamente altezza,
|
||||
// larghezza e profondità), lo spazio richiesto per il caching dei vertici già allocati passa da O(h*l*w)
|
||||
// a O(h*l).
|
||||
class Walker
|
||||
{
|
||||
private:
|
||||
typedef int VertexIndex;
|
||||
|
||||
public:
|
||||
Walker(const BoundingBox &bbox, const vcg::Point3i &resolution)
|
||||
{
|
||||
_bbox = bbox;
|
||||
_resolution = resolution;
|
||||
_cell_size.X() = _bbox.DimX()/_resolution.X();
|
||||
_cell_size.Y() = _bbox.DimY()/_resolution.Y();
|
||||
_cell_size.Z() = _bbox.DimZ()/_resolution.Z();
|
||||
_slice_dimension = resolution.X()*resolution.Z();
|
||||
|
||||
_x_cs = new VertexIndex[ _slice_dimension ];
|
||||
_y_cs = new VertexIndex[ _slice_dimension ];
|
||||
_z_cs = new VertexIndex[ _slice_dimension ];
|
||||
_x_ns = new VertexIndex[ _slice_dimension ];
|
||||
_z_ns = new VertexIndex[ _slice_dimension ];
|
||||
_v_cs = new float[_slice_dimension];
|
||||
_v_ns = new float[_slice_dimension];
|
||||
|
||||
};
|
||||
|
||||
~Walker()
|
||||
{}
|
||||
|
||||
template<class EXTRACTOR_TYPE>
|
||||
void BuildMesh(Mesh &mesh, Volume &volume, EXTRACTOR_TYPE &extractor)
|
||||
{
|
||||
_volume = &volume;
|
||||
_mesh = &mesh;
|
||||
_mesh->Clear();
|
||||
vcg::Point3i p1, p2;
|
||||
|
||||
Begin();
|
||||
extractor.Initialize();
|
||||
for (int j=_bbox.min.Y(); j<_bbox.max.Y()-_cell_size.Y(); j+=_cell_size.Y())
|
||||
{
|
||||
for (int i=_bbox.min.X(); i<_bbox.max.X()-_cell_size.X(); i+=_cell_size.X())
|
||||
{
|
||||
for (int k=_bbox.min.Z(); k<_bbox.max.Z()-_cell_size.Z(); k+=_cell_size.Z())
|
||||
{
|
||||
p1.X()=i; p1.Y()=j; p1.Z()=k;
|
||||
p2.X()=i+_cell_size.X(); p2.Y()=j+_cell_size.Y(); p2.Z()=k+_cell_size.Z();
|
||||
extractor.ProcessCell(p1, p2);
|
||||
}
|
||||
}
|
||||
NextSlice();
|
||||
}
|
||||
extractor.Finalize();
|
||||
_volume = NULL;
|
||||
_mesh = NULL;
|
||||
};
|
||||
|
||||
float V(int pi, int pj, int pk)
|
||||
{
|
||||
int i = pi - _bbox.min.X();
|
||||
int k = pk - _bbox.min.Z();
|
||||
return (pj==_current_slice) ? _v_cs[i+k*_resolution.X()] : _v_ns[i+k*_resolution.X()];
|
||||
}
|
||||
|
||||
bool Exist(const vcg::Point3i &p1, const vcg::Point3i &p2, VertexPointer &v)
|
||||
{
|
||||
int i_idx = p1.X()-_bbox.min.X();
|
||||
int k_idx = p2.Z()-_bbox.min.Z();
|
||||
int index = i_idx+k_idx*_resolution.X();
|
||||
if (p1.X()!=p2.X()) //intersezione della superficie con un Xedge
|
||||
return (p1.Y()==_current_slice)? _x_cs[index]!=-1 : _x_ns[index]!=-1;
|
||||
else if (p1.Y()!=p2.Y()) //intersezione della superficie con un Yedge
|
||||
return _y_cs[index]!=-1;
|
||||
else if (p1.Z()!=p2.Z()) //intersezione della superficie con un Zedge
|
||||
return (p1.Y()==_current_slice)? _z_cs[index]!=-1 : _z_ns[index]!=-1;
|
||||
}
|
||||
|
||||
void GetXIntercept(const vcg::Point3i &p1, const vcg::Point3i &p2, VertexPointer &v)
|
||||
{
|
||||
int i = p1.X() - _bbox.min.X();
|
||||
int z = p1.Z() - _bbox.min.Z();
|
||||
VertexIndex index = i+z*_resolution.X();
|
||||
VertexIndex pos;
|
||||
if (p1.Y()==_current_slice)
|
||||
{
|
||||
if ((pos=_x_cs[index])==-1)
|
||||
{
|
||||
_x_cs[index] = (VertexIndex) _mesh->vert.size();
|
||||
pos = _x_cs[index];
|
||||
Allocator::AddVertices( *_mesh, 1 );
|
||||
v = &_mesh->vert[pos];
|
||||
_volume->GetXIntercept(p1, p2, v);
|
||||
return;
|
||||
}
|
||||
}
|
||||
if (p1.Y()==_current_slice+_cell_size.Y())
|
||||
{
|
||||
if ((pos=_x_ns[index])==-1)
|
||||
{
|
||||
_x_ns[index] = (VertexIndex) _mesh->vert.size();
|
||||
pos = _x_ns[index];
|
||||
Allocator::AddVertices( *_mesh, 1 );
|
||||
v = &_mesh->vert[pos];
|
||||
_volume->GetXIntercept(p1, p2, v);
|
||||
return;
|
||||
}
|
||||
}
|
||||
v = &_mesh->vert[pos];
|
||||
}
|
||||
void GetYIntercept(const vcg::Point3i &p1, const vcg::Point3i &p2, VertexPointer &v)
|
||||
{
|
||||
int i = p1.X() - _bbox.min.X();
|
||||
int z = p1.Z() - _bbox.min.Z();
|
||||
VertexIndex index = i+z*_resolution.X();
|
||||
VertexIndex pos;
|
||||
if ((pos=_y_cs[index])==-1)
|
||||
{
|
||||
_y_cs[index] = (VertexIndex) _mesh->vert.size();
|
||||
pos = _y_cs[index];
|
||||
Allocator::AddVertices( *_mesh, 1);
|
||||
v = &_mesh->vert[ pos ];
|
||||
_volume->GetYIntercept(p1, p2, v);
|
||||
}
|
||||
v = &_mesh->vert[pos];
|
||||
}
|
||||
void GetZIntercept(const vcg::Point3i &p1, const vcg::Point3i &p2, VertexPointer &v)
|
||||
{
|
||||
int i = p1.X() - _bbox.min.X();
|
||||
int z = p1.Z() - _bbox.min.Z();
|
||||
VertexIndex index = i+z*_resolution.X();
|
||||
VertexIndex pos;
|
||||
if (p1.Y()==_current_slice)
|
||||
{
|
||||
if ((pos=_z_cs[index])==-1)
|
||||
{
|
||||
_z_cs[index] = (VertexIndex) _mesh->vert.size();
|
||||
pos = _z_cs[index];
|
||||
Allocator::AddVertices( *_mesh, 1 );
|
||||
v = &_mesh->vert[pos];
|
||||
_volume->GetZIntercept(p1, p2, v);
|
||||
return;
|
||||
}
|
||||
}
|
||||
if (p1.Y()==_current_slice+_cell_size.Y())
|
||||
{
|
||||
if ((pos=_z_ns[index])==-1)
|
||||
{
|
||||
_z_ns[index] = (VertexIndex) _mesh->vert.size();
|
||||
pos = _z_ns[index];
|
||||
Allocator::AddVertices( *_mesh, 1 );
|
||||
v = &_mesh->vert[pos];
|
||||
_volume->GetZIntercept(p1, p2, v);
|
||||
return;
|
||||
}
|
||||
}
|
||||
v = &_mesh->vert[pos];
|
||||
}
|
||||
|
||||
protected:
|
||||
BoundingBox _bbox;
|
||||
vcg::Point3i _resolution;
|
||||
vcg::Point3i _cell_size;
|
||||
|
||||
int _slice_dimension;
|
||||
int _current_slice;
|
||||
|
||||
float *_v_cs; // il valore del campo campionato nella fetta di volumecorrente
|
||||
float *_v_ns; // il valore del campo campionato nella prossima fetta di volume
|
||||
|
||||
VertexIndex *_x_cs; // indici dell'intersezioni della superficie lungo gli Xedge della fetta corrente
|
||||
VertexIndex *_y_cs; // indici dell'intersezioni della superficie lungo gli Yedge della fetta corrente
|
||||
VertexIndex *_z_cs; // indici dell'intersezioni della superficie lungo gli Zedge della fetta corrente
|
||||
VertexIndex *_x_ns; // indici dell'intersezioni della superficie lungo gli Xedge della prossima fetta
|
||||
VertexIndex *_z_ns; // indici dell'intersezioni della superficie lungo gli Zedge della prossima fetta
|
||||
|
||||
Mesh *_mesh;
|
||||
Volume *_volume;
|
||||
|
||||
void NextSlice()
|
||||
{
|
||||
memset(_x_cs, -1, _slice_dimension*sizeof(VertexIndex));
|
||||
memset(_y_cs, -1, _slice_dimension*sizeof(VertexIndex));
|
||||
memset(_z_cs, -1, _slice_dimension*sizeof(VertexIndex));
|
||||
|
||||
std::swap(_x_cs, _x_ns);
|
||||
std::swap(_z_cs, _z_ns);
|
||||
std::swap(_v_cs, _v_ns);
|
||||
|
||||
_current_slice += _cell_size.Y();
|
||||
int j = _current_slice + _cell_size.Y();
|
||||
int k_idx, i_idx, index;
|
||||
for (int i=_bbox.min.X(); i<_bbox.max.X(); i+=_cell_size.X())
|
||||
{
|
||||
i_idx = i-_bbox.min.X();
|
||||
for (int k=_bbox.min.Z(); k<_bbox.max.Z(); k+=_cell_size.Z())
|
||||
{
|
||||
k_idx = k-_bbox.min.Z();
|
||||
index = i_idx+k_idx*_resolution.X();
|
||||
_v_ns[ index ] = _volume->V(i, j, k);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
void Begin()
|
||||
{
|
||||
_current_slice = _bbox.min.Y();
|
||||
|
||||
memset(_x_cs, -1, _slice_dimension*sizeof(VertexIndex));
|
||||
memset(_y_cs, -1, _slice_dimension*sizeof(VertexIndex));
|
||||
memset(_z_cs, -1, _slice_dimension*sizeof(VertexIndex));
|
||||
memset(_x_ns, -1, _slice_dimension*sizeof(VertexIndex));
|
||||
memset(_z_ns, -1, _slice_dimension*sizeof(VertexIndex));
|
||||
|
||||
int index;
|
||||
int j = _current_slice;
|
||||
int i_idx, k_idx;
|
||||
for (int i=_bbox.min.X(); i<_bbox.max.X(); i+=_cell_size.X())
|
||||
{
|
||||
i_idx = i-_bbox.min.X();
|
||||
for (int k=_bbox.min.Z(); k<_bbox.max.Z(); k+=_cell_size.Z())
|
||||
{
|
||||
k_idx = k-_bbox.min.Z();
|
||||
index = i_idx+k_idx*_resolution.X();
|
||||
_v_cs[index] = _volume->V(i, j, k);
|
||||
_v_ns[index] = _volume->V(i, j+_cell_size.Y(), k);
|
||||
}
|
||||
}
|
||||
}
|
||||
};
|
||||
#endif // __VCGTEST_WALKER
|
|
@ -0,0 +1,155 @@
|
|||
<?xml version="1.0" encoding="Windows-1252"?>
|
||||
<VisualStudioProject
|
||||
ProjectType="Visual C++"
|
||||
Version="7.10"
|
||||
Name="extractor"
|
||||
ProjectGUID="{53BF323F-2659-4227-9D32-6FC9624028F3}"
|
||||
Keyword="Win32Proj">
|
||||
<Platforms>
|
||||
<Platform
|
||||
Name="Win32"/>
|
||||
</Platforms>
|
||||
<Configurations>
|
||||
<Configuration
|
||||
Name="Debug|Win32"
|
||||
OutputDirectory="Debug"
|
||||
IntermediateDirectory="Debug"
|
||||
ConfigurationType="1"
|
||||
CharacterSet="2">
|
||||
<Tool
|
||||
Name="VCCLCompilerTool"
|
||||
Optimization="0"
|
||||
AdditionalIncludeDirectories="c:\Libraries\vcg\"
|
||||
PreprocessorDefinitions="WIN32;_DEBUG;_CONSOLE"
|
||||
MinimalRebuild="TRUE"
|
||||
BasicRuntimeChecks="3"
|
||||
RuntimeLibrary="5"
|
||||
UsePrecompiledHeader="0"
|
||||
WarningLevel="3"
|
||||
Detect64BitPortabilityProblems="TRUE"
|
||||
DebugInformationFormat="4"/>
|
||||
<Tool
|
||||
Name="VCCustomBuildTool"/>
|
||||
<Tool
|
||||
Name="VCLinkerTool"
|
||||
OutputFile="$(OutDir)/extractor.exe"
|
||||
LinkIncremental="2"
|
||||
GenerateDebugInformation="TRUE"
|
||||
ProgramDatabaseFile="$(OutDir)/extractor.pdb"
|
||||
SubSystem="1"
|
||||
TargetMachine="1"/>
|
||||
<Tool
|
||||
Name="VCMIDLTool"/>
|
||||
<Tool
|
||||
Name="VCPostBuildEventTool"/>
|
||||
<Tool
|
||||
Name="VCPreBuildEventTool"/>
|
||||
<Tool
|
||||
Name="VCPreLinkEventTool"/>
|
||||
<Tool
|
||||
Name="VCResourceCompilerTool"/>
|
||||
<Tool
|
||||
Name="VCWebServiceProxyGeneratorTool"/>
|
||||
<Tool
|
||||
Name="VCXMLDataGeneratorTool"/>
|
||||
<Tool
|
||||
Name="VCWebDeploymentTool"/>
|
||||
<Tool
|
||||
Name="VCManagedWrapperGeneratorTool"/>
|
||||
<Tool
|
||||
Name="VCAuxiliaryManagedWrapperGeneratorTool"/>
|
||||
</Configuration>
|
||||
<Configuration
|
||||
Name="Release|Win32"
|
||||
OutputDirectory="Release"
|
||||
IntermediateDirectory="Release"
|
||||
ConfigurationType="1"
|
||||
CharacterSet="2">
|
||||
<Tool
|
||||
Name="VCCLCompilerTool"
|
||||
AdditionalIncludeDirectories="c:\Libraries\vcg"
|
||||
PreprocessorDefinitions="WIN32;NDEBUG;_CONSOLE"
|
||||
RuntimeLibrary="4"
|
||||
UsePrecompiledHeader="0"
|
||||
WarningLevel="3"
|
||||
Detect64BitPortabilityProblems="TRUE"
|
||||
DebugInformationFormat="3"/>
|
||||
<Tool
|
||||
Name="VCCustomBuildTool"/>
|
||||
<Tool
|
||||
Name="VCLinkerTool"
|
||||
OutputFile="$(OutDir)/extractor.exe"
|
||||
LinkIncremental="1"
|
||||
GenerateDebugInformation="TRUE"
|
||||
SubSystem="1"
|
||||
OptimizeReferences="2"
|
||||
EnableCOMDATFolding="2"
|
||||
TargetMachine="1"/>
|
||||
<Tool
|
||||
Name="VCMIDLTool"/>
|
||||
<Tool
|
||||
Name="VCPostBuildEventTool"/>
|
||||
<Tool
|
||||
Name="VCPreBuildEventTool"/>
|
||||
<Tool
|
||||
Name="VCPreLinkEventTool"/>
|
||||
<Tool
|
||||
Name="VCResourceCompilerTool"/>
|
||||
<Tool
|
||||
Name="VCWebServiceProxyGeneratorTool"/>
|
||||
<Tool
|
||||
Name="VCXMLDataGeneratorTool"/>
|
||||
<Tool
|
||||
Name="VCWebDeploymentTool"/>
|
||||
<Tool
|
||||
Name="VCManagedWrapperGeneratorTool"/>
|
||||
<Tool
|
||||
Name="VCAuxiliaryManagedWrapperGeneratorTool"/>
|
||||
</Configuration>
|
||||
</Configurations>
|
||||
<References>
|
||||
</References>
|
||||
<Files>
|
||||
<Filter
|
||||
Name="Source Files"
|
||||
Filter="cpp;c;cxx;def;odl;idl;hpj;bat;asm;asmx"
|
||||
UniqueIdentifier="{4FC737F1-C7A5-4376-A066-2A32D752A2FF}">
|
||||
<File
|
||||
RelativePath=".\main.cpp">
|
||||
</File>
|
||||
<File
|
||||
RelativePath="..\..\..\..\wrap\ply\plylib.cpp">
|
||||
</File>
|
||||
</Filter>
|
||||
<Filter
|
||||
Name="Header Files"
|
||||
Filter="h;hpp;hxx;hm;inl;inc;xsd"
|
||||
UniqueIdentifier="{93995380-89BD-4b04-88EB-625FBE52EBFB}">
|
||||
<File
|
||||
RelativePath=".\Definitions.h">
|
||||
</File>
|
||||
<File
|
||||
RelativePath=".\ImplicitSphere.h">
|
||||
</File>
|
||||
<File
|
||||
RelativePath=".\SphereDifference.h">
|
||||
</File>
|
||||
<File
|
||||
RelativePath=".\SphereUnion.h">
|
||||
</File>
|
||||
<File
|
||||
RelativePath=".\Volume.h">
|
||||
</File>
|
||||
<File
|
||||
RelativePath=".\Walker.h">
|
||||
</File>
|
||||
</Filter>
|
||||
<Filter
|
||||
Name="Resource Files"
|
||||
Filter="rc;ico;cur;bmp;dlg;rc2;rct;bin;rgs;gif;jpg;jpeg;jpe;resx"
|
||||
UniqueIdentifier="{67DA6AB6-F800-4c08-8B7A-83BB121AAD01}">
|
||||
</Filter>
|
||||
</Files>
|
||||
<Globals>
|
||||
</Globals>
|
||||
</VisualStudioProject>
|
|
@ -0,0 +1,33 @@
|
|||
#include <stdio.h>
|
||||
#include <wrap/io_trimesh/export_ply.h>
|
||||
#include "Definitions.h"
|
||||
#include "Volume.h"
|
||||
#include "Walker.h"
|
||||
#include <vcg/complex/trimesh/create/marching_cubes.h>
|
||||
#include <vcg/complex/trimesh/create/extended_marching_cubes.h>
|
||||
|
||||
int main(int argc, char *argv[])
|
||||
{
|
||||
BoundingBox bbox(vcg::Point3i(-20, -20, -20), vcg::Point3i(20, 20, 20));
|
||||
vcg::Point3i resolution(40, 40, 40);
|
||||
|
||||
Volume volume;
|
||||
Walker walker(bbox, resolution);
|
||||
|
||||
typedef vcg::tri::MarchingCubes<Mesh, Walker> MarchingCubes;
|
||||
typedef vcg::tri::ExtendedMarchingCubes<Mesh, Walker> ExtendedMarchingCubes;
|
||||
|
||||
|
||||
|
||||
// MARCHING CUBES
|
||||
Mesh mc_mesh;
|
||||
MarchingCubes mc(mc_mesh, walker);
|
||||
walker.BuildMesh<MarchingCubes>(mc_mesh, volume, mc);
|
||||
vcg::tri::io::ExporterPLY<Mesh>::Save( mc_mesh, "marching_cubes.ply");
|
||||
|
||||
// EXTENDED MARCHING CUBES
|
||||
Mesh emc_mesh;
|
||||
ExtendedMarchingCubes emc(emc_mesh, walker, 30);
|
||||
walker.BuildMesh<ExtendedMarchingCubes>(emc_mesh, volume, emc);
|
||||
vcg::tri::io::ExporterPLY<Mesh>::Save( emc_mesh, "extended_marching_cubes.ply");
|
||||
};
|
Binary file not shown.
|
@ -0,0 +1,21 @@
|
|||
Microsoft Visual Studio Solution File, Format Version 8.00
|
||||
Project("{8BC9CEB8-8B4A-11D0-8D11-00A0C91BC942}") = "extractor", "extractor\extractor.vcproj", "{53BF323F-2659-4227-9D32-6FC9624028F3}"
|
||||
ProjectSection(ProjectDependencies) = postProject
|
||||
EndProjectSection
|
||||
EndProject
|
||||
Global
|
||||
GlobalSection(SolutionConfiguration) = preSolution
|
||||
Debug = Debug
|
||||
Release = Release
|
||||
EndGlobalSection
|
||||
GlobalSection(ProjectConfiguration) = postSolution
|
||||
{53BF323F-2659-4227-9D32-6FC9624028F3}.Debug.ActiveCfg = Debug|Win32
|
||||
{53BF323F-2659-4227-9D32-6FC9624028F3}.Debug.Build.0 = Debug|Win32
|
||||
{53BF323F-2659-4227-9D32-6FC9624028F3}.Release.ActiveCfg = Release|Win32
|
||||
{53BF323F-2659-4227-9D32-6FC9624028F3}.Release.Build.0 = Release|Win32
|
||||
EndGlobalSection
|
||||
GlobalSection(ExtensibilityGlobals) = postSolution
|
||||
EndGlobalSection
|
||||
GlobalSection(ExtensibilityAddIns) = postSolution
|
||||
EndGlobalSection
|
||||
EndGlobal
|
Binary file not shown.
Loading…
Reference in New Issue