PolyFEM
Loading...
Searching...
No Matches
PointBasedProblem.cpp
Go to the documentation of this file.
3
4#include <iostream>
5#include <string>
6
7namespace polyfem
8{
9 using namespace io;
10 using namespace utils;
11
12 namespace problem
13 {
15 {
16 Eigen::Matrix<bool, 3, 1> dd;
17 dd.setConstant(true);
18 bool all_dimensions_dirichlet = true;
19
20 if (data.is_array())
21 {
22 assert(data.size() == 3);
23 init((double)data[0], (double)data[1], (double)data[2], dd);
24 // TODO add dimension
25 }
26 else if (data.is_object())
27 {
28 Eigen::MatrixXd fun, pts;
29 Eigen::MatrixXi tri;
30
31 std::string rbf = "multiquadric";
32 double eps = 0.1;
33
34 read_matrix(data["function"], fun);
35 read_matrix(data["points"], pts);
36
37 int coord = -1;
38 if (data.contains("coordinate"))
39 coord = data["coordinate"];
40
41 if (coord >= 0)
42 pts = pts.block(0, 0, pts.rows(), 2).eval();
43
44 is_tri = data.contains("triangles");
45 if (is_tri)
46 read_matrix(data["triangles"], tri);
47 else
48 {
49 if (data.contains("rbf"))
50 {
51 const std::string tmp = data["rbf"];
52 rbf = tmp;
53 }
54 if (data.contains("epsilon"))
55 eps = data["epsilon"];
56 }
57
58 if (data.contains("dimension"))
59 {
61 auto &tmp = data["dimension"];
62 assert(tmp.is_array());
63 for (size_t k = 0; k < tmp.size(); ++k)
64 dd(k) = tmp[k];
65 }
66
67 if (is_tri)
68 init(pts, tri, fun, coord, dd);
69 else
70 init(pts, fun, rbf, eps, coord, dd);
71 }
72 else
73 {
74 init(0, 0, 0, dd);
75 }
76
78 }
79
80 Eigen::RowVector3d PointBasedTensorProblem::BCValue::operator()(const Eigen::RowVector3d &pt) const
81 {
82 if (is_val)
83 {
84 return val.transpose();
85 }
86 Eigen::RowVector3d res;
87
88 if (coordiante_0 >= 0)
89 {
90 Eigen::RowVector2d pt2;
91 pt2 << pt(coordiante_0), pt(coordiante_1);
92
93 if (is_tri)
94 res = tri_func.interpolate(pt2);
95 else
96 res = rbf_func.interpolate(pt2);
97 }
98 else
99 {
100 if (is_tri)
101 res = tri_func.interpolate(pt);
102 else
103 res = rbf_func.interpolate(pt);
104 }
105
106 return res;
107 }
108
110 : Problem(name), rhs_(0), scaling_(1)
111 {
112 translation_.setZero();
113 }
114
115 void PointBasedTensorProblem::rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
116 {
117 val = Eigen::MatrixXd::Constant(pts.rows(), pts.cols(), rhs_);
118 }
119
120 void PointBasedTensorProblem::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
121 {
122 val = Eigen::MatrixXd::Zero(pts.rows(), mesh.dimension());
123
124 for (long i = 0; i < pts.rows(); ++i)
125 {
126 const int id = mesh.get_boundary_id(global_ids(i));
127 const auto &pt3d = (pts.row(i) + translation_.transpose()) / scaling_;
128
129 for (size_t b = 0; b < boundary_ids_.size(); ++b)
130 {
131 if (id == boundary_ids_[b])
132 {
133 val.row(i) = bc_[b](pt3d) * scaling_;
134 }
135 }
136 }
137 }
138
139 void PointBasedTensorProblem::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
140 {
141 val = Eigen::MatrixXd::Zero(pts.rows(), mesh.dimension());
142
143 for (long i = 0; i < pts.rows(); ++i)
144 {
145 const int id = mesh.get_boundary_id(global_ids(i));
146 const auto &pt3d = (pts.row(i) + translation_.transpose()) / scaling_;
147
148 for (size_t b = 0; b < neumann_boundary_ids_.size(); ++b)
149 {
150 if (id == neumann_boundary_ids_[b])
151 {
152 val.row(i) = neumann_bc_[b](pt3d) * scaling_;
153 }
154 }
155 }
156 }
157
158 void PointBasedTensorProblem::add_constant(const int bc_tag, const Eigen::Vector3d &value, const Eigen::Matrix<bool, 3, 1> &dd, const bool is_neumann)
159 {
160 if (is_neumann)
161 {
162 neumann_boundary_ids_.push_back(bc_tag);
163 neumann_bc_.emplace_back();
164 neumann_bc_.back().init(value, dd);
165 }
166 else
167 {
168 all_dimensions_dirichlet_ = all_dimensions_dirichlet_ && dd(0) && dd(1) && dd(2);
169 boundary_ids_.push_back(bc_tag);
170 bc_.emplace_back();
171 bc_.back().init(value, dd);
172 }
173 }
174
175 void PointBasedTensorProblem::add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const Eigen::MatrixXi &tri, const int coord, const Eigen::Matrix<bool, 3, 1> &dd, const bool is_neumann)
176 {
177 if (is_neumann)
178 {
179 neumann_boundary_ids_.push_back(bc_tag);
180 neumann_bc_.emplace_back();
181 neumann_bc_.back().init(pts.block(0, 0, pts.rows(), 2), tri, func, coord, dd);
182 }
183 else
184 {
185 all_dimensions_dirichlet_ = all_dimensions_dirichlet_ && dd(0) && dd(1) && dd(2);
186
187 boundary_ids_.push_back(bc_tag);
188 bc_.emplace_back();
189 bc_.back().init(pts.block(0, 0, pts.rows(), 2), tri, func, coord, dd);
190 }
191 }
192
193 void PointBasedTensorProblem::add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const std::string &rbf, const double eps, const int coord, const Eigen::Matrix<bool, 3, 1> &dd, const bool is_neumann)
194 {
195 if (is_neumann)
196 {
197 neumann_boundary_ids_.push_back(bc_tag);
198 neumann_bc_.emplace_back();
199 if (coord >= 0)
200 neumann_bc_.back().init(pts.block(0, 0, pts.rows(), 2), func, rbf, eps, coord, dd);
201 else
202 neumann_bc_.back().init(pts, func, rbf, eps, coord, dd);
203 }
204 else
205 {
206 all_dimensions_dirichlet_ = all_dimensions_dirichlet_ && dd(0) && dd(1) && dd(2);
207
208 boundary_ids_.push_back(bc_tag);
209 bc_.emplace_back();
210 if (coord >= 0)
211 bc_.back().init(pts.block(0, 0, pts.rows(), 2), func, rbf, eps, coord, dd);
212 else
213 bc_.back().init(pts, func, rbf, eps, coord, dd);
214 }
215 }
216
217 bool PointBasedTensorProblem::is_dimension_dirichet(const int tag, const int dim) const
218 {
220 return true;
221
222 for (size_t b = 0; b < boundary_ids_.size(); ++b)
223 {
224 if (tag == boundary_ids_[b])
225 {
226 return bc_[b].is_dirichet_dim(dim);
227 }
228 }
229
230 assert(false);
231 return true;
232 }
233
235 {
236 if (initialized_)
237 return;
238
239 if (params.contains("scaling"))
240 {
241 scaling_ = params["scaling"];
242 }
243
244 if (params.contains("rhs"))
245 {
246 rhs_ = params["rhs"];
247 }
248
249 if (params.contains("translation"))
250 {
251 const auto &j_translation = params["translation"];
252
253 assert(j_translation.is_array());
254 assert(j_translation.size() == 3);
255
256 for (int k = 0; k < 3; ++k)
257 translation_(k) = j_translation[k];
258 }
259
260 if (params.contains("boundary_ids"))
261 {
262 boundary_ids_.clear();
263 auto j_boundary = params["boundary_ids"];
264
265 boundary_ids_.resize(j_boundary.size());
266 bc_.resize(boundary_ids_.size());
267
268 for (size_t i = 0; i < boundary_ids_.size(); ++i)
269 {
270 boundary_ids_[i] = j_boundary[i]["id"];
271 const auto ff = j_boundary[i]["value"];
272 const bool all_d = bc_[i].init(ff);
274 }
275 }
276
277 if (params.contains("neumann_boundary_ids"))
278 {
279 neumann_boundary_ids_.clear();
280 auto j_boundary = params["neumann_boundary_ids"];
281
282 neumann_boundary_ids_.resize(j_boundary.size());
283 neumann_bc_.resize(neumann_boundary_ids_.size());
284
285 for (size_t i = 0; i < neumann_boundary_ids_.size(); ++i)
286 {
287 neumann_boundary_ids_[i] = j_boundary[i]["id"];
288 const auto ff = j_boundary[i]["value"];
289 neumann_bc_[i].init(ff);
290 }
291 }
292 }
293 } // namespace problem
294} // namespace polyfem
double val
Definition Assembler.cpp:86
const std::string & name() const
Definition Problem.hpp:23
std::vector< int > boundary_ids_
Definition Problem.hpp:78
std::vector< int > neumann_boundary_ids_
Definition Problem.hpp:79
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
virtual int get_boundary_id(const int primitive) const
Get the boundary selection of an element (face in 3d, edge in 2d)
Definition Mesh.hpp:475
int dimension() const
utily for dimension
Definition Mesh.hpp:151
Eigen::RowVector3d operator()(const Eigen::RowVector3d &pt) const
void set_parameters(const json &params) override
void add_constant(const int bc_tag, const Eigen::Vector3d &value)
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 add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const Eigen::MatrixXi &tri, const int coord)
bool is_dimension_dirichet(const int tag, const int dim) 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
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
bool read_matrix(const std::string &path, Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &mat)
Reads a matrix from a file. Determines the file format based on the path's extension.
Definition MatrixIO.cpp:18
nlohmann::json json
Definition Common.hpp:9