MeshLib C++ Docs
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#include <MRPch/MREigenSparseCore.h>
8
9namespace MR
10{
11
20class Laplacian
21{
22public:
23 enum class RememberShape
24 {
25 Yes, // true Laplacian mode when initial mesh shape is remembered and copied in apply
26 No // ignore initial mesh shape in the region and just position vertices smoothly in the region
27 };
28
29 MRMESH_API explicit Laplacian( Mesh & mesh );
30 Laplacian( const MeshTopology & topology, VertCoords & points ) : topology_( topology ), points_( points ) { }
31
34 MRMESH_API void init( const VertBitSet & freeVerts, EdgeWeights weights, VertexMass vmass = VertexMass::Unit,
35 RememberShape rem = Laplacian::RememberShape::Yes );
36
39 MRMESH_API void fixVertex( VertId v, bool smooth = true );
40
43 MRMESH_API void fixVertex( VertId v, const Vector3f & fixedPos, bool smooth = true );
44
47
49 MRMESH_API void apply();
50
52 MRMESH_API void applyToScalar( VertScalars & scalarField );
53
55 [[nodiscard]] const VertBitSet & region() const { return region_; }
56
58 [[nodiscard]] const VertBitSet & freeVerts() const { return freeVerts_; }
59
61 [[nodiscard]] const VertBitSet & firstLayerFixedVerts() const { assert( solverValid_ ); return firstLayerFixedVerts_; }
62
63private:
64 // updates solver_ only
65 void updateSolver_();
66
67 // updates rhs_ only
68 void updateRhs_();
69
70 template <typename I, typename G, typename S>
71 void prepareRhs_( I && iniRhs, G && g, S && s );
72
73 const MeshTopology & topology_;
74 VertCoords & points_;
75
76 // all initially free vertices and the first layer of vertices around 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
#define MRMESH_API
Definition MRMeshFwd.h:80
unsafe Laplacian(MR._ByValue_Laplacian _other)
unsafe void fixVertex(MR.VertId v, bool? smooth=null)
unsafe void updateSolver()
unsafe void init(MR.Const_VertBitSet freeVerts, MR.EdgeWeights weights, MR.VertexMass? vmass=null, MR.Laplacian.RememberShape? rem=null)
unsafe void applyToScalar(MR.VertScalars scalarField)
unsafe void apply()
Element
Definition MRImGuiMeasurementIndicators.h:99
Definition MRCameraOrientationPlugin.h:8