PolyFEM
Loading...
Searching...
No Matches
BarrierForms.hpp
Go to the documentation of this file.
1#pragma once
2
5
8
9namespace polyfem::solver
10{
12 {
13 public:
14 CollisionBarrierForm(const VariableToSimulationGroup &variable_to_simulation, const State &state, const double dhat, const double dmin = 0);
15
16 double value_unweighted(const Eigen::VectorXd &x) const override;
17
18 void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
19
20 void solution_changed(const Eigen::VectorXd &x) override;
21
22 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
23
24 std::string name() const override { return "collision barrier"; }
25
26 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
27
28 protected:
29 void build_collision_set(const Eigen::MatrixXd &displaced_surface);
30
31 Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const;
32
33 const State &state_;
34
35 Eigen::VectorXd X_init;
36
37 ipc::CollisionMesh collision_mesh_;
38 ipc::Collisions collision_set;
39 const double dhat_;
40 const double dmin_;
41 ipc::BroadPhaseMethod broad_phase_method_;
42
43 ipc::BarrierPotential barrier_potential_;
44 };
45
47 {
48 public:
49 LayerThicknessForm(const VariableToSimulationGroup &variable_to_simulations,
50 const State &state,
51 const std::vector<int> &boundary_ids,
52 const double dhat,
53 const bool use_log_barrier = false,
54 const double dmin = 0);
55
56 std::string name() const override { return "layer thickness"; }
57
58 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override { return 1.; }
59
60 protected:
62
63 std::vector<int> boundary_ids_;
64 std::map<int, std::set<int>> boundary_ids_to_dof_;
65 Eigen::MatrixXi can_collide_cache_;
66 };
67
69 {
70 public:
71 DeformedCollisionBarrierForm(const VariableToSimulationGroup &variable_to_simulation, const State &state, const double dhat);
72
73 std::string name() const override { return "deformed_collision_barrier"; }
74
75 double value_unweighted(const Eigen::VectorXd &x) const override;
76
77 void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
78
79 void solution_changed(const Eigen::VectorXd &x) override;
80
81 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
82
83 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
84
85 private:
86 void build_collision_set(const Eigen::MatrixXd &displaced_surface);
87
88 Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const;
89
90 const State &state_;
91
92 Eigen::VectorXd X_init;
93
94 ipc::CollisionMesh collision_mesh_;
95 ipc::Collisions collision_set;
96 const double dhat_;
97 ipc::BroadPhaseMethod broad_phase_method_;
98
99 const ipc::BarrierPotential barrier_potential_;
100 };
101} // namespace polyfem::solver
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
ipc::BarrierPotential barrier_potential_
std::string name() const override
void build_collision_set(const Eigen::MatrixXd &displaced_surface)
ipc::BroadPhaseMethod broad_phase_method_
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine the maximum step size allowable between the current and next solution.
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
void solution_changed(const Eigen::VectorXd &x) override
Update cached fields upon a change in the solution.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Checks if the step is collision free.
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine the maximum step size allowable between the current and next solution.
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
const ipc::BarrierPotential barrier_potential_
void build_collision_set(const Eigen::MatrixXd &displaced_surface)
Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const
void solution_changed(const Eigen::VectorXd &x) override
Update cached fields upon a change in the solution.
double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine the maximum step size allowable between the current and next solution.
std::map< int, std::set< int > > boundary_ids_to_dof_
std::string name() const override
A collection of VariableToSimulation.