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:
21 BodyForm(const int ndof,
22 const int n_pressure_bases,
23 const std::vector<int> &boundary_nodes,
24 const std::vector<mesh::LocalBoundary> &local_boundary,
25 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
26 const int n_boundary_samples,
27 const Eigen::MatrixXd &rhs,
28 const assembler::RhsAssembler &rhs_assembler,
29 const assembler::Density &density,
30 const bool apply_DBC,
31 const bool is_formulation_mixed,
32 const bool is_time_dependent);
33
34 std::string name() const override { return "body"; }
35
36 protected:
40 double value_unweighted(const Eigen::VectorXd &x) const override;
41
45 void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
46
50 void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override;
51
52 public:
56 void update_quantities(const double t, const Eigen::VectorXd &x) override;
57
58 bool get_apply_DBC() { return apply_DBC_; }
59 void set_apply_DBC(const Eigen::VectorXd &x, const bool val) override
60 {
61 if (val != apply_DBC_)
62 {
65 }
66 }
67
74 const int n_verts,
75 const double t,
76 const Eigen::MatrixXd &x,
77 const Eigen::MatrixXd &adjoint,
78 Eigen::VectorXd &term);
79
81 const Eigen::VectorXd &u_prev,
82 const double t,
83 StiffnessMatrix &hessian) const;
84
85 private:
86 double t_;
87 const int ndof_;
89
90 Eigen::MatrixXd x_prev_;
91
92 const std::vector<int> &boundary_nodes_;
93 const std::vector<mesh::LocalBoundary> &local_boundary_;
94 const std::vector<mesh::LocalBoundary> &local_neumann_boundary_;
96
97 const Eigen::MatrixXd &rhs_;
100
103
104 Eigen::MatrixXd current_rhs_;
105
107 void update_current_rhs(const Eigen::VectorXd &x);
108 };
109} // namespace polyfem::solver
double val
Definition Assembler.cpp:86
int x
Form representing body forces.
Definition BodyForm.hpp:15
std::string name() const override
Definition BodyForm.hpp:34
const assembler::RhsAssembler & rhs_assembler_
Reference to the RHS assembler.
Definition BodyForm.hpp:98
void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const override
Compute the second derivative of the value wrt x.
Definition BodyForm.cpp:75
bool is_formulation_mixed_
True if the formulation is mixed.
Definition BodyForm.hpp:102
const std::vector< int > & boundary_nodes_
Definition BodyForm.hpp:92
double t_
Current time.
Definition BodyForm.hpp:86
void update_quantities(const double t, const Eigen::VectorXd &x) override
Update time dependent quantities.
Definition BodyForm.cpp:80
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:115
void set_apply_DBC(const Eigen::VectorXd &x, const bool val) override
Set if the Dirichlet boundary conditions should be enforced.
Definition BodyForm.hpp:59
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:69
const std::vector< mesh::LocalBoundary > & local_neumann_boundary_
Definition BodyForm.hpp:94
const int ndof_
Number of degrees of freedom.
Definition BodyForm.hpp:87
const assembler::Density & density_
Definition BodyForm.hpp:99
Eigen::MatrixXd current_rhs_
Cached RHS for the current time.
Definition BodyForm.hpp:104
const std::vector< mesh::LocalBoundary > & local_boundary_
Definition BodyForm.hpp:93
void update_current_rhs(const Eigen::VectorXd &x)
Update current_rhs.
Definition BodyForm.cpp:87
bool apply_DBC_
If true, set the Dirichlet boundary conditions in the RHS.
Definition BodyForm.hpp:101
const Eigen::MatrixXd & rhs_
static RHS for the current time
Definition BodyForm.hpp:97
Eigen::MatrixXd x_prev_
Cached previous solution.
Definition BodyForm.hpp:90
void hessian_wrt_u_prev(const Eigen::VectorXd &u_prev, const double t, StiffnessMatrix &hessian) const
Definition BodyForm.cpp:341
double value_unweighted(const Eigen::VectorXd &x) const override
Compute the value of the body force form.
Definition BodyForm.cpp:64
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22