PolyFEM
Loading...
Searching...
No Matches
TargetForms.hpp
Go to the documentation of this file.
1#pragma once
2
4
5#include <igl/AABB.h>
8
9namespace polyfem::solver
10{
12 {
13 public:
14 TargetForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const json &args) : SpatialIntegralForm(variable_to_simulations, state, args)
15 {
17
18 auto tmp_ids = args["surface_selection"].get<std::vector<int>>();
19 ids_ = std::set(tmp_ids.begin(), tmp_ids.end());
20 }
21 ~TargetForm() = default;
22
23 virtual std::string name() const override { return "target"; }
24
25 void set_reference(const std::shared_ptr<const State> &target_state, const std::set<int> &reference_cached_body_ids); // target is another simulation solution
26 void set_reference(const Eigen::VectorXd &disp) { target_disp = disp; } // target is a constant displacement
27 void set_reference(const json &func, const json &grad_func); // target is a lambda function depending on deformed position
28 void set_active_dimension(const std::vector<bool> &mask) { active_dimension_mask = mask; }
29
30 protected:
32
33 std::shared_ptr<const State> target_state_;
34 std::map<int, int> e_to_ref_e_;
35
36 std::vector<bool> active_dimension_mask;
37 Eigen::VectorXd target_disp;
38
39 bool have_target_func = false;
41 std::array<utils::ExpressionValue, 3> target_func_grad;
42 };
43
45 {
46 public:
47 SDFTargetForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const json &args) : SpatialIntegralForm(variable_to_simulations, state, args)
48 {
50
51 auto tmp_ids = args["surface_selection"].get<std::vector<int>>();
52 ids_ = std::set(tmp_ids.begin(), tmp_ids.end());
53 }
54
55 virtual std::string name() const override { return "sdf-target"; }
56
57 void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x) override;
58 void set_bspline_target(const Eigen::MatrixXd &control_points, const Eigen::VectorXd &knots, const double delta);
59 void set_bspline_target(const Eigen::MatrixXd &control_points, const Eigen::VectorXd &knots_u, const Eigen::VectorXd &knots_v, const double delta);
60
61 protected:
63
64 private:
65 void compute_distance(const Eigen::MatrixXd &point, double &distance) const;
66
67 int dim;
68 double delta_;
69
70 Eigen::MatrixXd t_or_uv_sampling;
71 Eigen::MatrixXd point_sampling;
73
74 std::unique_ptr<LazyCubicInterpolator> interpolation_fn;
75 };
76
78 {
79 public:
80 MeshTargetForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const json &args) : SpatialIntegralForm(variable_to_simulations, state, args)
81 {
83
84 auto tmp_ids = args["surface_selection"].get<std::vector<int>>();
85 ids_ = std::set(tmp_ids.begin(), tmp_ids.end());
86 }
87
88 virtual std::string name() const override { return "mesh-target"; }
89
90 void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x) override;
91 void set_surface_mesh_target(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const double delta);
92
93 protected:
95
96 private:
97 int dim;
98 double delta_;
99
100 Eigen::MatrixXd V_;
101 Eigen::MatrixXi F_;
102 igl::AABB<Eigen::MatrixXd, 3> tree_;
103
104 std::unique_ptr<LazyCubicInterpolator> interpolation_fn;
105 };
106
108 {
109 public:
110 NodeTargetForm(const State &state, const VariableToSimulationGroup &variable_to_simulations, const json &args);
111 NodeTargetForm(const State &state, const VariableToSimulationGroup &variable_to_simulations, const std::vector<int> &active_nodes_, const Eigen::MatrixXd &target_vertex_positions_);
112 ~NodeTargetForm() = default;
113
114 std::string name() const override { return "node-target"; }
115
116 Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override;
117 double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override;
118 void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
119
120 protected:
121 const State &state_;
122
123 Eigen::MatrixXd target_vertex_positions;
124 std::vector<int> active_nodes;
125 };
126
128 {
129 public:
130 BarycenterTargetForm(const VariableToSimulationGroup &variable_to_simulations, const json &args, const std::shared_ptr<State> &state1, const std::shared_ptr<State> &state2);
131
132 Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override;
133 void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
134 double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override;
135
136 private:
137 std::vector<std::unique_ptr<PositionForm>> center1, center2;
138 int dim;
139 };
140
142 {
143 public:
144 MinTargetDistForm(const VariableToSimulationGroup &variable_to_simulations, const std::vector<int> &steps, const Eigen::VectorXd &target, const json &args, const std::shared_ptr<State> &state);
145 virtual ~MinTargetDistForm() = default;
146
147 Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const override;
148 void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
149
150 protected:
151 double value_unweighted(const Eigen::VectorXd &x) const override;
152
153 private:
154 static double eval1(Eigen::VectorXd &x)
155 {
156 return 1. / x.array().inverse().sum();
157 }
158 static Eigen::VectorXd eval1_grad(Eigen::VectorXd &x)
159 {
160 return x.array().pow(-2.0) * pow(eval1(x), 2);
161 }
162
163 double eval2(Eigen::VectorXd &x) const
164 {
165 assert(x.size() == dim+1);
166 double out = 0;
167 for (int d = 0; d < dim; d++)
168 out += pow(x(d) / x(dim) - target_(d), 2);
169 return out;
170 }
171 Eigen::VectorXd eval2_grad(Eigen::VectorXd &x) const
172 {
173 Eigen::VectorXd g = Eigen::VectorXd::Zero(dim + 1);
174 for (int d = 0; d < dim; d++)
175 {
176 const double tmp = 2 * (x(d) / x(dim) - target_(d));
177 g(d) += tmp / x(dim);
178 g(dim) += tmp * (-x(d) / x(dim) / x(dim));
179 }
180 return g;
181 }
182
183 const std::vector<int> steps_;
184 const Eigen::VectorXd target_;
185
186 int dim;
187 std::vector<std::unique_ptr<StaticForm>> objs;
188 };
189} // namespace polyfem::solver
int V
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
std::vector< std::unique_ptr< PositionForm > > center1
std::vector< std::unique_ptr< PositionForm > > center2
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override
void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x) override
IntegrableFunctional get_integral_functional() const override
igl::AABB< Eigen::MatrixXd, 3 > tree_
virtual std::string name() const override
void set_surface_mesh_target(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const double delta)
MeshTargetForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const json &args)
std::unique_ptr< LazyCubicInterpolator > interpolation_fn
static double eval1(Eigen::VectorXd &x)
double eval2(Eigen::VectorXd &x) const
std::vector< std::unique_ptr< StaticForm > > objs
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
const std::vector< int > steps_
static Eigen::VectorXd eval1_grad(Eigen::VectorXd &x)
virtual ~MinTargetDistForm()=default
Eigen::VectorXd eval2_grad(Eigen::VectorXd &x) const
Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const override
std::string name() const override
double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Eigen::MatrixXd target_vertex_positions
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override
void compute_distance(const Eigen::MatrixXd &point, double &distance) const
virtual std::string name() const override
void set_bspline_target(const Eigen::MatrixXd &control_points, const Eigen::VectorXd &knots, const double delta)
std::unique_ptr< LazyCubicInterpolator > interpolation_fn
void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x) override
IntegrableFunctional get_integral_functional() const override
SDFTargetForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const json &args)
void set_integral_type(const SpatialIntegralType type)
std::shared_ptr< const State > target_state_
void set_reference(const std::shared_ptr< const State > &target_state, const std::set< int > &reference_cached_body_ids)
TargetForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const json &args)
IntegrableFunctional get_integral_functional() const override
std::array< utils::ExpressionValue, 3 > target_func_grad
utils::ExpressionValue target_func
std::map< int, int > e_to_ref_e_
std::vector< bool > active_dimension_mask
void set_reference(const Eigen::VectorXd &disp)
virtual std::string name() const override
void set_active_dimension(const std::vector< bool > &mask)
A collection of VariableToSimulation.
nlohmann::json json
Definition Common.hpp:9