PolyFEM
Loading...
Searching...
No Matches
BSplineParametrization.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <map>
4#include <vector>
5#include <iostream>
6#include <Eigen/Dense>
7
8#include <nanospline/BSpline.h>
9#include <nanospline/BSplinePatch.h>
10
11namespace polyfem
12{
14 {
15 public:
16 BSplineParametrization(const Eigen::MatrixXd &V) { num_vertices = V.rows(); }
17 virtual ~BSplineParametrization() = default;
18
19 virtual int vertex_size() = 0;
20 virtual void reparametrize(const Eigen::MatrixXd &control_points, Eigen::MatrixXd &newV) = 0;
21 // virtual void get_parameters(const Eigen::MatrixXd &V, Eigen::MatrixXd &control_points) final;
22 virtual void get_parameters(const Eigen::MatrixXd &V, Eigen::MatrixXd &control_points, const bool mesh_changed) = 0;
23 virtual void derivative_wrt_params(const Eigen::VectorXd &grad_boundary, Eigen::VectorXd &grad_control_points) = 0;
24
25 protected:
27 };
28
30 {
31 public:
32 BSplineParametrization2D(const Eigen::MatrixXd &control_points, const Eigen::MatrixXd &knots, const Eigen::MatrixXd &V);
33
34 int vertex_size() override { return node_id_to_t_.size(); }
35
36 void reparametrize(const Eigen::MatrixXd &control_points, Eigen::MatrixXd &newV) override;
37
38 // Assume the connectivity has not changed. Does not work with remeshing
39 void get_parameters(const Eigen::MatrixXd &V, Eigen::MatrixXd &control_points, const bool mesh_changed) override;
40
41 void derivative_wrt_params(const Eigen::VectorXd &grad_boundary, Eigen::VectorXd &grad_control_points) override;
42
43 static void gradient(const Eigen::MatrixXd &point, const Eigen::MatrixXd &control_points, const double t_parameter, const double distance, Eigen::MatrixXd &grad);
44 static void eval(const Eigen::MatrixXd &control_points, const double t, Eigen::MatrixXd &val);
45 static void deriv(const Eigen::MatrixXd &control_points, const double t, Eigen::MatrixXd &val);
46
47 private:
48 std::vector<int> node_ids_;
49 std::map<int, double> node_id_to_t_;
50 const int dim;
51 nanospline::BSpline<double, 2, 3> curve;
52 };
53
55 {
56 public:
57 BSplineParametrization3D(const Eigen::MatrixXd &control_points, const Eigen::MatrixXd &knots_u, const Eigen::MatrixXd &knots_v, const Eigen::MatrixXd &V);
58
59 int vertex_size() override { return node_id_to_param_.size(); }
60
61 // Assume the connectivity has not changed. Does not work with remeshing
62 void get_parameters(const Eigen::MatrixXd &V, Eigen::MatrixXd &control_points, const bool mesh_changed) override;
63
64 void reparametrize(const Eigen::MatrixXd &control_points, Eigen::MatrixXd &newV) override;
65
66 void derivative_wrt_params(const Eigen::VectorXd &grad_boundary, Eigen::VectorXd &grad_control_points) override;
67
68 static void gradient(const Eigen::MatrixXd &point, const Eigen::MatrixXd &control_points, const Eigen::MatrixXd &uv_parameter, const double distance, Eigen::MatrixXd &grad);
69 static void eval(const Eigen::MatrixXd &control_points, const Eigen::MatrixXd &uv_parameter, Eigen::MatrixXd &val);
70 static void deriv(const Eigen::MatrixXd &control_points, const Eigen::MatrixXd &uv_parameter, Eigen::MatrixXd &deriv_u, Eigen::MatrixXd &deriv_v);
71
72 private:
73 std::vector<int> node_ids_;
74 std::map<int, Eigen::MatrixXd> node_id_to_param_;
75 const int dim;
76 nanospline::BSplinePatch<double, 3, 3, 3> patch;
77 };
78
79} // namespace polyfem
int V
double val
Definition Assembler.cpp:86
void derivative_wrt_params(const Eigen::VectorXd &grad_boundary, Eigen::VectorXd &grad_control_points) override
void reparametrize(const Eigen::MatrixXd &control_points, Eigen::MatrixXd &newV) override
static void eval(const Eigen::MatrixXd &control_points, const double t, Eigen::MatrixXd &val)
static void deriv(const Eigen::MatrixXd &control_points, const double t, Eigen::MatrixXd &val)
nanospline::BSpline< double, 2, 3 > curve
static void gradient(const Eigen::MatrixXd &point, const Eigen::MatrixXd &control_points, const double t_parameter, const double distance, Eigen::MatrixXd &grad)
void get_parameters(const Eigen::MatrixXd &V, Eigen::MatrixXd &control_points, const bool mesh_changed) override
static void gradient(const Eigen::MatrixXd &point, const Eigen::MatrixXd &control_points, const Eigen::MatrixXd &uv_parameter, const double distance, Eigen::MatrixXd &grad)
static void deriv(const Eigen::MatrixXd &control_points, const Eigen::MatrixXd &uv_parameter, Eigen::MatrixXd &deriv_u, Eigen::MatrixXd &deriv_v)
std::map< int, Eigen::MatrixXd > node_id_to_param_
void get_parameters(const Eigen::MatrixXd &V, Eigen::MatrixXd &control_points, const bool mesh_changed) override
static void eval(const Eigen::MatrixXd &control_points, const Eigen::MatrixXd &uv_parameter, Eigen::MatrixXd &val)
nanospline::BSplinePatch< double, 3, 3, 3 > patch
void derivative_wrt_params(const Eigen::VectorXd &grad_boundary, Eigen::VectorXd &grad_control_points) override
void reparametrize(const Eigen::MatrixXd &control_points, Eigen::MatrixXd &newV) override
virtual ~BSplineParametrization()=default
virtual void derivative_wrt_params(const Eigen::VectorXd &grad_boundary, Eigen::VectorXd &grad_control_points)=0
virtual void reparametrize(const Eigen::MatrixXd &control_points, Eigen::MatrixXd &newV)=0
BSplineParametrization(const Eigen::MatrixXd &V)
virtual void get_parameters(const Eigen::MatrixXd &V, Eigen::MatrixXd &control_points, const bool mesh_changed)=0