Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #ifndef __VCGLIB_HARMONIC_FIELD
00024 #define __VCGLIB_HARMONIC_FIELD
00025
00026 #include <vcg/complex/complex.h>
00027 #include <eigenlib/Eigen/Sparse>
00028
00029 namespace vcg {
00030 namespace tri {
00031
00032 template <class MeshType, typename Scalar = double>
00033 class Harmonic
00034 {
00035 public:
00036 typedef typename MeshType::VertexType VertexType;
00037 typedef typename MeshType::FaceType FaceType;
00038 typedef typename MeshType::CoordType CoordType;
00039 typedef typename MeshType::ScalarType ScalarType;
00040
00041 typedef double CoeffScalar;
00042
00043 typedef typename std::pair<VertexType *, Scalar> Constraint;
00044 typedef typename std::vector<Constraint> ConstraintVec;
00045 typedef typename ConstraintVec::const_iterator ConstraintIt;
00046
00057 template <typename ACCESSOR>
00058 static bool ComputeScalarField(MeshType & m, const ConstraintVec & constraints, ACCESSOR field, bool biharmonic = false)
00059 {
00060 typedef Eigen::SparseMatrix<CoeffScalar> SpMat;
00061 typedef Eigen::Triplet<CoeffScalar> Triple;
00062
00063 RequirePerVertexFlags(m);
00064 RequireCompactness(m);
00065 RequireFFAdjacency(m);
00066 MeshAssert<MeshType>::FFAdjacencyIsInitialized(m);
00067 MeshAssert<MeshType>::NoUnreferencedVertex(m);
00068
00069 if (constraints.empty())
00070 return false;
00071
00072 int n = m.VN();
00073
00074
00075 std::vector<Triple> coeffs;
00076 std::map<size_t,CoeffScalar> sums;
00077
00078 vcg::tri::UpdateFlags<MeshType>::FaceClearV(m);
00079 for (size_t i = 0; i < m.face.size(); ++i)
00080 {
00081 FaceType & f = m.face[i];
00082 assert(!f.IsV());
00083
00084 f.SetV();
00085
00086
00087 for (int edge = 0; edge < 3; ++edge)
00088 {
00089 CoeffScalar weight;
00090 WeightInfo res = CotangentWeightIfNotVisited(f, edge, weight);
00091
00092 if (res == EdgeAlreadyVisited) continue;
00093 assert(res == Success);
00094
00095
00096 size_t v0_idx = vcg::tri::Index(m, f.V0(edge));
00097 size_t v1_idx = vcg::tri::Index(m, f.V1(edge));
00098
00099 coeffs.push_back(Triple(v0_idx, v1_idx, -weight));
00100 coeffs.push_back(Triple(v1_idx, v0_idx, -weight));
00101
00102
00103 sums[v0_idx] += weight;
00104 sums[v1_idx] += weight;
00105 }
00106 }
00107
00108
00109 SpMat laplaceMat;
00110 laplaceMat.resize(n, n);
00111 laplaceMat.reserve(coeffs.size());
00112 for (std::map<size_t,CoeffScalar>::const_iterator it = sums.begin(); it != sums.end(); ++it)
00113 {
00114 coeffs.push_back(Triple(it->first, it->first, it->second));
00115 }
00116 laplaceMat.setFromTriplets(coeffs.begin(), coeffs.end());
00117
00118 if (biharmonic)
00119 {
00120 SpMat lap_t = laplaceMat;
00121 lap_t.transpose();
00122 laplaceMat = lap_t * laplaceMat;
00123 }
00124
00125
00126
00127 const CoeffScalar alpha = pow(10.0, 8.0);
00128
00129
00130 Eigen::Matrix<CoeffScalar, Eigen::Dynamic, 1> b, x;
00131 b.setZero(n);
00132
00133 for (ConstraintIt it=constraints.begin(); it!=constraints.end(); it++)
00134 {
00135 size_t v_idx = vcg::tri::Index(m, it->first);
00136 b(v_idx) = alpha * it->second;
00137 laplaceMat.coeffRef(v_idx, v_idx) += alpha;
00138 }
00139
00140
00141 Eigen::SimplicialLDLT<SpMat> solver;
00142 solver.compute(laplaceMat);
00143
00144 if(solver.info() != Eigen::Success)
00145 {
00146
00147 switch (solver.info())
00148 {
00149
00150 case Eigen::NumericalIssue :
00151 case Eigen::NoConvergence :
00152 case Eigen::InvalidInput :
00153 default : return false;
00154 }
00155 }
00156
00157
00158 x = solver.solve(b);
00159 if(solver.info() != Eigen::Success)
00160 {
00161 return false;
00162 }
00163
00164
00165 for (size_t i = 0; int(i) < n; ++i)
00166 {
00167 field[i] = Scalar(x[i]);
00168 }
00169
00170 return true;
00171 }
00172
00173 enum WeightInfo
00174 {
00175 Success = 0,
00176 EdgeAlreadyVisited
00177 };
00178
00179
00191 template <typename ScalarT>
00192 static WeightInfo CotangentWeightIfNotVisited(const FaceType &f, int edge, ScalarT & weight)
00193 {
00194 const FaceType * fp = f.cFFp(edge);
00195 if (fp != NULL && fp != &f)
00196 {
00197
00198 if (fp->IsV()) return EdgeAlreadyVisited;
00199 }
00200
00201 weight = CotangentWeight<ScalarT>(f, edge);
00202
00203 return Success;
00204 }
00205
00213 template <typename ScalarT>
00214 static ScalarT CotangentWeight(const FaceType &f, int edge)
00215 {
00216 assert(edge >= 0 && edge < 3);
00217
00218
00219 const FaceType * fp = f.cFFp(edge);
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232 ScalarT cotA = 0;
00233 ScalarT cotB = 0;
00234
00235
00236 VertexType * v0 = f.cV(edge);
00237 VertexType * v1 = f.cV((edge+1)%f.VN());
00238
00239 if (fp != NULL &&
00240 fp != &f)
00241 {
00242
00243 VertexType * vb = fp->cV((f.cFFi(edge)+2)%fp->VN());
00244 ScalarT angleB = ComputeAngle<ScalarT>(v0, vb, v1);
00245 cotB = vcg::math::Cos(angleB) / vcg::math::Sin(angleB);
00246 }
00247
00248 VertexType * va = f.cV((edge+2)%f.VN());
00249 ScalarT angleA = ComputeAngle<ScalarT>(v0, va, v1);
00250 cotA = vcg::math::Cos(angleA) / vcg::math::Sin(angleA);
00251
00252 return (cotA + cotB) / 2;
00253 }
00254
00255 template <typename ScalarT>
00256 static ScalarT ComputeAngle(const VertexType * a, const VertexType * b, const VertexType * c)
00257 {
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268 assert(a != NULL && b != NULL && c != NULL);
00269 Point3<ScalarT> A,B,C;
00270 A.Import(a->P());
00271 B.Import(b->P());
00272 C.Import(c->P());
00273 ScalarT angle = vcg::Angle(A - B, C - B);
00274 return angle;
00275 }
00276 };
00277
00278 }
00279 }
00280 #endif // __VCGLIB_HARMONIC_FIELD