28 #include <Eigen/Eigenvalues>
29 #include <vcg/complex/algorithms/update/normal.h>
53 template <
class MeshType>
56 typedef typename MeshType::VertexType VertexType;
57 typedef typename MeshType::VertexPointer VertexPointer;
58 typedef typename MeshType::VertexIterator VertexIterator;
59 typedef typename MeshType::ScalarType ScalarType;
60 typedef typename MeshType::FaceType FaceType;
61 typedef typename MeshType::FacePointer FacePointer;
62 typedef typename MeshType::FaceIterator FaceIterator;
63 typedef typename MeshType::ConstFaceIterator ConstFaceIterator;
64 typedef typename MeshType::FaceContainer FaceContainer;
65 typedef typename MeshType::CoordType CoordType;
69 inline ScalarType SQR(
const ScalarType &x)
const {
return x*x;}
70 inline ScalarType CUBE(
const ScalarType &x)
const {
return x*x*x;}
77 double P1, Pa, Pb, Paa, Pab, Pbb, Paaa, Paab, Pabb, Pbbb;
80 double Fa, Fb, Fc, Faa, Fbb, Fcc, Faaa, Fbbb, Fccc, Faab, Fbbc, Fcca;
83 double T0, T1[3], T2[3], TP[3];
94 void compProjectionIntegrals(
const FaceType &f)
98 double a0_2, a0_3, a0_4, b0_2, b0_3, b0_4;
99 double a1_2, a1_3, b1_2, b1_3;
100 double C1, Ca, Caa, Caaa, Cb, Cbb, Cbbb;
101 double Cab, Kab, Caab, Kaab, Cabb, Kabb;
104 P1 = Pa = Pb = Paa = Pab = Pbb = Paaa = Paab = Pabb = Pbbb = 0.0;
106 for (i = 0; i < 3; i++) {
109 a1 = f.V1(i)->P()[A];
110 b1 = f.V1(i)->P()[B];
113 a0_2 = a0 * a0; a0_3 = a0_2 * a0; a0_4 = a0_3 * a0;
114 b0_2 = b0 * b0; b0_3 = b0_2 * b0; b0_4 = b0_3 * b0;
115 a1_2 = a1 * a1; a1_3 = a1_2 * a1;
116 b1_2 = b1 * b1; b1_3 = b1_2 * b1;
119 Ca = a1*C1 + a0_2; Caa = a1*Ca + a0_3; Caaa = a1*Caa + a0_4;
120 Cb = b1*(b1 + b0) + b0_2; Cbb = b1*Cb + b0_3; Cbbb = b1*Cbb + b0_4;
121 Cab = 3*a1_2 + 2*a1*a0 + a0_2; Kab = a1_2 + 2*a1*a0 + 3*a0_2;
122 Caab = a0*Cab + 4*a1_3; Kaab = a1*Kab + 4*a0_3;
123 Cabb = 4*b1_3 + 3*b1_2*b0 + 2*b1*b0_2 + b0_3;
124 Kabb = b1_3 + 2*b1_2*b0 + 3*b1*b0_2 + 4*b0_3;
133 Pab += db*(b1*Cab + b0*Kab);
134 Paab += db*(b1*Caab + b0*Kaab);
135 Pabb += da*(a1*Cabb + a0*Kabb);
151 void CompFaceIntegrals(
const FaceType &f,
const Point3<ScalarType> &n)
154 double k1, k2, k3, k4;
156 compProjectionIntegrals(f);
159 k1 = 1 / n[C]; k2 = k1 * k1; k3 = k2 * k1; k4 = k3 * k1;
163 Fc = -k2 * (n[A]*Pa + n[B]*Pb + w*P1);
167 Fcc = k3 * (SQR(n[A])*Paa + 2*n[A]*n[B]*Pab + SQR(n[B])*Pbb
168 + w*(2*(n[A]*Pa + n[B]*Pb) + w*P1));
172 Fccc = -k4 * (CUBE(n[A])*Paaa + 3*SQR(n[A])*n[B]*Paab
173 + 3*n[A]*SQR(n[B])*Pabb + CUBE(n[B])*Pbbb
174 + 3*w*(SQR(n[A])*Paa + 2*n[A]*n[B]*Pab + SQR(n[B])*Pbb)
175 + w*w*(3*(n[A]*Pa + n[B]*Pb) + w*P1));
178 Fbbc = -k2 * (n[A]*Pabb + n[B]*Pbbb + w*Pbb);
179 Fcca = k3 * (SQR(n[A])*Paaa + 2*n[A]*n[B]*Paab + SQR(n[B])*Pabb
180 + w*(2*(n[A]*Paa + n[B]*Pab) + w*Pa));
193 T0 = T1[X] = T1[Y] = T1[Z]
194 = T2[X] = T2[Y] = T2[Z]
195 = TP[X] = TP[Y] = TP[Z] = 0;
196 for (
auto fi=m.face.begin(); fi!=m.face.end();++fi)
if(!(*fi).IsD() && vcg::DoubleArea(*fi)>std::numeric_limits<float>::min())
198 const FaceType &f=(*fi);
199 const auto fn = vcg::NormalizedTriangleNormal(f);
204 if (nx > ny && nx > nz) C = X;
205 else C = (ny > nz) ? Y : Z;
209 CompFaceIntegrals(f, fn);
211 T0 += fn[X] * ((A == X) ? Fa : ((B == X) ? Fb : Fc));
213 T1[A] += fn[A] * Faa;
214 T1[B] += fn[B] * Fbb;
215 T1[C] += fn[C] * Fcc;
216 T2[A] += fn[A] * Faaa;
217 T2[B] += fn[B] * Fbbb;
218 T2[C] += fn[C] * Fccc;
219 TP[A] += fn[A] * Faab;
220 TP[B] += fn[B] * Fbbc;
221 TP[C] += fn[C] * Fcca;
224 T1[X] /= 2; T1[Y] /= 2; T1[Z] /= 2;
225 T2[X] /= 3; T2[Y] /= 3; T2[Z] /= 3;
226 TP[X] /= 2; TP[Y] /= 2; TP[Z] /= 2;
235 return static_cast<ScalarType
>(T0);
251 void InertiaTensor(Matrix33<ScalarType> &J)
const
258 J[X][X] = (T2[Y] + T2[Z]);
259 J[Y][Y] = (T2[Z] + T2[X]);
260 J[Z][Z] = (T2[X] + T2[Y]);
261 J[X][Y] = J[Y][X] = - TP[X];
262 J[Y][Z] = J[Z][Y] = - TP[Y];
263 J[Z][X] = J[X][Z] = - TP[Z];
265 J[X][X] -= T0 * (r[Y]*r[Y] + r[Z]*r[Z]);
266 J[Y][Y] -= T0 * (r[Z]*r[Z] + r[X]*r[X]);
267 J[Z][Z] -= T0 * (r[X]*r[X] + r[Y]*r[Y]);
268 J[X][Y] = J[Y][X] += T0 * r[X] * r[Y];
269 J[Y][Z] = J[Z][Y] += T0 * r[Y] * r[Z];
270 J[Z][X] = J[X][Z] += T0 * r[Z] * r[X];
274 void InertiaTensor(Eigen::Matrix3d &J)
const
276 J=Eigen::Matrix3d::Identity();
282 J(X,X) = (T2[Y] + T2[Z]);
283 J(Y,Y) = (T2[Z] + T2[X]);
284 J(Z,Z) = (T2[X] + T2[Y]);
285 J(X,Y) = J(Y,X) = - TP[X];
286 J(Y,Z) = J(Z,Y) = - TP[Y];
287 J(Z,X) = J(X,Z) = - TP[Z];
289 J(X,X) -= T0 * (r[Y]*r[Y] + r[Z]*r[Z]);
290 J(Y,Y) -= T0 * (r[Z]*r[Z] + r[X]*r[X]);
291 J(Z,Z) -= T0 * (r[X]*r[X] + r[Y]*r[Y]);
292 J(X,Y) = J(Y,X) += T0 * r[X] * r[Y];
293 J(Y,Z) = J(Z,Y) += T0 * r[Y] * r[Z];
294 J(Z,X) = J(X,Z) += T0 * r[Z] * r[X];
307 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> eig(it);
308 Eigen::Vector3d c_val = eig.eigenvalues();
309 Eigen::Matrix3d c_vec = eig.eigenvectors();
310 EV.FromEigenMatrix(c_vec);
311 EV.transposeInPlace();
312 ev.FromEigenVector(c_val);
321 ConstFaceIterator fi;
322 ScalarType area = 0.0;
324 for(fi = m.face.begin(); fi != m.face.end(); ++fi)
327 bary += vcg::Barycenter( *fi )* vcg::DoubleArea(*fi);
328 area+=vcg::DoubleArea(*fi);
335 vcg::Matrix33<ScalarType> C0;
337 C0[0][0] = C0[1][1] = 2.0;
338 C0[0][1] = C0[1][0] = 1.0;
342 CoordType X(1/6.0,1/6.0,0);
343 vcg::Matrix33<ScalarType> A,
345 for(fi = m.face.begin(); fi != m.face.end(); ++fi)
348 const CoordType &P0 = (*fi).cP(0);
349 const CoordType &P1 = (*fi).cP(1);
350 const CoordType &P2 = (*fi).cP(2);
351 CoordType n = ((P1-P0)^(P2-P0));
352 const float da = n.Norm();
355 A.SetColumn(0, P1-P0);
356 A.SetColumn(1, P2-P0);
358 CoordType delta = P0 - bary;
365 DC+= A*C0*A.transpose();
366 vcg::Matrix33<ScalarType> tmp;
367 tmp.OuterProduct(A*X,delta);
368 DC += tmp + tmp.transpose();
370 tmp.OuterProduct(delta,delta);
Methods for computing Polyhedral Mass properties (like inertia tensor, volume, etc)
Definition: inertia.h:55
Point3< ScalarType > CenterOfMass(void) const
Return the Center of Mass (or barycenter) of the mesh.
Definition: inertia.h:242
Inertia(const MeshType &m)
Basic constructor.
Definition: inertia.h:91
void InertiaTensorEigen(Matrix33< ScalarType > &EV, Point3< ScalarType > &ev) const
Return the Inertia tensor the mesh.
Definition: inertia.h:303
void Compute(const MeshType &m)
Definition: inertia.h:189
static void Covariance(const MeshType &m, vcg::Point3< ScalarType > &bary, vcg::Matrix33< ScalarType > &C)
Definition: inertia.h:318
ScalarType Mass(void) const
Return the Volume (or mass) of the mesh.
Definition: inertia.h:233
Definition: namespaces.dox:6