MeshLib
 
Loading...
Searching...
No Matches
MRLaplacian.h
Go to the documentation of this file.
1#pragma once
2
3#include "MRBitSet.h"
4#include "MRVector.h"
5#include "MRVector3.h"
6#include "MREnums.h"
7
8#pragma warning(push)
9#pragma warning(disable: 4068) // unknown pragmas
10#pragma warning(disable: 4127) // conditional expression is constant
11#pragma warning(disable: 4464) // relative include path contains '..'
12#pragma warning(disable: 5054) // operator '|': deprecated between enumerations of different types
13#pragma clang diagnostic push
14#pragma clang diagnostic ignored "-Wdeprecated-anon-enum-enum-conversion"
15#pragma clang diagnostic ignored "-Wunknown-warning-option" // for next one
16#pragma clang diagnostic ignored "-Wunused-but-set-variable" // for newer clang
17#include <Eigen/SparseCore>
18#pragma clang diagnostic pop
19#pragma warning(pop)
20
21namespace MR
22{
23
24// Laplacian to smoothly deform a region preserving mesh fine details.
25// How to use:
26// 1. Initialize Laplacian for the region being deformed, here region properties are remembered.
27// 2. Change positions of some vertices within the region and call fixVertex for them.
28// 3. Optionally call updateSolver()
29// 4. Call apply() to change the remaining vertices within the region
30// Then steps 1-4 or 2-4 can be repeated.
32{
33public:
34 enum class RememberShape
35 {
36 Yes, // true Laplacian mode when initial mesh shape is remembered and copied in apply
37 No // ignore initial mesh shape in the region and just position vertices smoothly in the region
38 };
39
40 Laplacian( Mesh & mesh ) : mesh_( mesh ) { }
41 // initialize Laplacian for the region being deformed, here region properties are remembered and precomputed;
42 // \param freeVerts must not include all vertices of a mesh connected component
44 // notify Laplacian that given vertex has changed after init and must be fixed during apply;
45 // \param smooth whether to make the surface smooth in this vertex (sharp otherwise)
46 MRMESH_API void fixVertex( VertId v, bool smooth = true );
47 // sets position of given vertex after init and it must be fixed during apply (THIS METHOD CHANGES THE MESH);
48 // \param smooth whether to make the surface smooth in this vertex (sharp otherwise)
49 MRMESH_API void fixVertex( VertId v, const Vector3f & fixedPos, bool smooth = true );
50 // if you manually call this method after initialization and fixing vertices then next apply call will be much faster
52 // given fixed vertices, computes positions of remaining region vertices
54 // given a pre-resized scalar field with set values in fixed vertices, computes the values in free vertices
55 MRMESH_API void applyToScalar( VertScalars & scalarField );
56
57 // return all initially free vertices and the first layer around the them
58 const VertBitSet & region() const { return region_; }
59 // return currently free vertices
60 const VertBitSet & freeVerts() const { return freeVerts_; }
61 // return fixed vertices from the first layer around free vertices
62 VertBitSet firstLayerFixedVerts() const { assert( solverValid_ ); return firstLayerFixedVerts_; }
63
64 using EdgeWeights [[deprecated]] = MR::EdgeWeights;
65
66private:
67 // updates solver_ only
68 void updateSolver_();
69 // updates rhs_ only
70 void updateRhs_();
71 template <typename I, typename G, typename S>
72 void prepareRhs_( I && iniRhs, G && g, S && s );
73
74 Mesh & mesh_;
75
76 // all initially free vertices and the first layer around the them
77 VertBitSet region_;
78
79 // currently free vertices
80 VertBitSet freeVerts_;
81
82 // fixed vertices where no smoothness is required
83 VertBitSet fixedSharpVertices_;
84
85 // fixed vertices from the first layer around free vertices
86 VertBitSet firstLayerFixedVerts_;
87
88 // for all vertices in the region
89 struct Equation
90 {
91 Vector3d rhs; // equation right hand side
92 double centerCoeff = 0; // coefficient on matrix diagonal
93 int firstElem = 0; // index in nonZeroElements_
94 };
95 std::vector<Equation> equations_;
96
97 struct Element
98 {
99 double coeff = 0;
100 VertId neiVert;
101 };
102 std::vector<Element> nonZeroElements_;
103
104 // map from vertex index to matrix row/col
105 Vector< int, VertId > regionVert2id_;
106 Vector< int, VertId > freeVert2id_;
107
108 using SparseMatrix = Eigen::SparseMatrix<double,Eigen::RowMajor>;
109 SparseMatrix M_;
110
111 // if true then we do not need to recompute solver_ in the apply
112 bool solverValid_ = false;
113 using SparseMatrixColMajor = Eigen::SparseMatrix<double,Eigen::ColMajor>;
114
115 // interface needed to hide implementation headers
116 class Solver
117 {
118 public:
119 virtual ~Solver() = default;
120 virtual void compute( const SparseMatrixColMajor& A ) = 0;
121 virtual Eigen::VectorXd solve( const Eigen::VectorXd& rhs ) = 0;
122 };
123 std::unique_ptr<Solver> solver_;
124
125 // if true then we do not need to recompute rhs_ in the apply
126 bool rhsValid_ = false;
127 Eigen::VectorXd rhs_[3];
128};
129
130} //namespace MR
int VertId
Definition MRDotNet/MRMeshFwd.h:51
#define MRMESH_API
Definition MRMesh/MRMeshFwd.h:46
Definition MRDotNet/MRBitSet.h:39
Definition MRLaplacian.h:32
VertBitSet firstLayerFixedVerts() const
Definition MRLaplacian.h:62
const VertBitSet & freeVerts() const
Definition MRLaplacian.h:60
MRMESH_API void apply()
MRMESH_API void fixVertex(VertId v, const Vector3f &fixedPos, bool smooth=true)
MRMESH_API void updateSolver()
MRMESH_API void init(const VertBitSet &freeVerts, EdgeWeights weights, RememberShape rem=Laplacian::RememberShape::Yes)
RememberShape
Definition MRLaplacian.h:35
Laplacian(Mesh &mesh)
Definition MRLaplacian.h:40
const VertBitSet & region() const
Definition MRLaplacian.h:58
MRMESH_API void fixVertex(VertId v, bool smooth=true)
MRMESH_API void applyToScalar(VertScalars &scalarField)
represents a 3-dimentional float-typed vector
Definition MRDotNet/MRVector3.h:8
Definition MRCameraOrientationPlugin.h:7
EdgeWeights
determines the weight of each edge in applications like Laplacian
Definition MREnums.h:8
I
Definition MRMesh/MRMeshFwd.h:88
Definition MRMesh/MRMesh.h:23