PolyFEM
Loading...
Searching...
No Matches
GenericProblem.cpp
Go to the documentation of this file.
1#include "GenericProblem.hpp"
2
7
8namespace polyfem
9{
10 using namespace utils;
11
12 namespace assembler
13 {
14 namespace
15 {
16 std::vector<json> flatten_ids(const json &p_j_boundary_tmp)
17 {
18 const std::vector<json> j_boundary_tmp = utils::json_as_array(p_j_boundary_tmp);
19
20 std::vector<json> j_boundary;
21
22 for (size_t i = 0; i < j_boundary_tmp.size(); ++i)
23 {
24 const auto &tmp = j_boundary_tmp[i];
25
26 if (tmp.is_string())
27 {
28 j_boundary.push_back(tmp);
29 continue;
30 }
31 if (!tmp.contains("id"))
32 continue;
33
34 if (tmp["id"].is_array())
35 {
36 for (size_t j = 0; j < tmp["id"].size(); ++j)
37 {
38 json newj = tmp;
39 newj["id"] = tmp["id"][j].get<int>();
40 j_boundary.push_back(newj);
41 }
42 }
43 else
44 j_boundary.push_back(tmp);
45 }
46
47 return j_boundary;
48 }
49 } // namespace
50
51 double TensorBCValue::eval(const RowVectorNd &pts, const int dim, const double t, const int el_id) const
52 {
53 double x = pts(0);
54 double y = pts(1);
55 double z = pts.size() == 2 ? 0 : pts(2);
56
57 double val = value[dim](x, y, z, t, el_id);
58
59 if (interpolation.empty())
60 {
61 }
62 else if (interpolation.size() == 1)
63 val *= interpolation[0]->eval(t);
64 else
65 {
66 assert(dim < interpolation.size());
67 val *= interpolation[dim]->eval(t);
68 }
69
70 return val;
71 }
72
73 double ScalarBCValue::eval(const RowVectorNd &pts, const double t) const
74 {
75 assert(pts.size() == 2 || pts.size() == 3);
76 double x = pts(0), y = pts(1), z = pts.size() == 3 ? pts(2) : 0.0;
77 return value(x, y, z, t) * interpolation->eval(t);
78 }
79
81 : Problem(name), is_all_(false)
82 {
83 }
84
85 void GenericTensorProblem::set_units(const assembler::Assembler &assembler, const Units &units)
86 {
87 if (assembler.is_fluid())
88 {
89 // TODO
90 assert(false);
91 }
92 else
93 {
94 for (int i = 0; i < 3; ++i)
95 {
96 rhs_[i].set_unit_type(units.acceleration());
97 exact_[i].set_unit_type(units.length());
98 }
99 for (int i = 0; i < 3; ++i)
100 exact_grad_[i].set_unit_type("");
101
102 for (auto &v : displacements_)
103 v.set_unit_type(units.length());
104
105 for (auto &v : forces_)
106 v.set_unit_type(units.force());
107
108 for (auto &v : normal_aligned_forces_)
109 v.set_unit_type(units.pressure());
110
111 for (auto &v : pressures_)
112 v.set_unit_type(units.pressure());
113
114 for (auto &v : cavity_pressures_)
115 v.second.set_unit_type(units.pressure());
116
117 for (auto &v : initial_position_)
118 for (int i = 0; i < 3; ++i)
119 v.second[i].set_unit_type(units.length());
120
121 for (auto &v : initial_velocity_)
122 for (int i = 0; i < 3; ++i)
123 v.second[i].set_unit_type(units.velocity());
124
125 for (auto &v : initial_acceleration_)
126 for (int i = 0; i < 3; ++i)
127 v.second[i].set_unit_type(units.acceleration());
128
129 for (auto &v : nodal_dirichlet_)
130 v.second.set_unit_type(units.length());
131
132 for (auto &v : nodal_neumann_)
133 v.second.set_unit_type(units.force());
134 }
135 }
136
137 void GenericTensorProblem::rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
138 {
139 val.resize(pts.rows(), pts.cols());
140
141 if (is_rhs_zero())
142 {
143 val.setZero();
144 return;
145 }
146
147 const bool planar = pts.cols() == 2;
148 for (int i = 0; i < pts.rows(); ++i)
149 {
150 for (int j = 0; j < pts.cols(); ++j)
151 {
152 double x = pts(i, 0), y = pts(i, 1), z = planar ? 0 : pts(i, 2);
153 val(i, j) = rhs_[j](x, y, z, t);
154 }
155 }
156 }
157
158 bool GenericTensorProblem::is_dimension_dirichet(const int tag, const int dim) const
159 {
161 return true;
162
163 if (is_all_)
164 {
165 assert(displacements_.size() == 1);
166 return displacements_[0].dirichlet_dimension[dim];
167 }
168
169 for (size_t b = 0; b < boundary_ids_.size(); ++b)
170 {
171 if (tag == boundary_ids_[b])
172 {
173 auto &tmp = displacements_[b].dirichlet_dimension;
174 return tmp[dim];
175 }
176 }
177
178 assert(false);
179 return true;
180 }
181
182 void GenericTensorProblem::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
183 {
184 val = Eigen::MatrixXd::Zero(pts.rows(), mesh.dimension());
185
186 for (long i = 0; i < pts.rows(); ++i)
187 {
188 if (is_all_)
189 {
190 assert(displacements_.size() == 1);
191 for (int d = 0; d < val.cols(); ++d)
192 {
193 val(i, d) = displacements_[0].eval(pts.row(i), d, t);
194 }
195 }
196 else
197 {
198 const int id = mesh.get_boundary_id(global_ids(i));
199 for (size_t b = 0; b < boundary_ids_.size(); ++b)
200 {
201 if (id == boundary_ids_[b])
202 {
203 for (int d = 0; d < val.cols(); ++d)
204 {
205 val(i, d) = displacements_[b].eval(pts.row(i), d, t);
206 }
207
208 break;
209 }
210 }
211 }
212 }
213 }
214
215 void GenericTensorProblem::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
216 {
217 val = Eigen::MatrixXd::Zero(pts.rows(), mesh.dimension());
218
219 for (long i = 0; i < pts.rows(); ++i)
220 {
221 const int id = mesh.get_boundary_id(global_ids(i));
222
223 for (size_t b = 0; b < neumann_boundary_ids_.size(); ++b)
224 {
225 if (id == neumann_boundary_ids_[b])
226 {
227 for (int d = 0; d < val.cols(); ++d)
228 {
229 val(i, d) = forces_[b].eval(pts.row(i), d, t);
230 }
231
232 break;
233 }
234 }
235
236 for (size_t b = 0; b < normal_aligned_neumann_boundary_ids_.size(); ++b)
237 {
239 {
240 for (int d = 0; d < val.cols(); ++d)
241 {
242 val(i, d) = normal_aligned_forces_[b].eval(pts.row(i), t) * normals(i, d);
243 }
244 break;
245 }
246 }
247 }
248 }
249
250 void GenericTensorProblem::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
251 {
252 val = Eigen::MatrixXd::Zero(pts.rows(), 1);
253
254 for (long i = 0; i < pts.rows(); ++i)
255 {
256 const int id = mesh.get_boundary_id(global_ids(i));
257
258 for (size_t b = 0; b < pressure_boundary_ids_.size(); ++b)
259 {
260 if (id == pressure_boundary_ids_[b])
261 {
262 val(i) = pressures_[b].eval(pts.row(i), t);
263 break;
264 }
265 }
266 }
267 }
268
269 double GenericTensorProblem::pressure_cavity_bc(const int boundary_id, const double t) const
270 {
271 Eigen::VectorXd pt;
272 pt.setZero(3);
273 return cavity_pressures_.at(boundary_id).eval(pt, t);
274 }
275
276 void GenericTensorProblem::exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
277 {
278 assert(has_exact_sol());
279 const bool planar = pts.cols() == 2;
280 val.resize(pts.rows(), pts.cols());
281
282 for (int i = 0; i < pts.rows(); ++i)
283 {
284 for (int j = 0; j < pts.cols(); ++j)
285 {
286 double x = pts(i, 0), y = pts(i, 1), z = planar ? 0 : pts(i, 2);
287 val(i, j) = exact_[j](x, y, z, t);
288 }
289 }
290 }
291
292 void GenericTensorProblem::exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
293 {
294 const int size = pts.cols();
295 val.resize(pts.rows(), pts.cols() * size);
296 if (!has_exact_grad_)
297 return;
298
299 const bool planar = size == 2;
300 for (int i = 0; i < pts.rows(); ++i)
301 {
302 for (int j = 0; j < pts.cols() * size; ++j)
303 {
304 double x = pts(i, 0), y = pts(i, 1), z = planar ? 0 : pts(i, 2);
305 val(i, j) = exact_grad_[j](x, y, z, t);
306 }
307 }
308 }
309
310 void GenericTensorProblem::add_dirichlet_boundary(const int id, const Eigen::RowVector3d &val, const bool isx, const bool isy, const bool isz, const std::shared_ptr<Interpolation> &interp)
311 {
312 boundary_ids_.push_back(id);
313
314 displacements_.emplace_back();
315 for (size_t k = 0; k < val.size(); ++k)
316 displacements_.back().value[k].init(val[k]);
317
318 displacements_.back().dirichlet_dimension << isx, isy, isz;
319 displacements_.back().interpolation.push_back(interp);
320
321 if (!isx || !isy || !isz)
323 }
324
325 void GenericTensorProblem::update_dirichlet_boundary(const int id, const Eigen::RowVector3d &val, const bool isx, const bool isy, const bool isz, const std::shared_ptr<Interpolation> &interp)
326 {
327 int index = -1;
328 for (int i = 0; i < boundary_ids_.size(); ++i)
329 {
330 if (boundary_ids_[i] == id)
331 {
332 index = i;
333 break;
334 }
335 }
336 if (index == -1)
337 {
338 throw "Invalid boundary id";
339 }
340
341 displacements_[index].interpolation.clear();
342 displacements_[index].interpolation.push_back(interp);
343
344 for (size_t k = 0; k < val.size(); ++k)
345 displacements_[index].value[k].init(val[k]);
346
347 displacements_[index].dirichlet_dimension << isx, isy, isz;
348
349 if (!isx || !isy || !isz)
351 }
352
353 void GenericTensorProblem::add_neumann_boundary(const int id, const Eigen::RowVector3d &val, const std::shared_ptr<Interpolation> &interp)
354 {
355 neumann_boundary_ids_.push_back(id);
356 forces_.emplace_back();
357 for (size_t k = 0; k < val.size(); ++k)
358 forces_.back().value[k].init(val[k]);
359
360 forces_.back().interpolation.push_back(interp);
361 }
362
363 void GenericTensorProblem::update_neumann_boundary(const int id, const Eigen::RowVector3d &val, const std::shared_ptr<Interpolation> &interp)
364 {
365 int index = -1;
366 for (int i = 0; i < neumann_boundary_ids_.size(); ++i)
367 {
368 if (neumann_boundary_ids_[i] == id)
369 {
370 index = i;
371 break;
372 }
373 }
374 if (index == -1)
375 {
376 throw "Invalid boundary id";
377 }
378
379 forces_[index].interpolation.clear();
380 forces_[index].interpolation.push_back(interp);
381 for (size_t k = 0; k < val.size(); ++k)
382 forces_[index].value[k].init(val[k]);
383 }
384
385 void GenericTensorProblem::add_pressure_boundary(const int id, const double val, const std::shared_ptr<Interpolation> &interp)
386 {
387 pressure_boundary_ids_.push_back(id);
388 pressures_.emplace_back();
389 pressures_.back().value.init(val);
390 pressures_.back().interpolation = interp;
391 }
392
393 void GenericTensorProblem::update_pressure_boundary(const int id, const double val, const std::shared_ptr<Interpolation> &interp)
394 {
395 int index = -1;
396 for (int i = 0; i < pressure_boundary_ids_.size(); ++i)
397 {
398 if (pressure_boundary_ids_[i] == id)
399 {
400 index = i;
401 break;
402 }
403 }
404 if (index == -1)
405 {
406 throw "Invalid boundary id";
407 }
408 pressures_[index].interpolation = interp;
409 pressures_[index].value.init(val);
410 }
411
412 void GenericTensorProblem::add_dirichlet_boundary(const int id, const std::function<Eigen::MatrixXd(double x, double y, double z, double t)> &func, const bool isx, const bool isy, const bool isz, const std::shared_ptr<Interpolation> &interp)
413 {
414 boundary_ids_.push_back(id);
415 displacements_.emplace_back();
416 displacements_.back().interpolation.push_back(interp);
417 for (size_t k = 0; k < displacements_.back().value.size(); ++k)
418 displacements_.back().value[k].init(func, k);
419
420 displacements_.back().dirichlet_dimension << isx, isy, isz;
421
422 if (!isx || !isy || !isz)
424 }
425
426 void GenericTensorProblem::update_dirichlet_boundary(const int id, const std::function<Eigen::MatrixXd(double x, double y, double z, double t)> &func, const bool isx, const bool isy, const bool isz, const std::shared_ptr<Interpolation> &interp)
427 {
428 int index = -1;
429 for (int i = 0; i < boundary_ids_.size(); ++i)
430 {
431 if (boundary_ids_[i] == id)
432 {
433 index = i;
434 break;
435 }
436 }
437 if (index == -1)
438 {
439 throw "Invalid boundary id";
440 }
441 displacements_[index].interpolation.clear();
442 displacements_[index].interpolation.push_back(interp);
443 for (size_t k = 0; k < displacements_.back().value.size(); ++k)
444 displacements_[index].value[k].init(func, k);
445
446 displacements_[index].dirichlet_dimension << isx, isy, isz;
447
448 if (!isx || !isy || !isz)
450 }
451
452 void GenericTensorProblem::add_neumann_boundary(const int id, const std::function<Eigen::MatrixXd(double x, double y, double z, double t)> &func, const std::shared_ptr<Interpolation> &interp)
453 {
454 neumann_boundary_ids_.push_back(id);
455 forces_.emplace_back();
456 forces_.back().interpolation.push_back(interp);
457 for (size_t k = 0; k < forces_.back().value.size(); ++k)
458 forces_.back().value[k].init(func, k);
459 }
460
461 void GenericTensorProblem::update_neumann_boundary(const int id, const std::function<Eigen::MatrixXd(double x, double y, double z, double t)> &func, const std::shared_ptr<Interpolation> &interp)
462 {
463 int index = -1;
464 for (int i = 0; i < neumann_boundary_ids_.size(); ++i)
465 {
466 if (neumann_boundary_ids_[i] == id)
467 {
468 index = i;
469 break;
470 }
471 }
472 if (index == -1)
473 {
474 throw "Invalid boundary id";
475 }
476
477 forces_[index].interpolation.clear();
478 forces_[index].interpolation.push_back(interp);
479 for (size_t k = 0; k < forces_.back().value.size(); ++k)
480 forces_[index].value[k].init(func, k);
481 }
482
483 void GenericTensorProblem::add_pressure_boundary(const int id, const std::function<double(double x, double y, double z, double t)> &func, const std::shared_ptr<Interpolation> &interp)
484 {
485 pressure_boundary_ids_.push_back(id);
486 pressures_.emplace_back();
487 pressures_.back().value.init(func);
488 pressures_.back().interpolation = interp;
489 }
490
491 void GenericTensorProblem::update_pressure_boundary(const int id, const std::function<double(double x, double y, double z, double t)> &func, const std::shared_ptr<Interpolation> &interp)
492 {
493 int index = -1;
494 for (int i = 0; i < pressure_boundary_ids_.size(); ++i)
495 {
496 if (pressure_boundary_ids_[i] == id)
497 {
498 index = i;
499 break;
500 }
501 }
502 if (index == -1)
503 {
504 throw "Invalid boundary id";
505 }
506
507 pressures_[index].interpolation = interp;
508 pressures_[index].value.init(func);
509 }
510
511 void GenericTensorProblem::add_dirichlet_boundary(const int id, const json &val, const bool isx, const bool isy, const bool isz, const std::string &interpolation)
512 {
513 if (!val.is_array())
514 throw "Val must be an array";
515
516 boundary_ids_.push_back(id);
517 displacements_.emplace_back();
518 if (!interpolation.empty())
519 displacements_.back().interpolation.push_back(Interpolation::build(interpolation));
520
521 for (size_t k = 0; k < val.size(); ++k)
522 displacements_.back().value[k].init(val[k]);
523
524 displacements_.back().dirichlet_dimension << isx, isy, isz;
525
526 if (!isx || !isy || !isz)
528 }
529
530 void GenericTensorProblem::add_neumann_boundary(const int id, const json &val, const std::string &interpolation)
531 {
532 if (!val.is_array())
533 throw "Val must be an array";
534
535 neumann_boundary_ids_.push_back(id);
536
537 forces_.emplace_back();
538 if (!interpolation.empty())
539 forces_.back().interpolation.push_back(Interpolation::build(interpolation));
540
541 for (size_t k = 0; k < val.size(); ++k)
542 forces_.back().value[k].init(val[k]);
543 }
544
545 void GenericTensorProblem::add_pressure_boundary(const int id, json val, const std::string &interpolation)
546 {
547 pressure_boundary_ids_.push_back(id);
548 pressures_.emplace_back();
549
550 if (interpolation.empty())
551 pressures_.back().interpolation = std::make_shared<NoInterpolation>();
552 else
553 pressures_.back().interpolation = Interpolation::build(interpolation);
554 pressures_.back().value.init(val);
555 }
556
557 void GenericTensorProblem::update_dirichlet_boundary(const int id, const json &val, const bool isx, const bool isy, const bool isz, const std::string &interpolation)
558 {
559 if (!val.is_array())
560 throw "Val must be an array";
561 int index = -1;
562 for (int i = 0; i < boundary_ids_.size(); ++i)
563 {
564 if (boundary_ids_[i] == id)
565 {
566 index = i;
567 break;
568 }
569 }
570 if (index == -1)
571 {
572 throw "Invalid boundary id";
573 }
574
575 displacements_[index].interpolation.clear();
576 if (!interpolation.empty())
577 displacements_[index].interpolation.push_back(Interpolation::build(interpolation));
578
579 for (size_t k = 0; k < val.size(); ++k)
580 displacements_[index].value[k].init(val[k]);
581
582 displacements_[index].dirichlet_dimension << isx, isy, isz;
583
584 if (!isx || !isy || !isz)
586 }
587
588 void GenericTensorProblem::update_neumann_boundary(const int id, const json &val, const std::string &interpolation)
589 {
590 if (!val.is_array())
591 throw "Val must be an array";
592
593 int index = -1;
594 for (int i = 0; i < neumann_boundary_ids_.size(); ++i)
595 {
596 if (neumann_boundary_ids_[i] == id)
597 {
598 index = i;
599 break;
600 }
601 }
602 if (index == -1)
603 {
604 throw "Invalid boundary id";
605 }
606 forces_[index].interpolation.clear();
607 if (!interpolation.empty())
608 forces_[index].interpolation.push_back(Interpolation::build(interpolation));
609
610 for (size_t k = 0; k < val.size(); ++k)
611 forces_[index].value[k].init(val[k]);
612 }
613
614 void GenericTensorProblem::update_pressure_boundary(const int id, json val, const std::string &interpolation)
615 {
616 int index = -1;
617 for (int i = 0; i < pressure_boundary_ids_.size(); ++i)
618 {
619 if (pressure_boundary_ids_[i] == id)
620 {
621 index = i;
622 break;
623 }
624 }
625 if (index == -1)
626 {
627 throw "Invalid boundary id";
628 }
629
630 if (interpolation.empty())
631 pressures_[index].interpolation = std::make_shared<NoInterpolation>();
632 else
633 pressures_[index].interpolation = Interpolation::build(interpolation);
634 pressures_[index].value.init(val);
635 }
636
637 void GenericTensorProblem::set_rhs(double x, double y, double z)
638 {
639 rhs_[0].init(x);
640 rhs_[1].init(y);
641 rhs_[2].init(z);
642 }
643
644 void GenericTensorProblem::dirichlet_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const double t, Eigen::MatrixXd &val) const
645 {
646 val = Eigen::MatrixXd::Zero(1, mesh.dimension());
647 const int tag = mesh.get_node_id(node_id);
648
649 if (is_all_)
650 {
651 assert(nodal_dirichlet_.size() == 1);
652 const auto &tmp = nodal_dirichlet_.begin()->second;
653
654 for (int d = 0; d < val.cols(); ++d)
655 {
656 val(d) = tmp.eval(pt, d, t);
657 }
658
659 return;
660 }
661
662 const auto it = nodal_dirichlet_.find(tag);
663 if (it != nodal_dirichlet_.end())
664 {
665 for (int d = 0; d < val.cols(); ++d)
666 {
667 val(d) = it->second.eval(pt, d, t);
668 }
669
670 return;
671 }
672
673 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
674 {
675 for (int i = 0; i < n_dirichlet.rows(); ++i)
676 {
677 if (n_dirichlet(i, 0) == node_id)
678 {
679 for (int d = 0; d < val.cols(); ++d)
680 {
681 val(d) = n_dirichlet(i, d + 1);
682 }
683
684 return;
685 }
686 }
687 }
688
689 assert(false);
690 }
691
692 void GenericTensorProblem::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
693 {
694 // TODO implement me;
695 log_and_throw_error("Nodal neumann not implemented");
696 }
697
698 bool GenericTensorProblem::is_nodal_dirichlet_boundary(const int n_id, const int tag)
699 {
700 if (nodal_dirichlet_.find(tag) != nodal_dirichlet_.end())
701 return true;
702
703 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
704 {
705 for (int i = 0; i < n_dirichlet.rows(); ++i)
706 {
707 if (n_dirichlet(i, 0) == n_id)
708 return true;
709 }
710 }
711
712 return false;
713 }
714
715 bool GenericTensorProblem::is_nodal_neumann_boundary(const int n_id, const int tag)
716 {
717 return nodal_neumann_.find(tag) != nodal_neumann_.end();
718 }
719
721 {
722 return nodal_dirichlet_mat_.size() > 0;
723 }
724
726 {
727 return false; // nodal_neumann_.size() > 0;
728 }
729
730 bool GenericTensorProblem::is_nodal_dimension_dirichlet(const int n_id, const int tag, const int dim) const
731 {
732 const auto it = nodal_dirichlet_.find(tag);
733 if (it != nodal_dirichlet_.end())
734 {
735 return it->second.dirichlet_dimension(dim);
736 }
737
738 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
739 {
740 for (int i = 0; i < n_dirichlet.rows(); ++i)
741 {
742 if (n_dirichlet(i, 0) == n_id)
743 {
744 return !std::isnan(n_dirichlet(i, dim + 1));
745 }
746 }
747 }
748
749 assert(false);
750 return true;
751 }
752
753 void GenericTensorProblem::update_nodes(const Eigen::VectorXi &in_node_to_node)
754 {
755 for (auto &n_dirichlet : nodal_dirichlet_mat_)
756 {
757 for (int n = 0; n < n_dirichlet.rows(); ++n)
758 {
759 const int node_id = in_node_to_node[n_dirichlet(n, 0)];
760 n_dirichlet(n, 0) = node_id;
761 }
762 }
763 }
764
766 {
767 if (is_param_valid(params, "is_time_dependent"))
768 {
769 is_time_dept_ = params["is_time_dependent"];
770 }
771
772 if (is_param_valid(params, "rhs"))
773 {
774 auto rr = params["rhs"];
775 if (rr.is_array())
776 {
777 for (size_t k = 0; k < rr.size(); ++k)
778 rhs_[k].init(rr[k]);
779 }
780 else
781 {
782 logger().warn("Invalid problem rhs: should be an array.");
783 assert(false);
784 }
785 }
786
787 if (is_param_valid(params, "reference") && is_param_valid(params["reference"], "solution"))
788 {
789 auto ex = params["reference"]["solution"];
790 has_exact_ = ex.size() > 0;
791 if (ex.is_array())
792 {
793 for (size_t k = 0; k < ex.size(); ++k)
794 exact_[k].init(ex[k]);
795 }
796 else
797 {
798 assert(false);
799 }
800 }
801
802 if (is_param_valid(params, "reference") && is_param_valid(params["reference"], "gradient"))
803 {
804 auto ex = params["reference"]["gradient"];
805 has_exact_grad_ = ex.size() > 0;
806 if (ex.is_array())
807 {
808 for (size_t k = 0; k < ex.size(); ++k)
809 exact_grad_[k].init(ex[k]);
810 }
811 else
812 {
813 assert(false);
814 }
815 }
816
817 if (is_param_valid(params, "dirichlet_boundary"))
818 {
819 // boundary_ids_.clear();
820 int offset = boundary_ids_.size();
821 std::vector<json> j_boundary = flatten_ids(params["dirichlet_boundary"]);
822
823 boundary_ids_.resize(offset + j_boundary.size());
824 displacements_.resize(offset + j_boundary.size());
825
826 for (size_t i = offset; i < boundary_ids_.size(); ++i)
827 {
828 if (j_boundary[i - offset].is_string())
829 {
830 const std::string path = resolve_path(j_boundary[i - offset], params["root_path"]);
831 if (!std::filesystem::is_regular_file(path))
832 log_and_throw_error("unable to open {} file", path);
833
834 Eigen::MatrixXd tmp;
835 io::read_matrix(path, tmp);
836 nodal_dirichlet_mat_.emplace_back(tmp);
837
838 continue;
839 }
840
841 int current_id = -1;
842
843 if (j_boundary[i - offset]["id"] == "all")
844 {
845 assert(boundary_ids_.size() == 1);
846 boundary_ids_.clear();
847 is_all_ = true;
848 nodal_dirichlet_[current_id] = TensorBCValue();
849 }
850 else
851 {
852 boundary_ids_[i] = j_boundary[i - offset]["id"];
853 current_id = boundary_ids_[i];
854 nodal_dirichlet_[current_id] = TensorBCValue();
855 }
856
857 auto ff = j_boundary[i - offset]["value"];
858 if (ff.is_array())
859 {
860 for (size_t k = 0; k < ff.size(); ++k)
861 {
862 displacements_[i].value[k].init(ff[k]);
863 if (j_boundary[i - offset].contains("time_reference") && j_boundary[i - offset]["time_reference"].size() > 0)
864 displacements_[i].value[k].set_t(j_boundary[i - offset]["time_reference"]);
865 nodal_dirichlet_[current_id].value[k].init(ff[k]);
866 }
867 }
868 else
869 {
870 assert(false);
871 displacements_[i].value[0].init(0);
872 displacements_[i].value[1].init(0);
873 displacements_[i].value[2].init(0);
874 }
875
876 displacements_[i].dirichlet_dimension.setConstant(true);
877 nodal_dirichlet_[current_id].dirichlet_dimension.setConstant(true);
878 if (j_boundary[i - offset].contains("dimension"))
879 {
881 auto &tmp = j_boundary[i - offset]["dimension"];
882 assert(tmp.is_array());
883 for (size_t k = 0; k < tmp.size(); ++k)
884 {
885 displacements_[i].dirichlet_dimension[k] = tmp[k];
886 nodal_dirichlet_[current_id].dirichlet_dimension[k] = tmp[k];
887 }
888 }
889
890 if (j_boundary[i - offset]["interpolation"].is_array())
891 {
892 for (int ii = 0; ii < j_boundary[i - offset]["interpolation"].size(); ++ii)
893 displacements_[i].interpolation.push_back(Interpolation::build(j_boundary[i - offset]["interpolation"][ii]));
894 }
895 else
896 displacements_[i].interpolation.push_back(Interpolation::build(j_boundary[i - offset]["interpolation"]));
897
898 nodal_dirichlet_[current_id].interpolation = displacements_[i].interpolation;
899 }
900 }
901
902 if (is_param_valid(params, "neumann_boundary"))
903 {
904 // neumann_boundary_ids_.clear();
905 const int offset = neumann_boundary_ids_.size();
906
907 auto j_boundary_tmp = params["neumann_boundary"];
908 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
909
910 neumann_boundary_ids_.resize(offset + j_boundary.size());
911 forces_.resize(offset + j_boundary.size());
912
913 for (size_t i = offset; i < neumann_boundary_ids_.size(); ++i)
914 {
915 neumann_boundary_ids_[i] = j_boundary[i - offset]["id"];
916
917 auto ff = j_boundary[i - offset]["value"];
918 assert(ff.is_array());
919
920 for (size_t k = 0; k < ff.size(); ++k)
921 forces_[i].value[k].init(ff[k]);
922
923 if (j_boundary[i - offset]["interpolation"].is_array())
924 {
925 for (int ii = 0; ii < j_boundary[i - offset]["interpolation"].size(); ++ii)
926 forces_[i].interpolation.push_back(Interpolation::build(j_boundary[i - offset]["interpolation"][ii]));
927 }
928 else
929 {
930 forces_[i].interpolation.push_back(Interpolation::build(j_boundary[i - offset]["interpolation"]));
931 }
932 }
933 }
934
935 if (is_param_valid(params, "normal_aligned_neumann_boundary"))
936 {
937 const int offset = normal_aligned_neumann_boundary_ids_.size();
938
939 auto j_boundary_tmp = params["normal_aligned_neumann_boundary"];
940 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
941
942 normal_aligned_neumann_boundary_ids_.resize(offset + j_boundary.size());
943 normal_aligned_forces_.resize(offset + j_boundary.size());
944
945 for (size_t i = offset; i < normal_aligned_neumann_boundary_ids_.size(); ++i)
946 {
947 normal_aligned_neumann_boundary_ids_[i] = j_boundary[i - offset]["id"];
948
949 auto ff = j_boundary[i - offset]["value"];
950 normal_aligned_forces_[i].value.init(ff);
951
952 if (j_boundary[i - offset].contains("interpolation"))
953 normal_aligned_forces_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"]);
954 else
955 normal_aligned_forces_[i].interpolation = std::make_shared<NoInterpolation>();
956 }
957 }
958
959 if (is_param_valid(params, "pressure_boundary"))
960 {
961 // pressure_boundary_ids_.clear();
962 const int offset = pressure_boundary_ids_.size();
963
964 auto j_boundary_tmp = params["pressure_boundary"];
965 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
966
967 pressure_boundary_ids_.resize(offset + j_boundary.size());
968 pressures_.resize(offset + j_boundary.size());
969
970 for (size_t i = offset; i < pressure_boundary_ids_.size(); ++i)
971 {
972 pressure_boundary_ids_[i] = j_boundary[i - offset]["id"];
973
974 auto ff = j_boundary[i - offset]["value"];
975 pressures_[i].value.init(ff);
976
977 if (j_boundary[i - offset].contains("interpolation"))
978 pressures_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"]);
979 else
980 pressures_[i].interpolation = std::make_shared<NoInterpolation>();
981 }
982 }
983
984 if (is_param_valid(params, "pressure_cavity"))
985 {
986 const int offset = pressure_cavity_ids_.size();
987
988 auto j_boundary_tmp = params["pressure_cavity"];
989 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
990
991 pressure_cavity_ids_.resize(offset + j_boundary.size());
992
993 for (size_t i = offset; i < pressure_cavity_ids_.size(); ++i)
994 {
995 int boundary_id = j_boundary[i - offset]["id"];
996 pressure_cavity_ids_[i] = boundary_id;
997
998 if (cavity_pressures_.find(boundary_id) == cavity_pressures_.end())
999 {
1000 cavity_pressures_[boundary_id] = ScalarBCValue();
1001
1002 auto ff = j_boundary[i - offset]["value"];
1003 cavity_pressures_[boundary_id].value.init(ff);
1004 }
1005 }
1006 }
1007
1008 if (is_param_valid(params, "solution"))
1009 {
1010 auto rr = params["solution"];
1011 initial_position_.resize(rr.size());
1012 assert(rr.is_array());
1013
1014 for (size_t k = 0; k < rr.size(); ++k)
1015 {
1016 initial_position_[k].first = rr[k]["id"];
1017 const auto v = rr[k]["value"];
1018 for (size_t d = 0; d < v.size(); ++d)
1019 initial_position_[k].second[d].init(v[d]);
1020 }
1021 }
1022
1023 if (is_param_valid(params, "velocity"))
1024 {
1025 auto rr = params["velocity"];
1026 initial_velocity_.resize(rr.size());
1027 assert(rr.is_array());
1028
1029 for (size_t k = 0; k < rr.size(); ++k)
1030 {
1031 initial_velocity_[k].first = rr[k]["id"];
1032 const auto v = rr[k]["value"];
1033 for (size_t d = 0; d < v.size(); ++d)
1034 initial_velocity_[k].second[d].init(v[d]);
1035 }
1036 }
1037
1038 if (is_param_valid(params, "acceleration"))
1039 {
1040 auto rr = params["acceleration"];
1041 initial_acceleration_.resize(rr.size());
1042 assert(rr.is_array());
1043
1044 for (size_t k = 0; k < rr.size(); ++k)
1045 {
1046 initial_acceleration_[k].first = rr[k]["id"];
1047 const auto v = rr[k]["value"];
1048 for (size_t d = 0; d < v.size(); ++d)
1049 initial_acceleration_[k].second[d].init(v[d]);
1050 }
1051 }
1052 }
1053
1054 void GenericTensorProblem::initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const
1055 {
1056 val.resize(pts.rows(), pts.cols());
1057 if (initial_position_.empty())
1058 {
1059 val.setZero();
1060 return;
1061 }
1062
1063 const bool planar = pts.cols() == 2;
1064 for (int i = 0; i < pts.rows(); ++i)
1065 {
1066 const int id = mesh.get_body_id(global_ids(i));
1067 int index = -1;
1068 for (int j = 0; j < initial_position_.size(); ++j)
1069 {
1070 if (initial_position_[j].first == id)
1071 {
1072 index = j;
1073 break;
1074 }
1075 }
1076 if (index < 0)
1077 {
1078 val.row(i).setZero();
1079 continue;
1080 }
1081
1082 for (int j = 0; j < pts.cols(); ++j)
1083 val(i, j) = planar ? initial_position_[index].second[j](pts(i, 0), pts(i, 1)) : initial_position_[index].second[j](pts(i, 0), pts(i, 1), pts(i, 2));
1084 }
1085 }
1086
1087 void GenericTensorProblem::initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const
1088 {
1089 val.resize(pts.rows(), pts.cols());
1090 if (initial_velocity_.empty())
1091 {
1092 val.setZero();
1093 return;
1094 }
1095
1096 const bool planar = pts.cols() == 2;
1097 for (int i = 0; i < pts.rows(); ++i)
1098 {
1099 const int id = mesh.get_body_id(global_ids(i));
1100 int index = -1;
1101 for (int j = 0; j < initial_velocity_.size(); ++j)
1102 {
1103 if (initial_velocity_[j].first == id)
1104 {
1105 index = j;
1106 break;
1107 }
1108 }
1109 if (index < 0)
1110 {
1111 val.row(i).setZero();
1112 continue;
1113 }
1114
1115 for (int j = 0; j < pts.cols(); ++j)
1116 val(i, j) = planar ? initial_velocity_[index].second[j](pts(i, 0), pts(i, 1)) : initial_velocity_[index].second[j](pts(i, 0), pts(i, 1), pts(i, 2));
1117 }
1118 }
1119
1120 void GenericTensorProblem::initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const
1121 {
1122 val.resize(pts.rows(), pts.cols());
1123 if (initial_acceleration_.empty())
1124 {
1125 val.setZero();
1126 return;
1127 }
1128
1129 const bool planar = pts.cols() == 2;
1130 for (int i = 0; i < pts.rows(); ++i)
1131 {
1132 const int id = mesh.get_body_id(global_ids(i));
1133 int index = -1;
1134 for (int j = 0; j < initial_acceleration_.size(); ++j)
1135 {
1136 if (initial_acceleration_[j].first == id)
1137 {
1138 index = j;
1139 break;
1140 }
1141 }
1142 if (index < 0)
1143 {
1144 val.row(i).setZero();
1145 continue;
1146 }
1147
1148 for (int j = 0; j < pts.cols(); ++j)
1149 val(i, j) = planar ? initial_acceleration_[index].second[j](pts(i, 0), pts(i, 1)) : initial_acceleration_[index].second[j](pts(i, 0), pts(i, 1), pts(i, 2));
1150 }
1151 }
1152
1154 {
1156 has_exact_ = false;
1157 has_exact_grad_ = false;
1158 is_time_dept_ = false;
1159
1160 forces_.clear();
1161 displacements_.clear();
1162 normal_aligned_forces_.clear();
1163 pressures_.clear();
1164 cavity_pressures_.clear();
1165
1166 nodal_dirichlet_.clear();
1167 nodal_neumann_.clear();
1168 nodal_dirichlet_mat_.clear();
1169
1170 initial_position_.clear();
1171 initial_velocity_.clear();
1172 initial_acceleration_.clear();
1173
1174 for (int i = 0; i < rhs_.size(); ++i)
1175 rhs_[i].clear();
1176 for (int i = 0; i < exact_.size(); ++i)
1177 exact_[i].clear();
1178 for (int i = 0; i < exact_grad_.size(); ++i)
1179 exact_grad_[i].clear();
1180 is_all_ = false;
1181 }
1182
1184 : Problem(name), is_all_(false)
1185 {
1186 }
1187
1189 {
1190 // TODO?
1191
1192 for (auto &v : neumann_)
1193 v.set_unit_type("");
1194
1195 for (auto &v : dirichlet_)
1196 v.set_unit_type("");
1197
1198 for (auto &v : initial_solution_)
1199 v.second.set_unit_type("");
1200
1201 rhs_.set_unit_type("");
1203
1204 for (int i = 0; i < 3; ++i)
1205 exact_grad_[i].set_unit_type("");
1206 }
1207
1208 void GenericScalarProblem::rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
1209 {
1210 val.resize(pts.rows(), 1);
1211 if (is_rhs_zero())
1212 {
1213 val.setZero();
1214 return;
1215 }
1216 const bool planar = pts.cols() == 2;
1217 for (int i = 0; i < pts.rows(); ++i)
1218 {
1219 double x = pts(i, 0), y = pts(i, 1), z = planar ? 0 : pts(i, 2);
1220 val(i) = rhs_(x, y, z, t);
1221 }
1222 }
1223
1224 void GenericScalarProblem::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
1225 {
1226 val = Eigen::MatrixXd::Zero(pts.rows(), 1);
1227
1228 for (long i = 0; i < pts.rows(); ++i)
1229 {
1230 const int id = mesh.get_boundary_id(global_ids(i));
1231 if (is_all_)
1232 {
1233 assert(dirichlet_.size() == 1);
1234 val(i) = dirichlet_[0].eval(pts.row(i), t);
1235 }
1236 else
1237 {
1238 for (size_t b = 0; b < boundary_ids_.size(); ++b)
1239 {
1240 if (id == boundary_ids_[b])
1241 {
1242 val(i) = dirichlet_[b].eval(pts.row(i), t);
1243 break;
1244 }
1245 }
1246 }
1247 }
1248 }
1249
1250 void GenericScalarProblem::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
1251 {
1252 val = Eigen::MatrixXd::Zero(pts.rows(), 1);
1253
1254 for (long i = 0; i < pts.rows(); ++i)
1255 {
1256 const int id = mesh.get_boundary_id(global_ids(i));
1257
1258 for (size_t b = 0; b < neumann_boundary_ids_.size(); ++b)
1259 {
1260 if (id == neumann_boundary_ids_[b])
1261 {
1262 double x = pts(i, 0), y = pts(i, 1), z = pts.cols() == 2 ? 0 : pts(i, 2);
1263 val(i) = neumann_[b].eval(pts.row(i), t);
1264 break;
1265 }
1266 }
1267 }
1268 }
1269
1270 void GenericScalarProblem::initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const
1271 {
1272 val.resize(pts.rows(), 1);
1273 if (initial_solution_.empty())
1274 {
1275 val.setZero();
1276 return;
1277 }
1278
1279 const bool planar = pts.cols() == 2;
1280 for (int i = 0; i < pts.rows(); ++i)
1281 {
1282 const int id = mesh.get_body_id(global_ids(i));
1283 int index = -1;
1284 for (int j = 0; j < initial_solution_.size(); ++j)
1285 {
1286 if (initial_solution_[j].first == id)
1287 {
1288 index = j;
1289 break;
1290 }
1291 }
1292 if (index < 0)
1293 {
1294 val(i) = 0;
1295 continue;
1296 }
1297
1298 val(i) = planar ? initial_solution_[index].second(pts(i, 0), pts(i, 1)) : initial_solution_[index].second(pts(i, 0), pts(i, 1), pts(i, 2));
1299 }
1300 }
1301
1302 void GenericScalarProblem::exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
1303 {
1304 assert(has_exact_sol());
1305 const bool planar = pts.cols() == 2;
1306 val.resize(pts.rows(), 1);
1307
1308 for (int i = 0; i < pts.rows(); ++i)
1309 {
1310 double x = pts(i, 0), y = pts(i, 1), z = pts.cols() == 2 ? 0 : pts(i, 2);
1311 val(i) = exact_(x, y, z, t);
1312 }
1313 }
1314
1315 void GenericScalarProblem::exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
1316 {
1317 val.resize(pts.rows(), pts.cols());
1318 if (!has_exact_grad_)
1319 return;
1320
1321 const bool planar = pts.cols() == 2;
1322 for (int i = 0; i < pts.rows(); ++i)
1323 {
1324 for (int j = 0; j < pts.cols(); ++j)
1325 {
1326 double x = pts(i, 0), y = pts(i, 1), z = pts.cols() == 2 ? 0 : pts(i, 2);
1327 val(i, j) = exact_grad_[j](x, y, z, t);
1328 }
1329 }
1330 }
1331 void GenericScalarProblem::dirichlet_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const double t, Eigen::MatrixXd &val) const
1332 {
1333 val = Eigen::MatrixXd::Zero(1, 1);
1334 const int tag = mesh.get_node_id(node_id);
1335
1336 if (is_all_)
1337 {
1338 assert(nodal_dirichlet_.size() == 1);
1339 const auto &tmp = nodal_dirichlet_.begin()->second;
1340
1341 val(0) = tmp.eval(pt, t);
1342 return;
1343 }
1344
1345 const auto it = nodal_dirichlet_.find(tag);
1346 if (it != nodal_dirichlet_.end())
1347 {
1348 val(0) = it->second.eval(pt, t);
1349 return;
1350 }
1351
1352 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
1353 {
1354 for (int i = 0; i < n_dirichlet.rows(); ++i)
1355 {
1356 if (n_dirichlet(i, 0) == node_id)
1357 {
1358 val(0) = n_dirichlet(i, 1);
1359 return;
1360 }
1361 }
1362 }
1363
1364 assert(false);
1365 }
1366
1367 void GenericScalarProblem::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
1368 {
1369 // TODO implement me;
1370 log_and_throw_error("Nodal neumann not implemented");
1371 }
1372
1373 bool GenericScalarProblem::is_nodal_dirichlet_boundary(const int n_id, const int tag)
1374 {
1375 if (nodal_dirichlet_.find(tag) != nodal_dirichlet_.end())
1376 return true;
1377
1378 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
1379 {
1380 for (int i = 0; i < n_dirichlet.rows(); ++i)
1381 {
1382 if (n_dirichlet(i, 0) == n_id)
1383 return true;
1384 }
1385 }
1386
1387 return false;
1388 }
1389
1390 bool GenericScalarProblem::is_nodal_neumann_boundary(const int n_id, const int tag)
1391 {
1392 return nodal_neumann_.find(tag) != nodal_neumann_.end();
1393 }
1394
1396 {
1397 return nodal_dirichlet_mat_.size() > 0;
1398 }
1399
1401 {
1402 return false; // nodal_neumann_.size() > 0;
1403 }
1404
1405 void GenericScalarProblem::update_nodes(const Eigen::VectorXi &in_node_to_node)
1406 {
1407 for (auto &n_dirichlet : nodal_dirichlet_mat_)
1408 {
1409 for (int n = 0; n < n_dirichlet.rows(); ++n)
1410 {
1411 const int node_id = in_node_to_node[n_dirichlet(n, 0)];
1412 n_dirichlet(n, 0) = node_id;
1413 }
1414 }
1415 }
1416
1418 {
1419 if (is_param_valid(params, "is_time_dependent"))
1420 {
1421 is_time_dept_ = params["is_time_dependent"];
1422 }
1423
1424 if (is_param_valid(params, "rhs"))
1425 {
1426 rhs_.init(params["rhs"]);
1427 }
1428
1429 if (is_param_valid(params, "reference") && is_param_valid(params["reference"], "solution"))
1430 {
1431 has_exact_ = !params["reference"]["solution"].empty();
1432 exact_.init(params["reference"]["solution"]);
1433 }
1434
1435 if (is_param_valid(params, "reference") && is_param_valid(params["reference"], "gradient"))
1436 {
1437 auto ex = params["reference"]["gradient"];
1438 has_exact_grad_ = ex.size() > 0;
1439 if (ex.is_array())
1440 {
1441 for (size_t k = 0; k < ex.size(); ++k)
1442 exact_grad_[k].init(ex[k]);
1443 }
1444 else
1445 {
1446 assert(false);
1447 }
1448 }
1449
1450 if (is_param_valid(params, "dirichlet_boundary"))
1451 {
1452 // boundary_ids_.clear();
1453 const int offset = boundary_ids_.size();
1454 std::vector<json> j_boundary = flatten_ids(params["dirichlet_boundary"]);
1455
1456 boundary_ids_.resize(offset + j_boundary.size());
1457 dirichlet_.resize(offset + j_boundary.size());
1458
1459 for (size_t i = offset; i < boundary_ids_.size(); ++i)
1460 {
1461 if (j_boundary[i - offset].is_string())
1462 {
1463 const std::string path = resolve_path(j_boundary[i - offset], params["root_path"]);
1464 if (!std::filesystem::is_regular_file(path))
1465 log_and_throw_error(fmt::format("unable to open {} file", path));
1466
1467 Eigen::MatrixXd tmp;
1468 io::read_matrix(path, tmp);
1469 nodal_dirichlet_mat_.emplace_back(tmp);
1470
1471 continue;
1472 }
1473
1474 int current_id = -1;
1475
1476 if (j_boundary[i - offset]["id"] == "all")
1477 {
1478 assert(boundary_ids_.size() == 1);
1479
1480 is_all_ = true;
1481 boundary_ids_.clear();
1482 nodal_dirichlet_[current_id] = ScalarBCValue();
1483 }
1484 else
1485 {
1486 boundary_ids_[i] = j_boundary[i - offset]["id"];
1487 current_id = boundary_ids_[i];
1488 nodal_dirichlet_[current_id] = ScalarBCValue();
1489 }
1490
1491 auto ff = j_boundary[i - offset]["value"];
1492 dirichlet_[i].value.init(ff);
1493 nodal_dirichlet_[current_id].value.init(ff);
1494
1495 if (j_boundary[i - offset]["interpolation"].is_array())
1496 {
1497 if (j_boundary[i - offset]["interpolation"].size() == 0)
1498 dirichlet_[i].interpolation = std::make_shared<NoInterpolation>();
1499 else if (j_boundary[i - offset]["interpolation"].size() == 1)
1500 dirichlet_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"][0]);
1501 else
1502 log_and_throw_error("Only one Dirichlet interpolation supported");
1503 }
1504 else
1505 dirichlet_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"]);
1506
1507 nodal_dirichlet_[current_id].interpolation = dirichlet_[i].interpolation;
1508 }
1509 }
1510
1511 if (is_param_valid(params, "neumann_boundary"))
1512 {
1513 // neumann_boundary_ids_.clear();
1514 const int offset = neumann_boundary_ids_.size();
1515 auto j_boundary_tmp = params["neumann_boundary"];
1516 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
1517
1518 neumann_boundary_ids_.resize(offset + j_boundary.size());
1519 neumann_.resize(offset + j_boundary.size());
1520
1521 for (size_t i = offset; i < neumann_boundary_ids_.size(); ++i)
1522 {
1523 neumann_boundary_ids_[i] = j_boundary[i - offset]["id"];
1524
1525 auto ff = j_boundary[i - offset]["value"];
1526 neumann_[i].value.init(ff);
1527
1528 if (j_boundary[i - offset]["interpolation"].is_array())
1529 {
1530 if (j_boundary[i - offset]["interpolation"].size() == 0)
1531 neumann_[i].interpolation = std::make_shared<NoInterpolation>();
1532 else if (j_boundary[i - offset]["interpolation"].size() == 1)
1533 neumann_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"][0]);
1534 else
1535 log_and_throw_error("Only one Neumann interpolation supported");
1536 }
1537 else
1538 neumann_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"]);
1539 }
1540 }
1541
1542 if (is_param_valid(params, "solution"))
1543 {
1544 auto rr = params["solution"];
1545 initial_solution_.resize(rr.size());
1546 assert(rr.is_array());
1547
1548 for (size_t k = 0; k < rr.size(); ++k)
1549 {
1550 initial_solution_[k].first = rr[k]["id"];
1551 initial_solution_[k].second.init(rr[k]["value"]);
1552 }
1553 }
1554 }
1555
1556 void GenericScalarProblem::add_dirichlet_boundary(const int id, const double val, const std::shared_ptr<Interpolation> &interp)
1557 {
1558 boundary_ids_.push_back(id);
1559 dirichlet_.emplace_back();
1560 dirichlet_.back().value.init(val);
1561 dirichlet_.back().interpolation = interp;
1562 }
1563
1564 void GenericScalarProblem::update_dirichlet_boundary(const int id, const double val, const std::shared_ptr<Interpolation> &interp)
1565 {
1566 int index = -1;
1567 for (int i = 0; i < boundary_ids_.size(); ++i)
1568 {
1569 if (boundary_ids_[i] == id)
1570 {
1571 index = i;
1572 break;
1573 }
1574 }
1575 if (index == -1)
1576 {
1577 throw "Invalid boundary id";
1578 }
1579
1580 dirichlet_[index].value.init(val);
1581 dirichlet_[index].interpolation = interp;
1582 }
1583
1584 void GenericScalarProblem::add_neumann_boundary(const int id, const double val, const std::shared_ptr<Interpolation> &interp)
1585 {
1586 neumann_boundary_ids_.push_back(id);
1587 neumann_.emplace_back();
1588 neumann_.back().value.init(val);
1589 neumann_.back().interpolation = interp;
1590 }
1591
1592 void GenericScalarProblem::update_neumann_boundary(const int id, const double val, const std::shared_ptr<Interpolation> &interp)
1593 {
1594 int index = -1;
1595 for (int i = 0; i < neumann_boundary_ids_.size(); ++i)
1596 {
1597 if (neumann_boundary_ids_[i] == id)
1598 {
1599 index = i;
1600 break;
1601 }
1602 }
1603 if (index == -1)
1604 {
1605 throw "Invalid boundary id";
1606 }
1607
1608 neumann_[index].value.init(val);
1609 neumann_[index].interpolation = interp;
1610 }
1611
1612 void GenericScalarProblem::add_dirichlet_boundary(const int id, const std::function<double(double x, double y, double z, double t)> &func, const std::shared_ptr<Interpolation> &interp)
1613 {
1614 boundary_ids_.push_back(id);
1615 dirichlet_.emplace_back();
1616 dirichlet_.back().value.init(func);
1617 dirichlet_.back().interpolation = interp;
1618 }
1619
1620 void GenericScalarProblem::update_dirichlet_boundary(const int id, const std::function<double(double x, double y, double z, double t)> &func, const std::shared_ptr<Interpolation> &interp)
1621 {
1622 int index = -1;
1623 for (int i = 0; i < boundary_ids_.size(); ++i)
1624 {
1625 if (boundary_ids_[i] == id)
1626 {
1627 index = i;
1628 break;
1629 }
1630 }
1631 if (index == -1)
1632 {
1633 throw "Invalid boundary id";
1634 }
1635 dirichlet_[index].value.init(func);
1636 dirichlet_[index].interpolation = interp;
1637 }
1638
1639 void GenericScalarProblem::add_neumann_boundary(const int id, const std::function<double(double x, double y, double z, double t)> &func, const std::shared_ptr<Interpolation> &interp)
1640 {
1641 neumann_boundary_ids_.push_back(id);
1642 neumann_.emplace_back();
1643 neumann_.back().value.init(func);
1644 neumann_.back().interpolation = interp;
1645 }
1646
1647 void GenericScalarProblem::update_neumann_boundary(const int id, const std::function<double(double x, double y, double z, double t)> &func, const std::shared_ptr<Interpolation> &interp)
1648 {
1649 int index = -1;
1650 for (int i = 0; i < neumann_boundary_ids_.size(); ++i)
1651 {
1652 if (neumann_boundary_ids_[i] == id)
1653 {
1654 index = i;
1655 break;
1656 }
1657 }
1658 if (index == -1)
1659 {
1660 throw "Invalid boundary id";
1661 }
1662 neumann_[index].value.init(func);
1663 neumann_[index].interpolation = interp;
1664 }
1665
1666 void GenericScalarProblem::add_dirichlet_boundary(const int id, const json &val, const std::string &interp)
1667 {
1668 boundary_ids_.push_back(id);
1669 dirichlet_.emplace_back();
1670 dirichlet_.back().value.init(val);
1671 if (interp.empty())
1672 dirichlet_.back().interpolation = std::make_shared<NoInterpolation>();
1673 else
1674 dirichlet_.back().interpolation = Interpolation::build(interp);
1675 }
1676
1677 void GenericScalarProblem::add_neumann_boundary(const int id, const json &val, const std::string &interp)
1678 {
1679 neumann_boundary_ids_.push_back(id);
1680 neumann_.emplace_back();
1681 neumann_.back().value.init(val);
1682 if (interp.empty())
1683 neumann_.back().interpolation = std::make_shared<NoInterpolation>();
1684 else
1685 neumann_.back().interpolation = Interpolation::build(interp);
1686 }
1687
1688 void GenericScalarProblem::update_dirichlet_boundary(const int id, const json &val, const std::string &interp)
1689 {
1690 int index = -1;
1691 for (int i = 0; i < boundary_ids_.size(); ++i)
1692 {
1693 if (boundary_ids_[i] == id)
1694 {
1695 index = i;
1696 break;
1697 }
1698 }
1699 if (index == -1)
1700 {
1701 throw "Invalid boundary id";
1702 }
1703
1704 dirichlet_[index].value.init(val);
1705 if (interp.empty())
1706 dirichlet_[index].interpolation = std::make_shared<NoInterpolation>();
1707 else
1708 dirichlet_[index].interpolation = Interpolation::build(interp);
1709 }
1710
1711 void GenericScalarProblem::update_neumann_boundary(const int id, const json &val, const std::string &interp)
1712 {
1713 int index = -1;
1714 for (int i = 0; i < neumann_boundary_ids_.size(); ++i)
1715 {
1716 if (neumann_boundary_ids_[i] == id)
1717 {
1718 index = i;
1719 break;
1720 }
1721 }
1722 if (index == -1)
1723 {
1724 throw "Invalid boundary id";
1725 }
1726 neumann_[index].value.init(val);
1727
1728 if (interp.empty())
1729 neumann_[index].interpolation = std::make_shared<NoInterpolation>();
1730 else
1731 neumann_[index].interpolation = Interpolation::build(interp);
1732 }
1733
1735 {
1736 neumann_.clear();
1737 dirichlet_.clear();
1738
1739 nodal_dirichlet_.clear();
1740 nodal_neumann_.clear();
1741 nodal_dirichlet_mat_.clear();
1742
1743 rhs_.clear();
1744 exact_.clear();
1745 for (int i = 0; i < exact_grad_.size(); ++i)
1746 exact_grad_[i].clear();
1747 is_all_ = false;
1748 has_exact_ = false;
1749 has_exact_grad_ = false;
1750 is_time_dept_ = false;
1751 }
1752 } // namespace assembler
1753} // namespace polyfem
double val
Definition Assembler.cpp:86
int y
int z
int x
std::string force() const
Definition Units.hpp:29
std::string pressure() const
Definition Units.hpp:30
std::string velocity() const
Definition Units.hpp:27
std::string acceleration() const
Definition Units.hpp:28
const std::string & length() const
Definition Units.hpp:19
virtual bool is_fluid() const
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 add_dirichlet_boundary(const int id, const double val, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
void exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void add_neumann_boundary(const int id, const double val, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
std::map< int, ScalarBCValue > nodal_dirichlet_
void update_dirichlet_boundary(const int id, const double val, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
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
void update_neumann_boundary(const int id, const double val, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
std::vector< Eigen::MatrixXd > nodal_dirichlet_mat_
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 add_dirichlet_boundary(const int id, const Eigen::RowVector3d &val, const bool isx, const bool isy, const bool isz, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
void set_parameters(const json &params) override
void update_nodes(const Eigen::VectorXi &in_node_to_node) override
GenericTensorProblem(const std::string &name)
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_
void set_rhs(double x, double y, double z)
void update_pressure_boundary(const int id, const double val, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
std::vector< std::pair< int, std::array< utils::ExpressionValue, 3 > > > initial_velocity_
std::vector< TensorBCValue > displacements_
std::vector< ScalarBCValue > normal_aligned_forces_
bool is_nodal_dirichlet_boundary(const int n_id, const int tag) override
void update_dirichlet_boundary(const int id, const Eigen::RowVector3d &val, const bool isx, const bool isy, const bool isz, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
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 update_neumann_boundary(const int id, const Eigen::RowVector3d &val, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
void add_pressure_boundary(const int id, const double val, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
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 add_neumann_boundary(const int id, const Eigen::RowVector3d &val, const std::shared_ptr< utils::Interpolation > &interp=std::make_shared< utils::NoInterpolation >())
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_
std::vector< int > normal_aligned_neumann_boundary_ids_
Definition Problem.hpp:80
std::vector< int > pressure_boundary_ids_
Definition Problem.hpp:81
std::vector< int > boundary_ids_
Definition Problem.hpp:78
virtual void init(const mesh::Mesh &mesh)
Definition Problem.hpp:22
std::vector< int > neumann_boundary_ids_
Definition Problem.hpp:79
std::vector< int > pressure_cavity_ids_
Definition Problem.hpp:82
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
virtual int get_body_id(const int primitive) const
Get the volume selection of an element (cell in 3d, face in 2d)
Definition Mesh.hpp:504
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:478
int dimension() const
utily for dimension
Definition Mesh.hpp:148
virtual int get_node_id(const int node_id) const
Get the boundary selection of a node.
Definition Mesh.hpp:487
void set_unit_type(const std::string &unit_type)
static std::shared_ptr< Interpolation > build(const json &params)
list tmp
Definition p_bases.py:339
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
std::string resolve_path(const std::string &path, const std::string &input_file_path, const bool only_if_exists=false)
std::vector< T > json_as_array(const json &j)
Return the value of a json object as an array.
Definition JSONUtils.hpp:38
bool is_param_valid(const json &params, const std::string &key)
Determine if a key exists and is non-null in a json object.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:11
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
std::shared_ptr< utils::Interpolation > interpolation
double eval(const RowVectorNd &pts, const double t) const
std::vector< std::shared_ptr< utils::Interpolation > > interpolation
std::array< utils::ExpressionValue, 3 > value
double eval(const RowVectorNd &pts, const int dim, const double t, const int el_id=-1) const