PolyFEM
Loading...
Searching...
No Matches
BarrierForms.hpp
Go to the documentation of this file.
1#pragma once
2
5
7#include <ipc/potentials/barrier_potential.hpp>
8#include <ipc/smooth_contact/smooth_collisions.hpp>
9#include <ipc/smooth_contact/smooth_contact_potential.hpp>
11
12namespace polyfem::solver
13{
15 {
16 public:
17 CollisionBarrierForm(const VariableToSimulationGroup &variable_to_simulation, const State &state, const double dhat, const double dmin = 0);
18
19 double value_unweighted(const Eigen::VectorXd &x) const override;
20
21 void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
22
23 void solution_changed(const Eigen::VectorXd &x) override;
24
25 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
26
27 std::string name() const override { return "collision barrier"; }
28
29 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
30
31 protected:
32 void build_collision_set(const Eigen::MatrixXd &displaced_surface);
33
34 Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const;
35
36 const State &state_;
37
38 Eigen::VectorXd X_init;
39
40 ipc::CollisionMesh collision_mesh_;
41 ipc::NormalCollisions collision_set;
42 const double dhat_;
43 const double dmin_;
44 ipc::BroadPhaseMethod broad_phase_method_;
45
46 ipc::BarrierPotential barrier_potential_;
47 };
48
50 {
51 public:
52 LayerThicknessForm(const VariableToSimulationGroup &variable_to_simulations,
53 const State &state,
54 const std::vector<int> &boundary_ids,
55 const double dhat,
56 const bool use_log_barrier = false,
57 const double dmin = 0);
58
59 std::string name() const override { return "layer thickness"; }
60
61 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override { return 1.; }
62
63 protected:
65
66 std::vector<int> boundary_ids_;
67 std::map<int, std::set<int>> boundary_ids_to_dof_;
68 Eigen::MatrixXi can_collide_cache_;
69 };
70
72 {
73 public:
74 DeformedCollisionBarrierForm(const VariableToSimulationGroup &variable_to_simulation, const State &state, const double dhat);
75
76 std::string name() const override { return "deformed_collision_barrier"; }
77
78 double value_unweighted(const Eigen::VectorXd &x) const override;
79
80 void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
81
82 void solution_changed(const Eigen::VectorXd &x) override;
83
84 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
85
86 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
87
88 private:
89 void build_collision_set(const Eigen::MatrixXd &displaced_surface);
90
91 Eigen::VectorXd get_updated_mesh_nodes(const Eigen::VectorXd &x) const;
92
93 const State &state_;
94
95 Eigen::VectorXd X_init;
96
97 ipc::CollisionMesh collision_mesh_;
98 ipc::NormalCollisions collision_set;
99 const double dhat_;
100 ipc::BroadPhaseMethod broad_phase_method_;
101
102 const ipc::BarrierPotential barrier_potential_;
103 };
104
106 {
107 public:
109 const VariableToSimulationGroup &variable_to_simulations,
110 const State &state,
111 const json &args);
113
114 double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override;
115 Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override;
116 void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
117 void solution_changed_step(const int time_step, const Eigen::VectorXd &x) override;
118
119 protected:
121 ipc::SmoothCollisions get_smooth_collision_set(const Eigen::MatrixXd &displaced_surface);
122
123 const State &state_;
124 std::set<int> boundary_ids_;
125 std::map<int, std::set<int>> boundary_ids_to_dof_;
126
127 ipc::CollisionMesh collision_mesh_;
128 ipc::SmoothCollisions collisions_;
129 const ipc::ParameterType params_;
130 const double dmin_ = 0;
131
132 ipc::SmoothContactPotential potential_;
133 };
134} // 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
ipc::SmoothCollisions get_smooth_collision_set(const Eigen::MatrixXd &displaced_surface)
void solution_changed_step(const int time_step, const Eigen::VectorXd &x) override
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) 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
std::map< int, std::set< int > > boundary_ids_to_dof_
ipc::SmoothContactPotential potential_
A collection of VariableToSimulation.
nlohmann::json json
Definition Common.hpp:9