PolyFEM
Loading...
Searching...
No Matches
InversionBarrierForm.hpp
Go to the documentation of this file.
1
2#pragma once
3
4#include "Form.hpp"
5
7
8namespace polyfem::solver
9{
11 {
12 public:
14 const Eigen::MatrixXd &rest_positions, const Eigen::MatrixXi &elements, const int dim, const double vhat);
15
16 std::string name() const override { return "inversion_barrier"; }
17
18 protected:
22 double value_unweighted(const Eigen::VectorXd &x) const override;
23
27 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
28
32 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
33
34 public:
39 bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override;
40
41 protected:
42 static double element_volume(const Eigen::MatrixXd &element_vertices);
43 static Eigen::VectorXd element_volume_gradient(const Eigen::MatrixXd &element_vertices);
44 static Eigen::MatrixXd element_volume_hessian(const Eigen::MatrixXd &element_vertices);
45
46 private:
47 Eigen::MatrixXd rest_positions_;
48 Eigen::MatrixXi elements_;
49 int dim_;
50 double vhat_;
51 };
52} // namespace polyfem::solver
int x
bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) const override
Determine if a step from solution x0 to solution x1 is allowed.
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
static double element_volume(const Eigen::MatrixXd &element_vertices)
static Eigen::MatrixXd element_volume_hessian(const Eigen::MatrixXd &element_vertices)
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the form.
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
static Eigen::VectorXd element_volume_gradient(const Eigen::MatrixXd &element_vertices)
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22