Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GenericProblem.hpp
Go to the documentation of this file.
1#pragma once
2
6
7namespace polyfem
8{
9 namespace assembler
10 {
12 {
13 std::array<utils::ExpressionValue, 3> value;
14 std::vector<std::shared_ptr<utils::Interpolation>> interpolation;
15 Eigen::Matrix<bool, 1, 3> dirichlet_dimension;
16
17 void set_unit_type(const std::string &unit_type)
18 {
19 for (auto &v : value)
20 v.set_unit_type(unit_type);
21 }
22
23 double eval(const RowVectorNd &pts, const int dim, const double t, const int el_id = -1) const;
24 };
25
27 {
29 std::shared_ptr<utils::Interpolation> interpolation;
30
31 void set_unit_type(const std::string &unit_type)
32 {
33 value.set_unit_type(unit_type);
34 }
35
36 double eval(const RowVectorNd &pts, const double t) const;
37 };
38
40 {
41 public:
42 GenericTensorProblem(const std::string &name);
43 void set_units(const assembler::Assembler &assembler, const Units &units) override;
44
45 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
46 bool is_rhs_zero() const override
47 {
48 for (int i = 0; i < 3; ++i)
49 {
50 if (!rhs_[i].is_zero())
51 return false;
52 }
53 return true;
54 }
55
56 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;
57 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;
58 void pressure_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;
59 double pressure_cavity_bc(const int boundary_id, const double t) const override;
60
61 void dirichlet_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const double t, Eigen::MatrixXd &val) const override;
62 void neumann_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const Eigen::MatrixXd &normal, const double t, Eigen::MatrixXd &val) const override;
63 bool is_nodal_dirichlet_boundary(const int n_id, const int tag) override;
64 bool is_nodal_neumann_boundary(const int n_id, const int tag) override;
65 bool has_nodal_dirichlet() override;
66 bool has_nodal_neumann() override;
67 bool is_nodal_dimension_dirichlet(const int n_id, const int tag, const int dim) const override;
68 void update_nodes(const Eigen::VectorXi &in_node_to_node) override;
69
70 bool has_exact_sol() const override { return has_exact_; }
71 bool is_scalar() const override { return false; }
72 bool is_time_dependent() const override { return is_time_dept_; }
73 void set_time_dependent(const bool val) { is_time_dept_ = val; }
74 bool is_constant_in_time() const override { return !is_time_dept_; }
75 bool might_have_no_dirichlet() override { return !is_all_; }
76
77 void initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
78 void initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
79 void initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
80
81 void set_parameters(const json &params) override;
82
83 bool is_dimension_dirichet(const int tag, const int dim) const override;
84 bool all_dimensions_dirichlet() const override { return all_dimensions_dirichlet_; }
85
86 void exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
87 void exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
88
89 void clear() override;
90
91 void update_dirichlet_boundary(const int id, const int time_step, const Eigen::VectorXd &val);
92 void update_dirichlet_nodes(const Eigen::VectorXi &in_node_to_node, const Eigen::VectorXi &node_ids, const Eigen::MatrixXd &nodal_dirichlet);
93 void update_pressure_boundary(const int id, const int time_step, const double val);
94
95 private:
97 bool has_exact_ = false;
98 bool has_exact_grad_ = false;
99 bool is_time_dept_ = false;
100 // bool is_mixed_ = false;
101
102 std::vector<TensorBCValue> forces_;
103 std::vector<TensorBCValue> displacements_;
104 std::vector<ScalarBCValue> normal_aligned_forces_;
105 std::vector<ScalarBCValue> pressures_;
106 std::unordered_map<int, ScalarBCValue> cavity_pressures_;
107
108 std::vector<std::pair<int, std::array<utils::ExpressionValue, 3>>> initial_position_;
109 std::vector<std::pair<int, std::array<utils::ExpressionValue, 3>>> initial_velocity_;
110 std::vector<std::pair<int, std::array<utils::ExpressionValue, 3>>> initial_acceleration_;
111
112 std::array<utils::ExpressionValue, 3> rhs_;
113 std::array<utils::ExpressionValue, 3> exact_;
114 std::array<utils::ExpressionValue, 9> exact_grad_;
115
116 std::map<int, TensorBCValue> nodal_dirichlet_;
117 std::map<int, TensorBCValue> nodal_neumann_;
118 std::vector<Eigen::MatrixXd> nodal_dirichlet_mat_;
119
121 };
122
124 {
125 public:
126 GenericScalarProblem(const std::string &name);
127 void set_units(const assembler::Assembler &assembler, const Units &units) override;
128
129 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
130 bool is_rhs_zero() const override { return rhs_.is_zero(); }
131
132 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;
133 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;
134 void initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override;
135
136 void dirichlet_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const double t, Eigen::MatrixXd &val) const override;
137 void neumann_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const Eigen::MatrixXd &normal, const double t, Eigen::MatrixXd &val) const override;
138 bool is_nodal_dirichlet_boundary(const int n_id, const int tag) override;
139 bool is_nodal_neumann_boundary(const int n_id, const int tag) override;
140 bool has_nodal_dirichlet() override;
141 bool has_nodal_neumann() override;
142 void update_nodes(const Eigen::VectorXi &in_node_to_node) override;
143
144 bool has_exact_sol() const override { return has_exact_; }
145 bool is_scalar() const override { return true; }
146 bool is_time_dependent() const override { return is_time_dept_; }
147 void set_time_dependent(const bool val) { is_time_dept_ = val; }
148 bool is_constant_in_time() const override { return !is_time_dept_; }
149 bool might_have_no_dirichlet() override { return !is_all_; }
150
151 void set_parameters(const json &params) override;
152
153 void exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
154 void exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
155
156 void clear() override;
157
158 private:
159 std::vector<ScalarBCValue> neumann_;
160 std::vector<ScalarBCValue> dirichlet_;
161 std::vector<std::pair<int, utils::ExpressionValue>> initial_solution_;
162
163 std::map<int, ScalarBCValue> nodal_dirichlet_;
164 std::map<int, ScalarBCValue> nodal_neumann_;
165 std::vector<Eigen::MatrixXd> nodal_dirichlet_mat_;
166
169 std::array<utils::ExpressionValue, 3> exact_grad_;
171 bool has_exact_ = false;
172 bool has_exact_grad_ = false;
173 bool is_time_dept_ = false;
174 };
175 } // namespace assembler
176} // namespace polyfem
double val
Definition Assembler.cpp:86
bool is_nodal_neumann_boundary(const int n_id, const int tag) 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_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
std::vector< std::pair< int, utils::ExpressionValue > > initial_solution_
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
std::vector< ScalarBCValue > neumann_
void update_nodes(const Eigen::VectorXi &in_node_to_node) override
bool is_nodal_dirichlet_boundary(const int n_id, const int tag) override
void dirichlet_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const double t, Eigen::MatrixXd &val) const override
std::array< utils::ExpressionValue, 3 > exact_grad_
void exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
std::map< int, ScalarBCValue > nodal_dirichlet_
void set_units(const assembler::Assembler &assembler, const Units &units) override
std::map< int, ScalarBCValue > nodal_neumann_
void neumann_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const Eigen::MatrixXd &normal, 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
std::vector< ScalarBCValue > dirichlet_
void exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
std::vector< Eigen::MatrixXd > nodal_dirichlet_mat_
void update_pressure_boundary(const int id, const int time_step, const double val)
std::map< int, TensorBCValue > nodal_neumann_
std::unordered_map< int, ScalarBCValue > cavity_pressures_
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
std::vector< std::pair< int, std::array< utils::ExpressionValue, 3 > > > initial_position_
void set_parameters(const json &params) override
void update_nodes(const Eigen::VectorXi &in_node_to_node) override
void initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const override
std::vector< std::pair< int, std::array< utils::ExpressionValue, 3 > > > initial_acceleration_
std::vector< std::pair< int, std::array< utils::ExpressionValue, 3 > > > initial_velocity_
void update_dirichlet_nodes(const Eigen::VectorXi &in_node_to_node, const Eigen::VectorXi &node_ids, const Eigen::MatrixXd &nodal_dirichlet)
std::vector< TensorBCValue > displacements_
std::vector< ScalarBCValue > normal_aligned_forces_
bool is_nodal_dirichlet_boundary(const int n_id, const int tag) override
void neumann_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const Eigen::MatrixXd &normal, 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
void pressure_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
bool is_dimension_dirichet(const int tag, const int dim) const override
void exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
std::map< int, TensorBCValue > nodal_dirichlet_
std::array< utils::ExpressionValue, 3 > rhs_
std::array< utils::ExpressionValue, 9 > exact_grad_
void update_dirichlet_boundary(const int id, const int time_step, const Eigen::VectorXd &val)
double pressure_cavity_bc(const int boundary_id, const double t) const override
void dirichlet_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, 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
std::vector< ScalarBCValue > pressures_
bool is_nodal_neumann_boundary(const int n_id, const int tag) override
void exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
bool is_nodal_dimension_dirichlet(const int n_id, const int tag, const int dim) const override
void set_units(const assembler::Assembler &assembler, const Units &units) override
void initial_solution(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
std::array< utils::ExpressionValue, 3 > exact_
std::vector< Eigen::MatrixXd > nodal_dirichlet_mat_
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
void set_unit_type(const std::string &unit_type)
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
std::shared_ptr< utils::Interpolation > interpolation
void set_unit_type(const std::string &unit_type)
double eval(const RowVectorNd &pts, const double t) const
Eigen::Matrix< bool, 1, 3 > dirichlet_dimension
std::vector< std::shared_ptr< utils::Interpolation > > interpolation
void set_unit_type(const std::string &unit_type)
std::array< utils::ExpressionValue, 3 > value
double eval(const RowVectorNd &pts, const int dim, const double t, const int el_id=-1) const