PolyFEM
Loading...
Searching...
No Matches
BodyForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Form.hpp"
4
8
10
11namespace polyfem::solver
12{
14 class BodyForm : public Form
15 {
16 public:
20 BodyForm(const int ndof,
21 const int n_pressure_bases,
22 const std::vector<int> &boundary_nodes,
23 const std::vector<mesh::LocalBoundary> &local_boundary,
24 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
25 const int n_boundary_samples,
26 const Eigen::MatrixXd &rhs,
27 const assembler::RhsAssembler &rhs_assembler,
28 const assembler::Density &density,
29 const bool is_formulation_mixed,
30 const bool is_time_dependent);
31
32 std::string name() const override { return "body"; }
33
34 protected:
38 double value_unweighted(const Eigen::VectorXd &x) const override;
39
43 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
44
48 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
49
50 public:
54 void update_quantities(const double t, const Eigen::VectorXd &x) override;
55
62 const int n_verts,
63 const double t,
64 const Eigen::MatrixXd &x,
65 const Eigen::MatrixXd &adjoint,
66 Eigen::VectorXd &term);
67
69 const Eigen::VectorXd &u_prev,
70 const double t,
71 StiffnessMatrix &hessian) const;
72
73 private:
74 double t_;
75 const int ndof_;
77
78 Eigen::MatrixXd x_prev_;
79
80 const std::vector<int> &boundary_nodes_;
81 const std::vector<mesh::LocalBoundary> &local_boundary_;
82 const std::vector<mesh::LocalBoundary> &local_neumann_boundary_;
84
85 const Eigen::MatrixXd &rhs_;
88
90
91 Eigen::MatrixXd current_rhs_;
92
94 void update_current_rhs(const Eigen::VectorXd &x);
95 };
96} // namespace polyfem::solver
int x
Form representing body forces.
Definition BodyForm.hpp:15
std::string name() const override
Definition BodyForm.hpp:32
const assembler::RhsAssembler & rhs_assembler_
Reference to the RHS assembler.
Definition BodyForm.hpp:86
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
Definition BodyForm.cpp:73
bool is_formulation_mixed_
True if the formulation is mixed.
Definition BodyForm.hpp:89
const std::vector< int > & boundary_nodes_
Definition BodyForm.hpp:80
double t_
Current time.
Definition BodyForm.hpp:74
void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time dependent quantities.
Definition BodyForm.cpp:78
void force_shape_derivative(const int n_verts, const double t, const Eigen::MatrixXd &x, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &term)
Compute the derivative of the force wrt vertex positions, then multiply the resulting matrix with adj...
Definition BodyForm.cpp:106
void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
Compute the first derivative of the value wrt x.
Definition BodyForm.cpp:67
const std::vector< mesh::LocalBoundary > & local_neumann_boundary_
Definition BodyForm.hpp:82
const int ndof_
Number of degrees of freedom.
Definition BodyForm.hpp:75
const assembler::Density & density_
Definition BodyForm.hpp:87
Eigen::MatrixXd current_rhs_
Cached RHS for the current time.
Definition BodyForm.hpp:91
const std::vector< mesh::LocalBoundary > & local_boundary_
Definition BodyForm.hpp:81
void update_current_rhs(const Eigen::VectorXd &x)
Update current_rhs.
Definition BodyForm.cpp:85
const Eigen::MatrixXd & rhs_
static RHS for the current time
Definition BodyForm.hpp:85
Eigen::MatrixXd x_prev_
Cached previous solution.
Definition BodyForm.hpp:78
void hessian_wrt_u_prev(const Eigen::VectorXd &u_prev, const double t, StiffnessMatrix &hessian) const
Definition BodyForm.cpp:332
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the body force form.
Definition BodyForm.cpp:62
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22