PolyFEM
Loading...
Searching...
No Matches
Parametrizations.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Parametrization.hpp"
4#include <polyfem/Common.hpp>
5#include <map>
6
7namespace polyfem::mesh
8{
9 class Mesh;
10}
11
12namespace polyfem::basis
13{
14 class ElementBases;
15}
16
17namespace polyfem::solver
18{
20 {
21 private:
23
24 public:
25 static std::vector<std::shared_ptr<Parametrization>> build(const json &params, const int full_size);
26 };
27
29 {
30 public:
31 ExponentialMap(const int from = -1, const int to = -1);
32
33 int size(const int x_size) const override { return x_size; }
34 Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override;
35 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
36 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
37
38 private:
39 const int from_, to_;
40 };
41
42 class Scaling : public Parametrization
43 {
44 public:
45 Scaling(const double scale, const int from = -1, const int to = -1);
46
47 int size(const int x_size) const override { return x_size; }
48 Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override;
49 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
50 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
51
52 private:
53 const int from_, to_;
54 const double scale_;
55 };
56
58 {
59 public:
60 PowerMap(const double power = 1, const int from = -1, const int to = -1) : power_(power), from_(from), to_(to)
61 {
62 assert(from_ < to_ || from_ < 0);
63 assert(power_ > 0);
64 }
65
66 int size(const int x_size) const override { return x_size; }
67 Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override;
68 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
69 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
70
71 private:
72 const double power_;
73 const int from_, to_;
74 };
75
77 {
78 public:
79 ENu2LambdaMu(const bool is_volume);
80
81 int size(const int x_size) const override { return x_size; }
82
83 Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override;
84 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
85 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
86
87 private:
88 const bool is_volume_;
89 };
90
92 {
93 public:
94 PerBody2PerNode(const mesh::Mesh &mesh, const std::vector<basis::ElementBases> &bases, const int n_bases);
95
96 int size(const int x_size) const override;
97 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
98 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
99
100 private:
102 const std::vector<basis::ElementBases> &bases_;
105 Eigen::VectorXi node_id_to_body_id_;
106 };
107
109 {
110 public:
111 PerBody2PerElem(const mesh::Mesh &mesh);
112
113 int size(const int x_size) const override;
114 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
115 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
116
117 private:
121 std::map<int, std::array<int, 2>> body_id_map_; // from body_id to {elem_id, index}
122 };
123
125 {
126 public:
127 SliceMap(const int from = -1, const int to = -1, const int total = -1);
128
129 int size(const int x_size) const override { return to_ - from_; }
130
131 Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override;
132 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
133 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
134
135 private:
136 const int from_, to_, total_;
137 };
138
140 {
141 public:
142 InsertConstantMap(const int size = -1, const double val = 0, const int start_index = -1);
143 InsertConstantMap(const Eigen::VectorXd &values, const int start_index = -1);
144
145 int size(const int x_size) const override;
146 Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override;
147 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
148 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
149
150 private:
151 // const int size_;
152 // const double val_;
153 int start_index_ = -1;
154 Eigen::VectorXd values_;
155 };
156
158 {
159 public:
160 LinearFilter(const mesh::Mesh &mesh, const double radius);
161
162 int size(const int x_size) const override { return x_size; }
163 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
164 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
165
166 private:
167 Eigen::SparseMatrix<double> tt_radius_adjacency;
169 };
170
172 {
173 public:
174 ScalarVelocityParametrization(const double start_val, const double dt) : start_val_(start_val), dt_(dt) {}
175
176 int size(const int x_size) const override { return x_size; }
177 Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override;
178 Eigen::VectorXd eval(const Eigen::VectorXd &x) const override;
179 Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override;
180
181 private:
182 const double start_val_;
183 const double dt_;
184 };
185} // namespace polyfem::solver
double val
Definition Assembler.cpp:86
int y
int x
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
int size(const int x_size) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
int size(const int x_size) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
int size(const int x_size) const override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
Eigen::SparseMatrix< double > tt_radius_adjacency
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
int size(const int x_size) const override
static std::vector< std::shared_ptr< Parametrization > > build(const json &params, const int full_size)
This parameterize a function f : x -> y and provides the chain rule with respect to previous gradient...
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
int size(const int x_size) const override
std::map< int, std::array< int, 2 > > body_id_map_
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
int size(const int x_size) const override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
const std::vector< basis::ElementBases > & bases_
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
int size(const int x_size) const override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
PowerMap(const double power=1, const int from=-1, const int to=-1)
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
ScalarVelocityParametrization(const double start_val, const double dt)
int size(const int x_size) const override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
int size(const int x_size) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
Eigen::VectorXd apply_jacobian(const Eigen::VectorXd &grad, const Eigen::VectorXd &x) const override
Eigen::VectorXd inverse_eval(const Eigen::VectorXd &y) override
int size(const int x_size) const override
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
nlohmann::json json
Definition Common.hpp:9