PolyFEM
Loading...
Searching...
No Matches
RhsAssembler.hpp
Go to the documentation of this file.
1#pragma once
2
5
9
10namespace polyfem
11{
12 namespace assembler
13 {
14 // computes the rhs of a problem by \int \phi rho rhs
16 {
17 public:
18 // initialization with assembler factory mesh
19 // size of the problem, bases
20 // and solver used internally
22 const Assembler &assembler,
23 const mesh::Mesh &mesh,
24 const mesh::Obstacle &obstacle,
25 const std::vector<int> &dirichlet_nodes,
26 const std::vector<int> &neumann_nodes,
27 const std::vector<RowVectorNd> &dirichlet_nodes_position,
28 const std::vector<RowVectorNd> &neumann_nodes_position,
29 const int n_basis,
30 const int size,
31 const std::vector<basis::ElementBases> &bases,
32 const std::vector<basis::ElementBases> &gbases,
34 const Problem &problem,
35 const std::string bc_method,
36 const json &solver_params);
37
38 // computes the rhs of a problem by \int \phi rho rhs
39 void assemble(const Density &density, Eigen::MatrixXd &rhs, const double t = 1) const;
40
41 // computes the initial soltion for time dependent, calls time_bc
42 void initial_solution(Eigen::MatrixXd &sol) const;
43 // computes the initial velocity for time dependent, calls time_bc
44 void initial_velocity(Eigen::MatrixXd &sol) const;
45 // computes the initial acceleration for time dependent, calls time_bc
46 void initial_acceleration(Eigen::MatrixXd &sol) const;
47
48 // sets boundary conditions to rhs, the boundary conditions are projected (Dirichlet) integrated (Neumann) at resolution
49 // local boundary stores the mapping from elemment to nodes for Dirichlet nodes
50 // local local_neumann_boundary stores the mapping from elemment to nodes for Neumann nodes
51 // calls set_bc
52 void set_bc(
53 const std::vector<mesh::LocalBoundary> &local_boundary,
54 const std::vector<int> &bounday_nodes,
55 const int resolution,
56 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
57 Eigen::MatrixXd &rhs,
58 const Eigen::MatrixXd &displacement = Eigen::MatrixXd(),
59 const double t = 1) const;
60
61 // compute body energy
62 double compute_energy(
63 const Eigen::MatrixXd &displacement,
64 const Eigen::MatrixXd &displacement_prev,
65 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
66 const Density &density,
67 const int resolution,
68 const double t) const;
69 // compute body energy gradient, hessian is zero, rhs is a linear function
71 const std::vector<mesh::LocalBoundary> &local_boundary,
72 const std::vector<int> &bounday_nodes,
73 const Density &density,
74 const int resolution,
75 const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
76 const Eigen::MatrixXd &final_rhs,
77 const double t,
78 Eigen::MatrixXd &rhs) const;
79
80 // compute body hessian wrt to previous solution
81 void compute_energy_hess(const std::vector<int> &bounday_nodes, const int resolution, const std::vector<mesh::LocalBoundary> &local_neumann_boundary, const Eigen::MatrixXd &displacement, const double t, const bool project_to_psd, StiffnessMatrix &hess) const;
82
83 inline const Problem &problem() const { return problem_; }
84 inline const mesh::Mesh &mesh() const { return mesh_; }
85 inline const std::vector<basis::ElementBases> &bases() const { return bases_; }
86 inline const std::vector<basis::ElementBases> &gbases() const { return gbases_; }
87 inline const AssemblyValsCache &ass_vals_cache() const { return ass_vals_cache_; }
88 inline const Assembler &assembler() const { return assembler_; }
89
90 private:
91 // leastsquares fit bc
92 void lsq_bc(const std::function<void(const Eigen::MatrixXi &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &df,
93 const std::vector<mesh::LocalBoundary> &local_boundary, const std::vector<int> &bounday_nodes, const int resolution, Eigen::MatrixXd &rhs) const;
94
95 // sample bc at nodes
96 void sample_bc(const std::function<void(const Eigen::MatrixXi &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &df,
97 const std::vector<mesh::LocalBoundary> &local_boundary, const std::vector<int> &bounday_nodes, Eigen::MatrixXd &rhs) const;
98
99 // set boundary condition
100 // the 2 lambdas are callback to dirichlet df and neumann nf
101 // diriclet boundary condition are projected on the FEM bases, it inverts a linear system
102 void set_bc(
103 const std::function<void(const Eigen::MatrixXi &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &df,
104 const std::function<void(const Eigen::MatrixXi &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &nf,
105 const std::vector<mesh::LocalBoundary> &local_boundary, const std::vector<int> &bounday_nodes, const int resolution, const std::vector<mesh::LocalBoundary> &local_neumann_boundary,
106 const Eigen::MatrixXd &displacement, const double t, Eigen::MatrixXd &rhs) const;
107
108 // sets the time (initial) boundary condition
109 // the lambda depeneds if soltuion, velocity, or acceleration
110 // they are projected on the FEM bases, it inverts a linear system
111 void time_bc(const std::function<void(const mesh::Mesh &, const Eigen::MatrixXi &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &fun, Eigen::MatrixXd &sol) const;
112
116 const int n_basis_;
117 const int size_;
118 const std::vector<basis::ElementBases> &bases_;
119 const std::vector<basis::ElementBases> &gbases_;
122 const std::string bc_method_;
124 const std::vector<int> &dirichlet_nodes_;
125 const std::vector<RowVectorNd> &dirichlet_nodes_position_;
126 const std::vector<int> &neumann_nodes_;
127 const std::vector<RowVectorNd> &neumann_nodes_position_;
128 };
129 } // namespace assembler
130} // namespace polyfem
Caches basis evaluation and geometric mapping at every element.
void set_bc(const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bounday_nodes, const int resolution, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, Eigen::MatrixXd &rhs, const Eigen::MatrixXd &displacement=Eigen::MatrixXd(), const double t=1) const
void time_bc(const std::function< void(const mesh::Mesh &, const Eigen::MatrixXi &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &fun, Eigen::MatrixXd &sol) const
const AssemblyValsCache & ass_vals_cache() const
const mesh::Mesh & mesh() const
const std::vector< RowVectorNd > & neumann_nodes_position_
void set_bc(const std::function< void(const Eigen::MatrixXi &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &df, const std::function< void(const Eigen::MatrixXi &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &nf, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bounday_nodes, const int resolution, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const Eigen::MatrixXd &displacement, const double t, Eigen::MatrixXd &rhs) const
const std::vector< basis::ElementBases > & gbases_
basis functions associated with geometric mapping
const Problem & problem() const
const int size_
dimension of problem
const std::vector< int > & neumann_nodes_
void compute_energy_grad(const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bounday_nodes, const Density &density, const int resolution, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const Eigen::MatrixXd &final_rhs, const double t, Eigen::MatrixXd &rhs) const
const std::vector< RowVectorNd > & dirichlet_nodes_position_
double compute_energy(const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const Density &density, const int resolution, const double t) const
void compute_energy_hess(const std::vector< int > &bounday_nodes, const int resolution, const std::vector< mesh::LocalBoundary > &local_neumann_boundary, const Eigen::MatrixXd &displacement, const double t, const bool project_to_psd, StiffnessMatrix &hess) const
const Assembler & assembler() const
const AssemblyValsCache & ass_vals_cache_
void initial_solution(Eigen::MatrixXd &sol) const
const mesh::Obstacle & obstacle_
void lsq_bc(const std::function< void(const Eigen::MatrixXi &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &df, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bounday_nodes, const int resolution, Eigen::MatrixXd &rhs) const
const std::vector< basis::ElementBases > & bases_
basis functions associated with solution
const std::vector< basis::ElementBases > & gbases() const
void initial_acceleration(Eigen::MatrixXd &sol) const
void initial_velocity(Eigen::MatrixXd &sol) const
void assemble(const Density &density, Eigen::MatrixXd &rhs, const double t=1) const
const std::vector< int > & dirichlet_nodes_
const std::vector< basis::ElementBases > & bases() const
void sample_bc(const std::function< void(const Eigen::MatrixXi &, const Eigen::MatrixXd &, const Eigen::MatrixXd &, Eigen::MatrixXd &)> &df, const std::vector< mesh::LocalBoundary > &local_boundary, const std::vector< int > &bounday_nodes, Eigen::MatrixXd &rhs) const
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
nlohmann::json json
Definition Common.hpp:9
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22