PolyFEM
Loading...
Searching...
No Matches
ElasticProblem.hpp
Go to the documentation of this file.
1#pragma once
2
5
6#include <vector>
7#include <Eigen/Dense>
8
9namespace polyfem
10{
11 namespace problem
12 {
14 {
15 public:
16 ElasticProblem(const std::string &name);
17
18 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
19 bool is_rhs_zero() const override { return true; }
20
21 void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
22
23 bool has_exact_sol() const override { return false; }
24 bool is_scalar() const override { return false; }
25 };
26
28 {
29 public:
30 TorsionElasticProblem(const std::string &name);
31
32 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
33 bool is_rhs_zero() const override { return true; }
34
35 void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
36
37 bool has_exact_sol() const override { return false; }
38 bool is_scalar() const override { return false; }
39 bool is_constant_in_time() const override { return false; }
40
41 void set_parameters(const json &params) override;
42
43 private:
44 double n_turns_ = 0.5;
48 };
49
51 {
52 public:
53 DoubleTorsionElasticProblem(const std::string &name);
54
55 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
56 bool is_rhs_zero() const override { return true; }
57
58 void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
59
60 void initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
61 void initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
62 void initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
63
64 bool has_exact_sol() const override { return false; }
65 bool is_scalar() const override { return false; }
66 bool is_constant_in_time() const override { return false; }
67 bool is_time_dependent() const override { return true; }
68
69 void set_parameters(const json &params) override;
70
71 private:
72 double angular_v0_ = 0.5;
73 double angular_v1_ = -0.5;
74 std::array<int, 2> coordiante_0_ = {{0, 1}};
75 std::array<int, 2> coordiante_1_ = {{0, 1}};
78 };
79
81 {
82 public:
83 ElasticProblemZeroBC(const std::string &name);
84 bool is_rhs_zero() const override { return false; }
85
86 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
87 void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
88
89 bool has_exact_sol() const override { return false; }
90 bool is_scalar() const override { return false; }
91 };
92
94 {
95 public:
96 ElasticProblemExact(const std::string &name);
97
98 VectorNd eval_fun(const VectorNd &pt, const double t) const override;
99 AutodiffGradPt eval_fun(const AutodiffGradPt &pt, const double t) const override;
100 AutodiffHessianPt eval_fun(const AutodiffHessianPt &pt, const double t) const override;
101
102 bool is_scalar() const override { return false; }
103 };
104
106 {
107 public:
108 CompressionElasticProblemExact(const std::string &name);
109
110 VectorNd eval_fun(const VectorNd &pt, const double t) const override;
111 AutodiffGradPt eval_fun(const AutodiffGradPt &pt, const double t) const override;
112 AutodiffHessianPt eval_fun(const AutodiffHessianPt &pt, const double t) const override;
113
114 bool is_scalar() const override { return false; }
115 };
116
118 {
119 public:
120 QuadraticElasticProblemExact(const std::string &name);
121
122 VectorNd eval_fun(const VectorNd &pt, const double t) const override;
123 AutodiffGradPt eval_fun(const AutodiffGradPt &pt, const double t) const override;
124 AutodiffHessianPt eval_fun(const AutodiffHessianPt &pt, const double t) const override;
125
126 bool is_scalar() const override { return false; }
127 };
128
130 {
131 public:
132 LinearElasticProblemExact(const std::string &name);
133
134 VectorNd eval_fun(const VectorNd &pt, const double t) const override;
135 AutodiffGradPt eval_fun(const AutodiffGradPt &pt, const double t) const override;
136 AutodiffHessianPt eval_fun(const AutodiffHessianPt &pt, const double t) const override;
137
138 bool is_scalar() const override { return false; }
139 };
140
142 {
143 public:
144 GravityProblem(const std::string &name);
145
146 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
147 bool is_rhs_zero() const override { return false; }
148
149 void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
150
151 void initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
152 void initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
153 void initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
154
155 bool has_exact_sol() const override { return false; }
156 bool is_scalar() const override { return false; }
157 bool is_time_dependent() const override { return true; }
158
159 void set_parameters(const json &params) override;
160
161 private:
162 double force_ = 0.1;
163 };
164
166 {
167 public:
168 WalkProblem(const std::string &name);
169
170 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
171 bool is_rhs_zero() const override { return true; }
172
173 void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
174
175 void initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
176 void initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
177 void initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
178
179 bool has_exact_sol() const override { return false; }
180 bool is_scalar() const override { return false; }
181 bool is_time_dependent() const override { return true; }
182 bool is_constant_in_time() const override { return false; }
183 };
184
186 {
187 public:
188 ElasticCantileverExact(const std::string &name);
189
190 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
191 bool is_rhs_zero() const override { return false; }
192
193 void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
194 void neumann_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &normals, const double t, Eigen::MatrixXd &val) const override;
195
196 bool has_exact_sol() const override { return true; }
197 void exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
198 void exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
199 bool is_scalar() const override { return false; }
200
201 void set_parameters(const json &params) override;
202
203 private:
204 VectorNd eval_fun(const VectorNd &pt, const double t) const;
205 AutodiffGradPt eval_fun(const AutodiffGradPt &pt, const double t) const;
206 AutodiffHessianPt eval_fun(const AutodiffHessianPt &pt, const double t) const;
207 int size_for(const Eigen::MatrixXd &pts) const { return is_scalar() ? 1 : pts.cols(); }
208
210 double E;
211 double nu;
212 std::string formulation;
213 double length;
214 double width;
215 };
216 } // namespace problem
217} // namespace polyfem
double val
Definition Assembler.cpp:86
const std::string & name() const
Definition Problem.hpp:23
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
VectorNd eval_fun(const VectorNd &pt, const double t) const override
void initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
void initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
void set_parameters(const json &params) override
void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
int size_for(const Eigen::MatrixXd &pts) const
void set_parameters(const json &params) override
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void neumann_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &normals, const double t, Eigen::MatrixXd &val) const override
VectorNd eval_fun(const VectorNd &pt, const double t) const
void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
VectorNd eval_fun(const VectorNd &pt, const double t) const override
bool is_rhs_zero() const override
bool has_exact_sol() const override
void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void set_parameters(const json &params) override
void initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
bool has_exact_sol() const override
bool is_time_dependent() const override
void initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
VectorNd eval_fun(const VectorNd &pt, const double t) const override
VectorNd eval_fun(const VectorNd &pt, const double t) const override
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void set_parameters(const json &params) override
void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
bool is_constant_in_time() const override
bool is_rhs_zero() const override
void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
bool has_exact_sol() const override
void initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
bool is_scalar() const override
void initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
bool is_time_dependent() const override
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:11
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
Eigen::Matrix< AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffGradPt