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
42 // initialize Laplacian for the region being deformed, here region properties are remembered and precomputed;
43 // \param freeVerts must not include all vertices of a mesh connected component
45
46 // notify Laplacian that given vertex has changed after init and must be fixed during apply;
47 // \param smooth whether to make the surface smooth in this vertex (sharp otherwise)
48 MRMESH_API void fixVertex( VertId v, bool smooth = true );
49
50 // sets position of given vertex after init and it must be fixed during apply (THIS METHOD CHANGES THE MESH);
51 // \param smooth whether to make the surface smooth in this vertex (sharp otherwise)
52 MRMESH_API void fixVertex( VertId v, const Vector3f & fixedPos, bool smooth = true );
53
54 // if you manually call this method after initialization and fixing vertices then next apply call will be much faster
56
57 // given fixed vertices, computes positions of remaining region vertices
59
60 // given a pre-resized scalar field with set values in fixed vertices, computes the values in free vertices
61 MRMESH_API void applyToScalar( VertScalars & scalarField );
62
63 // return all initially free vertices and the first layer of vertices around them
64 const VertBitSet & region() const { return region_; }
65
66 // return currently free vertices
67 const VertBitSet & freeVerts() const { return freeVerts_; }
68
69 // return fixed vertices from the first layer around free vertices
70 VertBitSet firstLayerFixedVerts() const { assert( solverValid_ ); return firstLayerFixedVerts_; }
71
72 using EdgeWeights [[deprecated]] = MR::EdgeWeights;
73
74private:
75 // updates solver_ only
76 void updateSolver_();
77
78 // updates rhs_ only
79 void updateRhs_();
80
81 template <typename I, typename G, typename S>
82 void prepareRhs_( I && iniRhs, G && g, S && s );
83
84 Mesh & mesh_;
85
86 // all initially free vertices and the first layer of vertices around them
87 VertBitSet region_;
88
89 // currently free vertices
90 VertBitSet freeVerts_;
91
92 // fixed vertices where no smoothness is required
93 VertBitSet fixedSharpVertices_;
94
95 // fixed vertices from the first layer around free vertices
96 VertBitSet firstLayerFixedVerts_;
97
98 // for all vertices in the region
99 struct Equation
100 {
101 Vector3d rhs; // equation right hand side
102 double centerCoeff = 0; // coefficient on matrix diagonal
103 int firstElem = 0; // index in nonZeroElements_
104 };
105 std::vector<Equation> equations_;
106
107 struct Element
108 {
109 double coeff = 0;
110 VertId neiVert;
111 };
112 std::vector<Element> nonZeroElements_;
113
114 // map from vertex index to matrix row/col
115 Vector< int, VertId > regionVert2id_;
116 Vector< int, VertId > freeVert2id_;
117
118 using SparseMatrix = Eigen::SparseMatrix<double,Eigen::RowMajor>;
119 SparseMatrix M_;
120
121 // if true then we do not need to recompute solver_ in the apply
122 bool solverValid_ = false;
123 using SparseMatrixColMajor = Eigen::SparseMatrix<double,Eigen::ColMajor>;
124
125 // interface needed to hide implementation headers
126 class Solver
127 {
128 public:
129 virtual ~Solver() = default;
130 virtual void compute( const SparseMatrixColMajor& A ) = 0;
131 virtual Eigen::VectorXd solve( const Eigen::VectorXd& rhs ) = 0;
132 };
133 std::unique_ptr<Solver> solver_;
134
135 // if true then we do not need to recompute rhs_ in the apply
136 bool rhsValid_ = false;
137 Eigen::VectorXd rhs_[3];
138};
139
140} //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:70
const VertBitSet & freeVerts() const
Definition MRLaplacian.h:67
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:64
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:8
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