added Covariance to to compute the covariance of a generic mesh (without the use of divergence theorem)

This commit is contained in:
ganovelli 2008-08-19 09:43:24 +00:00
parent 1ac5c66e78
commit 4e81e65145
1 changed files with 378 additions and 308 deletions

View File

@ -1,308 +1,378 @@
/**************************************************************************** /****************************************************************************
* VCGLib o o * * VCGLib o o *
* Visual and Computer Graphics Library o o * * Visual and Computer Graphics Library o o *
* _ O _ * * _ O _ *
* Copyright(C) 2005 \/)\/ * * Copyright(C) 2005 \/)\/ *
* Visual Computing Lab /\/| * * Visual Computing Lab /\/| *
* ISTI - Italian National Research Council | * * ISTI - Italian National Research Council | *
* \ * * \ *
* All rights reserved. * * All rights reserved. *
* * * *
* This program is free software; you can redistribute it and/or modify * * This program is free software; you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by * * it under the terms of the GNU General Public License as published by *
* the Free Software Foundation; either version 2 of the License, or * * the Free Software Foundation; either version 2 of the License, or *
* (at your option) any later version. * * (at your option) any later version. *
* * * *
* This program is distributed in the hope that it will be useful, * * This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of * * but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License (http://www.gnu.org/licenses/gpl.txt) * * GNU General Public License (http://www.gnu.org/licenses/gpl.txt) *
* for more details. * * for more details. *
* * * *
****************************************************************************/ ****************************************************************************/
/**************************************************************************** /****************************************************************************
History History
$Log: not supported by cvs2svn $ $Log: not supported by cvs2svn $
Revision 1.3 2006/03/29 10:12:08 corsini Revision 1.3 2006/03/29 10:12:08 corsini
Add cast to avoid warning Add cast to avoid warning
Revision 1.2 2005/12/12 12:08:30 cignoni Revision 1.2 2005/12/12 12:08:30 cignoni
First working version First working version
Revision 1.1 2005/11/21 15:58:12 cignoni Revision 1.1 2005/11/21 15:58:12 cignoni
First Release (not working!) First Release (not working!)
Revision 1.13 2005/11/17 00:42:03 cignoni Revision 1.13 2005/11/17 00:42:03 cignoni
****************************************************************************/ ****************************************************************************/
/* #ifndef _VCG_INERTIA_
The algorithm is based on a three step reduction of the volume integrals #define _VCG_INERTIA_
to successively simpler integrals. The algorithm is designed to minimize
the numerical errors that can result from poorly conditioned alignment of /*
polyhedral faces. It is also designed for efficiency. All required volume The algorithm is based on a three step reduction of the volume integrals
integrals of a polyhedron are computed together during a single walk over to successively simpler integrals. The algorithm is designed to minimize
the boundary of the polyhedron; exploiting common subexpressions reduces the numerical errors that can result from poorly conditioned alignment of
floating point operations. polyhedral faces. It is also designed for efficiency. All required volume
integrals of a polyhedron are computed together during a single walk over
For more information, check out: the boundary of the polyhedron; exploiting common subexpressions reduces
floating point operations.
Brian Mirtich,
``Fast and Accurate Computation of Polyhedral Mass Properties,'' For more information, check out:
journal of graphics tools, volume 1, number 2, 1996
Brian Mirtich,
*/ ``Fast and Accurate Computation of Polyhedral Mass Properties,''
journal of graphics tools, volume 1, number 2, 1996
#include <vcg/math/matrix33.h>
#include <vcg/math/lin_algebra.h> */
#include <vcg/complex/trimesh/update/normal.h> #include <vcg/math/matrix33.h>
namespace vcg #include <vcg/math/lin_algebra.h>
{
namespace tri #include <vcg/complex/trimesh/update/normal.h>
{ namespace vcg
template <class InertiaMeshType> {
class Inertia namespace tri
{ {
typedef InertiaMeshType MeshType; template <class InertiaMeshType>
typedef typename MeshType::VertexType VertexType; class Inertia
typedef typename MeshType::VertexPointer VertexPointer; {
typedef typename MeshType::VertexIterator VertexIterator; typedef InertiaMeshType MeshType;
typedef typename MeshType::ScalarType ScalarType; typedef typename MeshType::VertexType VertexType;
typedef typename MeshType::FaceType FaceType; typedef typename MeshType::VertexPointer VertexPointer;
typedef typename MeshType::FacePointer FacePointer; typedef typename MeshType::VertexIterator VertexIterator;
typedef typename MeshType::FaceIterator FaceIterator; typedef typename MeshType::ScalarType ScalarType;
typedef typename MeshType::FaceContainer FaceContainer; typedef typename MeshType::FaceType FaceType;
typedef typename MeshType::FacePointer FacePointer;
private : typedef typename MeshType::FaceIterator FaceIterator;
enum {X=0,Y=1,Z=2}; typedef typename MeshType::FaceContainer FaceContainer;
inline ScalarType SQR(ScalarType &x) const { return x*x;} typedef typename MeshType::CoordType CoordType;
inline ScalarType CUBE(ScalarType &x) const { return x*x*x;} typedef typename MeshType::ScalarType ScalarType;
int A; /* alpha */ private :
int B; /* beta */ enum {X=0,Y=1,Z=2};
int C; /* gamma */ inline ScalarType SQR(ScalarType &x) const { return x*x;}
inline ScalarType CUBE(ScalarType &x) const { return x*x*x;}
/* projection integrals */
double P1, Pa, Pb, Paa, Pab, Pbb, Paaa, Paab, Pabb, Pbbb; int A; /* alpha */
int B; /* beta */
/* face integrals */ int C; /* gamma */
double Fa, Fb, Fc, Faa, Fbb, Fcc, Faaa, Fbbb, Fccc, Faab, Fbbc, Fcca;
/* projection integrals */
/* volume integrals */ double P1, Pa, Pb, Paa, Pab, Pbb, Paaa, Paab, Pabb, Pbbb;
double T0, T1[3], T2[3], TP[3];
/* face integrals */
public: double Fa, Fb, Fc, Faa, Fbb, Fcc, Faaa, Fbbb, Fccc, Faab, Fbbc, Fcca;
/* compute various integrations over projection of face */ /* volume integrals */
void compProjectionIntegrals(FaceType &f) double T0, T1[3], T2[3], TP[3];
{
double a0, a1, da; public:
double b0, b1, db;
double a0_2, a0_3, a0_4, b0_2, b0_3, b0_4; /* compute various integrations over projection of face */
double a1_2, a1_3, b1_2, b1_3; void compProjectionIntegrals(FaceType &f)
double C1, Ca, Caa, Caaa, Cb, Cbb, Cbbb; {
double Cab, Kab, Caab, Kaab, Cabb, Kabb; double a0, a1, da;
int i; double b0, b1, db;
double a0_2, a0_3, a0_4, b0_2, b0_3, b0_4;
P1 = Pa = Pb = Paa = Pab = Pbb = Paaa = Paab = Pabb = Pbbb = 0.0; double a1_2, a1_3, b1_2, b1_3;
double C1, Ca, Caa, Caaa, Cb, Cbb, Cbbb;
for (i = 0; i < 3; i++) { double Cab, Kab, Caab, Kaab, Cabb, Kabb;
a0 = f.V(i)->P()[A]; int i;
b0 = f.V(i)->P()[B];
a1 = f.V1(i)->P()[A]; P1 = Pa = Pb = Paa = Pab = Pbb = Paaa = Paab = Pabb = Pbbb = 0.0;
b1 = f.V1(i)->P()[B];
da = a1 - a0; for (i = 0; i < 3; i++) {
db = b1 - b0; a0 = f.V(i)->P()[A];
a0_2 = a0 * a0; a0_3 = a0_2 * a0; a0_4 = a0_3 * a0; b0 = f.V(i)->P()[B];
b0_2 = b0 * b0; b0_3 = b0_2 * b0; b0_4 = b0_3 * b0; a1 = f.V1(i)->P()[A];
a1_2 = a1 * a1; a1_3 = a1_2 * a1; b1 = f.V1(i)->P()[B];
b1_2 = b1 * b1; b1_3 = b1_2 * b1; da = a1 - a0;
db = b1 - b0;
C1 = a1 + a0; a0_2 = a0 * a0; a0_3 = a0_2 * a0; a0_4 = a0_3 * a0;
Ca = a1*C1 + a0_2; Caa = a1*Ca + a0_3; Caaa = a1*Caa + a0_4; b0_2 = b0 * b0; b0_3 = b0_2 * b0; b0_4 = b0_3 * b0;
Cb = b1*(b1 + b0) + b0_2; Cbb = b1*Cb + b0_3; Cbbb = b1*Cbb + b0_4; a1_2 = a1 * a1; a1_3 = a1_2 * a1;
Cab = 3*a1_2 + 2*a1*a0 + a0_2; Kab = a1_2 + 2*a1*a0 + 3*a0_2; b1_2 = b1 * b1; b1_3 = b1_2 * b1;
Caab = a0*Cab + 4*a1_3; Kaab = a1*Kab + 4*a0_3;
Cabb = 4*b1_3 + 3*b1_2*b0 + 2*b1*b0_2 + b0_3; C1 = a1 + a0;
Kabb = b1_3 + 2*b1_2*b0 + 3*b1*b0_2 + 4*b0_3; Ca = a1*C1 + a0_2; Caa = a1*Ca + a0_3; Caaa = a1*Caa + a0_4;
Cb = b1*(b1 + b0) + b0_2; Cbb = b1*Cb + b0_3; Cbbb = b1*Cbb + b0_4;
P1 += db*C1; Cab = 3*a1_2 + 2*a1*a0 + a0_2; Kab = a1_2 + 2*a1*a0 + 3*a0_2;
Pa += db*Ca; Caab = a0*Cab + 4*a1_3; Kaab = a1*Kab + 4*a0_3;
Paa += db*Caa; Cabb = 4*b1_3 + 3*b1_2*b0 + 2*b1*b0_2 + b0_3;
Paaa += db*Caaa; Kabb = b1_3 + 2*b1_2*b0 + 3*b1*b0_2 + 4*b0_3;
Pb += da*Cb;
Pbb += da*Cbb; P1 += db*C1;
Pbbb += da*Cbbb; Pa += db*Ca;
Pab += db*(b1*Cab + b0*Kab); Paa += db*Caa;
Paab += db*(b1*Caab + b0*Kaab); Paaa += db*Caaa;
Pabb += da*(a1*Cabb + a0*Kabb); Pb += da*Cb;
} Pbb += da*Cbb;
Pbbb += da*Cbbb;
P1 /= 2.0; Pab += db*(b1*Cab + b0*Kab);
Pa /= 6.0; Paab += db*(b1*Caab + b0*Kaab);
Paa /= 12.0; Pabb += da*(a1*Cabb + a0*Kabb);
Paaa /= 20.0; }
Pb /= -6.0;
Pbb /= -12.0; P1 /= 2.0;
Pbbb /= -20.0; Pa /= 6.0;
Pab /= 24.0; Paa /= 12.0;
Paab /= 60.0; Paaa /= 20.0;
Pabb /= -60.0; Pb /= -6.0;
} Pbb /= -12.0;
Pbbb /= -20.0;
Pab /= 24.0;
void CompFaceIntegrals(FaceType &f) Paab /= 60.0;
{ Pabb /= -60.0;
Point3<ScalarType> n; }
ScalarType w;
double k1, k2, k3, k4;
void CompFaceIntegrals(FaceType &f)
compProjectionIntegrals(f); {
Point3<ScalarType> n;
n = f.N(); ScalarType w;
w = -f.V(0)->P()*n; double k1, k2, k3, k4;
k1 = 1 / n[C]; k2 = k1 * k1; k3 = k2 * k1; k4 = k3 * k1;
compProjectionIntegrals(f);
Fa = k1 * Pa;
Fb = k1 * Pb; n = f.N();
Fc = -k2 * (n[A]*Pa + n[B]*Pb + w*P1); w = -f.V(0)->P()*n;
k1 = 1 / n[C]; k2 = k1 * k1; k3 = k2 * k1; k4 = k3 * k1;
Faa = k1 * Paa;
Fbb = k1 * Pbb; Fa = k1 * Pa;
Fcc = k3 * (SQR(n[A])*Paa + 2*n[A]*n[B]*Pab + SQR(n[B])*Pbb Fb = k1 * Pb;
+ w*(2*(n[A]*Pa + n[B]*Pb) + w*P1)); Fc = -k2 * (n[A]*Pa + n[B]*Pb + w*P1);
Faaa = k1 * Paaa; Faa = k1 * Paa;
Fbbb = k1 * Pbbb; Fbb = k1 * Pbb;
Fccc = -k4 * (CUBE(n[A])*Paaa + 3*SQR(n[A])*n[B]*Paab Fcc = k3 * (SQR(n[A])*Paa + 2*n[A]*n[B]*Pab + SQR(n[B])*Pbb
+ 3*n[A]*SQR(n[B])*Pabb + CUBE(n[B])*Pbbb + w*(2*(n[A]*Pa + n[B]*Pb) + w*P1));
+ 3*w*(SQR(n[A])*Paa + 2*n[A]*n[B]*Pab + SQR(n[B])*Pbb)
+ w*w*(3*(n[A]*Pa + n[B]*Pb) + w*P1)); Faaa = k1 * Paaa;
Fbbb = k1 * Pbbb;
Faab = k1 * Paab; Fccc = -k4 * (CUBE(n[A])*Paaa + 3*SQR(n[A])*n[B]*Paab
Fbbc = -k2 * (n[A]*Pabb + n[B]*Pbbb + w*Pbb); + 3*n[A]*SQR(n[B])*Pabb + CUBE(n[B])*Pbbb
Fcca = k3 * (SQR(n[A])*Paaa + 2*n[A]*n[B]*Paab + SQR(n[B])*Pabb + 3*w*(SQR(n[A])*Paa + 2*n[A]*n[B]*Pab + SQR(n[B])*Pbb)
+ w*(2*(n[A]*Paa + n[B]*Pab) + w*Pa)); + w*w*(3*(n[A]*Pa + n[B]*Pb) + w*P1));
}
Faab = k1 * Paab;
Fbbc = -k2 * (n[A]*Pabb + n[B]*Pbbb + w*Pbb);
void Compute(MeshType &m) Fcca = k3 * (SQR(n[A])*Paaa + 2*n[A]*n[B]*Paab + SQR(n[B])*Pabb
{ + w*(2*(n[A]*Paa + n[B]*Pab) + w*Pa));
tri::UpdateNormals<MeshType>::PerFaceNormalized(m); }
double nx, ny, nz;
T0 = T1[X] = T1[Y] = T1[Z] void Compute(MeshType &m)
= T2[X] = T2[Y] = T2[Z] {
= TP[X] = TP[Y] = TP[Z] = 0; tri::UpdateNormals<MeshType>::PerFaceNormalized(m);
FaceIterator fi; double nx, ny, nz;
for (fi=m.face.begin(); fi!=m.face.end();++fi) if(!(*fi).IsD()) {
FaceType &f=(*fi); T0 = T1[X] = T1[Y] = T1[Z]
= T2[X] = T2[Y] = T2[Z]
nx = fabs(f.N()[0]); = TP[X] = TP[Y] = TP[Z] = 0;
ny = fabs(f.N()[1]); FaceIterator fi;
nz = fabs(f.N()[2]); for (fi=m.face.begin(); fi!=m.face.end();++fi) if(!(*fi).IsD()) {
if (nx > ny && nx > nz) C = X; FaceType &f=(*fi);
else C = (ny > nz) ? Y : Z;
A = (C + 1) % 3; nx = fabs(f.N()[0]);
B = (A + 1) % 3; ny = fabs(f.N()[1]);
nz = fabs(f.N()[2]);
CompFaceIntegrals(f); if (nx > ny && nx > nz) C = X;
else C = (ny > nz) ? Y : Z;
T0 += f.N()[X] * ((A == X) ? Fa : ((B == X) ? Fb : Fc)); A = (C + 1) % 3;
B = (A + 1) % 3;
T1[A] += f.N()[A] * Faa;
T1[B] += f.N()[B] * Fbb; CompFaceIntegrals(f);
T1[C] += f.N()[C] * Fcc;
T2[A] += f.N()[A] * Faaa; T0 += f.N()[X] * ((A == X) ? Fa : ((B == X) ? Fb : Fc));
T2[B] += f.N()[B] * Fbbb;
T2[C] += f.N()[C] * Fccc; T1[A] += f.N()[A] * Faa;
TP[A] += f.N()[A] * Faab; T1[B] += f.N()[B] * Fbb;
TP[B] += f.N()[B] * Fbbc; T1[C] += f.N()[C] * Fcc;
TP[C] += f.N()[C] * Fcca; T2[A] += f.N()[A] * Faaa;
} T2[B] += f.N()[B] * Fbbb;
T2[C] += f.N()[C] * Fccc;
T1[X] /= 2; T1[Y] /= 2; T1[Z] /= 2; TP[A] += f.N()[A] * Faab;
T2[X] /= 3; T2[Y] /= 3; T2[Z] /= 3; TP[B] += f.N()[B] * Fbbc;
TP[X] /= 2; TP[Y] /= 2; TP[Z] /= 2; TP[C] += f.N()[C] * Fcca;
} }
ScalarType Mass() T1[X] /= 2; T1[Y] /= 2; T1[Z] /= 2;
{ T2[X] /= 3; T2[Y] /= 3; T2[Z] /= 3;
return static_cast<ScalarType>(T0); TP[X] /= 2; TP[Y] /= 2; TP[Z] /= 2;
} }
Point3<ScalarType> CenterOfMass() ScalarType Mass()
{ {
Point3<ScalarType> r; return static_cast<ScalarType>(T0);
r[X] = T1[X] / T0; }
r[Y] = T1[Y] / T0;
r[Z] = T1[Z] / T0; Point3<ScalarType> CenterOfMass()
return r; {
} Point3<ScalarType> r;
void InertiaTensor(Matrix33<ScalarType> &J ){ r[X] = T1[X] / T0;
Point3<ScalarType> r; r[Y] = T1[Y] / T0;
r[X] = T1[X] / T0; r[Z] = T1[Z] / T0;
r[Y] = T1[Y] / T0; return r;
r[Z] = T1[Z] / T0; }
/* compute inertia tensor */ void InertiaTensor(Matrix33<ScalarType> &J ){
J[X][X] = (T2[Y] + T2[Z]); Point3<ScalarType> r;
J[Y][Y] = (T2[Z] + T2[X]); r[X] = T1[X] / T0;
J[Z][Z] = (T2[X] + T2[Y]); r[Y] = T1[Y] / T0;
J[X][Y] = J[Y][X] = - TP[X]; r[Z] = T1[Z] / T0;
J[Y][Z] = J[Z][Y] = - TP[Y]; /* compute inertia tensor */
J[Z][X] = J[X][Z] = - TP[Z]; J[X][X] = (T2[Y] + T2[Z]);
J[Y][Y] = (T2[Z] + T2[X]);
J[X][X] -= T0 * (r[Y]*r[Y] + r[Z]*r[Z]); J[Z][Z] = (T2[X] + T2[Y]);
J[Y][Y] -= T0 * (r[Z]*r[Z] + r[X]*r[X]); J[X][Y] = J[Y][X] = - TP[X];
J[Z][Z] -= T0 * (r[X]*r[X] + r[Y]*r[Y]); J[Y][Z] = J[Z][Y] = - TP[Y];
J[X][Y] = J[Y][X] += T0 * r[X] * r[Y]; J[Z][X] = J[X][Z] = - TP[Z];
J[Y][Z] = J[Z][Y] += T0 * r[Y] * r[Z];
J[Z][X] = J[X][Z] += T0 * r[Z] * r[X]; J[X][X] -= T0 * (r[Y]*r[Y] + r[Z]*r[Z]);
} J[Y][Y] -= T0 * (r[Z]*r[Z] + r[X]*r[X]);
J[Z][Z] -= T0 * (r[X]*r[X] + r[Y]*r[Y]);
void InertiaTensor(Matrix44<ScalarType> &J ) J[X][Y] = J[Y][X] += T0 * r[X] * r[Y];
{ J[Y][Z] = J[Z][Y] += T0 * r[Y] * r[Z];
J.SetIdentity(); J[Z][X] = J[X][Z] += T0 * r[Z] * r[X];
Point3<ScalarType> r; }
r[X] = T1[X] / T0;
r[Y] = T1[Y] / T0; void InertiaTensor(Matrix44<ScalarType> &J )
r[Z] = T1[Z] / T0; {
/* compute inertia tensor */ J.SetIdentity();
J[X][X] = (T2[Y] + T2[Z]); Point3<ScalarType> r;
J[Y][Y] = (T2[Z] + T2[X]); r[X] = T1[X] / T0;
J[Z][Z] = (T2[X] + T2[Y]); r[Y] = T1[Y] / T0;
J[X][Y] = J[Y][X] = - TP[X]; r[Z] = T1[Z] / T0;
J[Y][Z] = J[Z][Y] = - TP[Y]; /* compute inertia tensor */
J[Z][X] = J[X][Z] = - TP[Z]; J[X][X] = (T2[Y] + T2[Z]);
J[Y][Y] = (T2[Z] + T2[X]);
J[X][X] -= T0 * (r[Y]*r[Y] + r[Z]*r[Z]); J[Z][Z] = (T2[X] + T2[Y]);
J[Y][Y] -= T0 * (r[Z]*r[Z] + r[X]*r[X]); J[X][Y] = J[Y][X] = - TP[X];
J[Z][Z] -= T0 * (r[X]*r[X] + r[Y]*r[Y]); J[Y][Z] = J[Z][Y] = - TP[Y];
J[X][Y] = J[Y][X] += T0 * r[X] * r[Y]; J[Z][X] = J[X][Z] = - TP[Z];
J[Y][Z] = J[Z][Y] += T0 * r[Y] * r[Z];
J[Z][X] = J[X][Z] += T0 * r[Z] * r[X]; J[X][X] -= T0 * (r[Y]*r[Y] + r[Z]*r[Z]);
} J[Y][Y] -= T0 * (r[Z]*r[Z] + r[X]*r[X]);
J[Z][Z] -= T0 * (r[X]*r[X] + r[Y]*r[Y]);
J[X][Y] = J[Y][X] += T0 * r[X] * r[Y];
// Calcola autovalori ed autovettori dell'inertia tensor. J[Y][Z] = J[Z][Y] += T0 * r[Y] * r[Z];
// Gli autovettori fanno una rotmatrix che se applicata mette l'oggetto secondo gli assi id minima/max inerzia. J[Z][X] = J[X][Z] += T0 * r[Z] * r[X];
}
void InertiaTensorEigen(Matrix44<ScalarType> &EV, Point4<ScalarType> &ev )
{
Matrix44<ScalarType> it; /** Compute eigenvalues and eigenvectors of inertia tensor.
InertiaTensor(it); The eigenvectors make a rotation matrix that aligns the mesh along the axes of min/max inertia
Matrix44d EVd,ITd;ITd.Import(it); */
Point4d evd; evd.Import(ev); void InertiaTensorEigen(Matrix44<ScalarType> &EV, Point4<ScalarType> &ev )
int n; {
Jacobi(ITd,evd,EVd,n); Matrix44<ScalarType> it;
EV.Import(EVd); InertiaTensor(it);
ev.Import(evd); Matrix44d EVd,ITd;ITd.Import(it);
} Point4d evd; evd.Import(ev);
int n;
}; // end class Inertia Jacobi(ITd,evd,EVd,n);
EV.Import(EVd);
} // end namespace tri ev.Import(evd);
} // end namespace vcg }
/** Compute covariance matrix of a mesh, i.e. the integral
int_{M} { (x-b)(x-b)^T }dx where b is the barycenter and x spans over the mesh M
*/
static void Covariance(MeshType m, vcg::Point3<ScalarType> & bary, vcg::Matrix33<ScalarType> &C){
// find the barycenter
FaceIterator fi;
ScalarType area = 0.0;
bary.Zero();
for(fi = m.face.begin(); fi != m.face.end(); ++fi)
if(!(*fi).IsD())
{
bary += vcg::Barycenter( *fi )* vcg::DoubleArea(*fi);
area+=vcg::DoubleArea(*fi);
}
bary/=area;
C.SetZero();
// C as covariance of triangle (0,0,0)(1,0,0)(0,1,0)
vcg::Matrix33<ScalarType> C0;
C0.SetZero();
C0[0][0] = C0[1][1] = 2.0;
C0[0][1] = C0[1][0] = 1.0;
C0*=1/24.0;
// integral of (x,y,0) in the same triangle
CoordType X(2/3.0,2/3.0,0);
vcg::Matrix33<ScalarType> A,At,DC; // matrix that bring the vertices to (v1-v0,v2-v0,n)
for(fi = m.face.begin(); fi != m.face.end(); ++fi)
if(!(*fi).IsD())
{
CoordType n = (*fi).N().Normalize();
const CoordType &P0 = (*fi).P(0);
const CoordType &P1 = (*fi).P(1);
const CoordType &P2 = (*fi).P(2);
A.SetColumn(0,P1-P0);
A.SetColumn(1,P2-P0);
A.SetColumn(2,n);
CoordType delta = P0 - bary;
At= A;
At.Transpose();
/* DC is calculated as integral of (A*x+delta) * (A*x+delta)^T over the triangle,
where delta = v0-bary
*/
DC.SetZero();
DC+= A*C0*At;
vcg::Matrix33<ScalarType> tmp;
tmp.OuterProduct(X,delta);
DC+= A * tmp;
tmp.Transpose();
DC+= tmp * At;
tmp.OuterProduct(delta,delta);
DC+=tmp*0.5;
DC*=fabs(A.Determinant()); // the determinant of A is the jacobian of the change of variables A*x+delta
C+=DC;
}
}
}; // end class Inertia
} // end namespace tri
} // end namespace vcg
#endif