Loading [MathJax]/jax/output/HTML-CSS/config.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 for (int i = 0; i < 3; ++i)
90 {
91 rhs_[i].set_unit_type(units.force());
92 exact_[i].set_unit_type(units.velocity());
93 }
94 for (int i = 0; i < 3; ++i)
95 exact_grad_[i].set_unit_type("");
96
97 for (auto &v : displacements_)
98 v.set_unit_type(units.velocity());
99
100 for (auto &v : forces_)
101 v.set_unit_type(units.force());
102
103 for (auto &v : normal_aligned_forces_)
104 v.set_unit_type(units.pressure());
105
106 for (auto &v : pressures_)
107 v.set_unit_type(units.pressure());
108
109 for (auto &v : cavity_pressures_)
110 v.second.set_unit_type(units.pressure());
111
112 for (auto &v : initial_position_)
113 for (int i = 0; i < 3; ++i)
114 v.second[i].set_unit_type(units.velocity());
115
116 for (auto &v : initial_velocity_)
117 for (int i = 0; i < 3; ++i)
118 v.second[i].set_unit_type(units.velocity());
119
120 for (auto &v : initial_acceleration_)
121 for (int i = 0; i < 3; ++i)
122 v.second[i].set_unit_type(units.acceleration());
123
124 for (auto &v : nodal_dirichlet_)
125 v.second.set_unit_type(units.velocity());
126
127 for (auto &v : nodal_neumann_)
128 v.second.set_unit_type(units.force());
129 }
130 else
131 {
132 for (int i = 0; i < 3; ++i)
133 {
134 rhs_[i].set_unit_type(units.acceleration());
135 exact_[i].set_unit_type(units.length());
136 }
137 for (int i = 0; i < 3; ++i)
138 exact_grad_[i].set_unit_type("");
139
140 for (auto &v : displacements_)
141 v.set_unit_type(units.length());
142
143 for (auto &v : forces_)
144 v.set_unit_type(units.force());
145
146 for (auto &v : normal_aligned_forces_)
147 v.set_unit_type(units.pressure());
148
149 for (auto &v : pressures_)
150 v.set_unit_type(units.pressure());
151
152 for (auto &v : cavity_pressures_)
153 v.second.set_unit_type(units.pressure());
154
155 for (auto &v : initial_position_)
156 for (int i = 0; i < 3; ++i)
157 v.second[i].set_unit_type(units.length());
158
159 for (auto &v : initial_velocity_)
160 for (int i = 0; i < 3; ++i)
161 v.second[i].set_unit_type(units.velocity());
162
163 for (auto &v : initial_acceleration_)
164 for (int i = 0; i < 3; ++i)
165 v.second[i].set_unit_type(units.acceleration());
166
167 for (auto &v : nodal_dirichlet_)
168 v.second.set_unit_type(units.length());
169
170 for (auto &v : nodal_neumann_)
171 v.second.set_unit_type(units.force());
172 }
173 }
174
175 void GenericTensorProblem::rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
176 {
177 val.resize(pts.rows(), pts.cols());
178
179 if (is_rhs_zero())
180 {
181 val.setZero();
182 return;
183 }
184
185 const bool planar = pts.cols() == 2;
186 for (int i = 0; i < pts.rows(); ++i)
187 {
188 for (int j = 0; j < pts.cols(); ++j)
189 {
190 double x = pts(i, 0), y = pts(i, 1), z = planar ? 0 : pts(i, 2);
191 val(i, j) = rhs_[j](x, y, z, t);
192 }
193 }
194 }
195
196 bool GenericTensorProblem::is_dimension_dirichet(const int tag, const int dim) const
197 {
199 return true;
200
201 if (is_all_)
202 {
203 assert(displacements_.size() == 1);
204 return displacements_[0].dirichlet_dimension[dim];
205 }
206
207 for (size_t b = 0; b < boundary_ids_.size(); ++b)
208 {
209 if (tag == boundary_ids_[b])
210 {
211 auto &tmp = displacements_[b].dirichlet_dimension;
212 return tmp[dim];
213 }
214 }
215
216 assert(false);
217 return true;
218 }
219
220 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
221 {
222 val = Eigen::MatrixXd::Zero(pts.rows(), mesh.dimension());
223
224 for (long i = 0; i < pts.rows(); ++i)
225 {
226 if (is_all_)
227 {
228 assert(displacements_.size() == 1);
229 for (int d = 0; d < val.cols(); ++d)
230 {
231 val(i, d) = displacements_[0].eval(pts.row(i), d, t);
232 }
233 }
234 else
235 {
236 const int id = mesh.get_boundary_id(global_ids(i));
237 for (size_t b = 0; b < boundary_ids_.size(); ++b)
238 {
239 if (id == boundary_ids_[b])
240 {
241 for (int d = 0; d < val.cols(); ++d)
242 {
243 val(i, d) = displacements_[b].eval(pts.row(i), d, t);
244 }
245
246 break;
247 }
248 }
249 }
250 }
251 }
252
253 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
254 {
255 val = Eigen::MatrixXd::Zero(pts.rows(), mesh.dimension());
256
257 for (long i = 0; i < pts.rows(); ++i)
258 {
259 const int id = mesh.get_boundary_id(global_ids(i));
260
261 for (size_t b = 0; b < neumann_boundary_ids_.size(); ++b)
262 {
263 if (id == neumann_boundary_ids_[b])
264 {
265 for (int d = 0; d < val.cols(); ++d)
266 {
267 val(i, d) = forces_[b].eval(pts.row(i), d, t);
268 }
269
270 break;
271 }
272 }
273
274 for (size_t b = 0; b < normal_aligned_neumann_boundary_ids_.size(); ++b)
275 {
277 {
278 for (int d = 0; d < val.cols(); ++d)
279 {
280 val(i, d) = normal_aligned_forces_[b].eval(pts.row(i), t) * normals(i, d);
281 }
282 break;
283 }
284 }
285 }
286 }
287
288 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
289 {
290 val = Eigen::MatrixXd::Zero(pts.rows(), 1);
291
292 for (long i = 0; i < pts.rows(); ++i)
293 {
294 const int id = mesh.get_boundary_id(global_ids(i));
295
296 for (size_t b = 0; b < pressure_boundary_ids_.size(); ++b)
297 {
298 if (id == pressure_boundary_ids_[b])
299 {
300 val(i) = pressures_[b].eval(pts.row(i), t);
301 break;
302 }
303 }
304 }
305 }
306
307 double GenericTensorProblem::pressure_cavity_bc(const int boundary_id, const double t) const
308 {
309 Eigen::VectorXd pt;
310 pt.setZero(3);
311 return cavity_pressures_.at(boundary_id).eval(pt, t);
312 }
313
314 void GenericTensorProblem::exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
315 {
316 assert(has_exact_sol());
317 const bool planar = pts.cols() == 2;
318 val.resize(pts.rows(), pts.cols());
319
320 for (int i = 0; i < pts.rows(); ++i)
321 {
322 for (int j = 0; j < pts.cols(); ++j)
323 {
324 double x = pts(i, 0), y = pts(i, 1), z = planar ? 0 : pts(i, 2);
325 val(i, j) = exact_[j](x, y, z, t);
326 }
327 }
328 }
329
330 void GenericTensorProblem::exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
331 {
332 const int size = pts.cols();
333 val.resize(pts.rows(), pts.cols() * size);
334 if (!has_exact_grad_)
335 return;
336
337 const bool planar = size == 2;
338 for (int i = 0; i < pts.rows(); ++i)
339 {
340 for (int j = 0; j < pts.cols() * size; ++j)
341 {
342 double x = pts(i, 0), y = pts(i, 1), z = planar ? 0 : pts(i, 2);
343 val(i, j) = exact_grad_[j](x, y, z, t);
344 }
345 }
346 }
347
348 void GenericTensorProblem::dirichlet_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const double t, Eigen::MatrixXd &val) const
349 {
350 val = Eigen::MatrixXd::Zero(1, mesh.dimension());
351 const int tag = mesh.get_node_id(node_id);
352
353 if (is_all_)
354 {
355 assert(nodal_dirichlet_.size() == 1);
356 const auto &tmp = nodal_dirichlet_.begin()->second;
357
358 for (int d = 0; d < val.cols(); ++d)
359 {
360 val(d) = tmp.eval(pt, d, t);
361 }
362
363 return;
364 }
365
366 const auto it = nodal_dirichlet_.find(tag);
367 if (it != nodal_dirichlet_.end())
368 {
369 for (int d = 0; d < val.cols(); ++d)
370 {
371 val(d) = it->second.eval(pt, d, t);
372 }
373
374 return;
375 }
376
377 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
378 {
379 for (int i = 0; i < n_dirichlet.rows(); ++i)
380 {
381 if (n_dirichlet(i, 0) == node_id)
382 {
383 for (int d = 0; d < val.cols(); ++d)
384 {
385 val(d) = n_dirichlet(i, d + 1);
386 }
387
388 return;
389 }
390 }
391 }
392
393 assert(false);
394 }
395
396 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
397 {
398 // TODO implement me;
399 log_and_throw_error("Nodal neumann not implemented");
400 }
401
402 bool GenericTensorProblem::is_nodal_dirichlet_boundary(const int n_id, const int tag)
403 {
404 if (nodal_dirichlet_.find(tag) != nodal_dirichlet_.end())
405 return true;
406
407 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
408 {
409 for (int i = 0; i < n_dirichlet.rows(); ++i)
410 {
411 if (n_dirichlet(i, 0) == n_id)
412 return true;
413 }
414 }
415
416 return false;
417 }
418
419 bool GenericTensorProblem::is_nodal_neumann_boundary(const int n_id, const int tag)
420 {
421 return nodal_neumann_.find(tag) != nodal_neumann_.end();
422 }
423
425 {
426 return nodal_dirichlet_mat_.size() > 0;
427 }
428
430 {
431 return false; // nodal_neumann_.size() > 0;
432 }
433
434 bool GenericTensorProblem::is_nodal_dimension_dirichlet(const int n_id, const int tag, const int dim) const
435 {
436 const auto it = nodal_dirichlet_.find(tag);
437 if (it != nodal_dirichlet_.end())
438 {
439 return it->second.dirichlet_dimension(dim);
440 }
441
442 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
443 {
444 for (int i = 0; i < n_dirichlet.rows(); ++i)
445 {
446 if (n_dirichlet(i, 0) == n_id)
447 {
448 return !std::isnan(n_dirichlet(i, dim + 1));
449 }
450 }
451 }
452
453 assert(false);
454 return true;
455 }
456
457 void GenericTensorProblem::update_nodes(const Eigen::VectorXi &in_node_to_node)
458 {
460 {
461 if (nodal_dirichlet_mat_.size() > 0)
462 logger().debug("Skipping updating in nodes to nodes in problem, already done once...");
463 return;
464 }
465 for (auto &n_dirichlet : nodal_dirichlet_mat_)
466 {
467 for (int n = 0; n < n_dirichlet.rows(); ++n)
468 {
469 const int node_id = in_node_to_node[n_dirichlet(n, 0)];
470 n_dirichlet(n, 0) = node_id;
471 }
472 }
474 }
475
477 {
478 if (is_param_valid(params, "is_time_dependent"))
479 {
480 is_time_dept_ = params["is_time_dependent"];
481 }
482
483 if (is_param_valid(params, "rhs"))
484 {
485 auto rr = params["rhs"];
486 if (rr.is_array())
487 {
488 for (size_t k = 0; k < rr.size(); ++k)
489 rhs_[k].init(rr[k]);
490 }
491 else
492 {
493 logger().warn("Invalid problem rhs: should be an array.");
494 assert(false);
495 }
496 }
497
498 if (is_param_valid(params, "reference") && is_param_valid(params["reference"], "solution"))
499 {
500 auto ex = params["reference"]["solution"];
501 has_exact_ = ex.size() > 0;
502 if (ex.is_array())
503 {
504 for (size_t k = 0; k < ex.size(); ++k)
505 exact_[k].init(ex[k]);
506 }
507 else
508 {
509 assert(false);
510 }
511 }
512
513 if (is_param_valid(params, "reference") && is_param_valid(params["reference"], "gradient"))
514 {
515 auto ex = params["reference"]["gradient"];
516 has_exact_grad_ = ex.size() > 0;
517 if (ex.is_array())
518 {
519 for (size_t k = 0; k < ex.size(); ++k)
520 exact_grad_[k].init(ex[k]);
521 }
522 else
523 {
524 assert(false);
525 }
526 }
527
528 if (is_param_valid(params, "dirichlet_boundary"))
529 {
530 // boundary_ids_.clear();
531 int offset = boundary_ids_.size();
532 std::vector<json> j_boundary = flatten_ids(params["dirichlet_boundary"]);
533
534 boundary_ids_.resize(offset + j_boundary.size());
535 displacements_.resize(offset + j_boundary.size());
536
537 for (size_t i = offset; i < boundary_ids_.size(); ++i)
538 {
539 if (j_boundary[i - offset].is_string())
540 {
541 const std::string path = resolve_path(j_boundary[i - offset], params["root_path"]);
542 if (!std::filesystem::is_regular_file(path))
543 log_and_throw_error("unable to open {} file", path);
544
545 Eigen::MatrixXd tmp;
546 io::read_matrix(path, tmp);
547 nodal_dirichlet_mat_.emplace_back(tmp);
548
549 continue;
550 }
551
552 int current_id = -1;
553
554 if (j_boundary[i - offset]["id"] == "all")
555 {
556 assert(boundary_ids_.size() == 1);
557 boundary_ids_.clear();
558 is_all_ = true;
559 nodal_dirichlet_[current_id] = TensorBCValue();
560 }
561 else
562 {
563 boundary_ids_[i] = j_boundary[i - offset]["id"];
564 current_id = boundary_ids_[i];
565 nodal_dirichlet_[current_id] = TensorBCValue();
566 }
567
568 auto ff = j_boundary[i - offset]["value"];
569 if (ff.is_array())
570 {
571 for (size_t k = 0; k < ff.size(); ++k)
572 {
573 displacements_[i].value[k].init(ff[k]);
574 if (j_boundary[i - offset].contains("time_reference") && j_boundary[i - offset]["time_reference"].size() > 0)
575 displacements_[i].value[k].set_t(j_boundary[i - offset]["time_reference"]);
576 nodal_dirichlet_[current_id].value[k].init(ff[k]);
577 }
578 }
579 else
580 {
581 assert(false);
582 displacements_[i].value[0].init(0);
583 displacements_[i].value[1].init(0);
584 displacements_[i].value[2].init(0);
585 }
586
587 displacements_[i].dirichlet_dimension.setConstant(true);
588 nodal_dirichlet_[current_id].dirichlet_dimension.setConstant(true);
589 if (j_boundary[i - offset].contains("dimension"))
590 {
592 auto &tmp = j_boundary[i - offset]["dimension"];
593 assert(tmp.is_array());
594 for (size_t k = 0; k < tmp.size(); ++k)
595 {
596 displacements_[i].dirichlet_dimension[k] = tmp[k];
597 nodal_dirichlet_[current_id].dirichlet_dimension[k] = tmp[k];
598 }
599 }
600
601 if (j_boundary[i - offset]["interpolation"].is_array())
602 {
603 for (int ii = 0; ii < j_boundary[i - offset]["interpolation"].size(); ++ii)
604 displacements_[i].interpolation.push_back(Interpolation::build(j_boundary[i - offset]["interpolation"][ii]));
605 }
606 else
607 displacements_[i].interpolation.push_back(Interpolation::build(j_boundary[i - offset]["interpolation"]));
608
609 nodal_dirichlet_[current_id].interpolation = displacements_[i].interpolation;
610 }
611 }
612
613 if (is_param_valid(params, "neumann_boundary"))
614 {
615 // neumann_boundary_ids_.clear();
616 const int offset = neumann_boundary_ids_.size();
617
618 auto j_boundary_tmp = params["neumann_boundary"];
619 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
620
621 neumann_boundary_ids_.resize(offset + j_boundary.size());
622 forces_.resize(offset + j_boundary.size());
623
624 for (size_t i = offset; i < neumann_boundary_ids_.size(); ++i)
625 {
626 neumann_boundary_ids_[i] = j_boundary[i - offset]["id"];
627
628 auto ff = j_boundary[i - offset]["value"];
629 assert(ff.is_array());
630
631 for (size_t k = 0; k < ff.size(); ++k)
632 forces_[i].value[k].init(ff[k]);
633
634 if (j_boundary[i - offset]["interpolation"].is_array())
635 {
636 for (int ii = 0; ii < j_boundary[i - offset]["interpolation"].size(); ++ii)
637 forces_[i].interpolation.push_back(Interpolation::build(j_boundary[i - offset]["interpolation"][ii]));
638 }
639 else
640 {
641 forces_[i].interpolation.push_back(Interpolation::build(j_boundary[i - offset]["interpolation"]));
642 }
643 }
644 }
645
646 if (is_param_valid(params, "normal_aligned_neumann_boundary"))
647 {
648 const int offset = normal_aligned_neumann_boundary_ids_.size();
649
650 auto j_boundary_tmp = params["normal_aligned_neumann_boundary"];
651 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
652
653 normal_aligned_neumann_boundary_ids_.resize(offset + j_boundary.size());
654 normal_aligned_forces_.resize(offset + j_boundary.size());
655
656 for (size_t i = offset; i < normal_aligned_neumann_boundary_ids_.size(); ++i)
657 {
658 normal_aligned_neumann_boundary_ids_[i] = j_boundary[i - offset]["id"];
659
660 auto ff = j_boundary[i - offset]["value"];
661 normal_aligned_forces_[i].value.init(ff);
662
663 if (j_boundary[i - offset].contains("interpolation"))
664 normal_aligned_forces_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"]);
665 else
666 normal_aligned_forces_[i].interpolation = std::make_shared<NoInterpolation>();
667 }
668 }
669
670 if (is_param_valid(params, "pressure_boundary"))
671 {
672 // pressure_boundary_ids_.clear();
673 const int offset = pressure_boundary_ids_.size();
674
675 auto j_boundary_tmp = params["pressure_boundary"];
676 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
677
678 pressure_boundary_ids_.resize(offset + j_boundary.size());
679 pressures_.resize(offset + j_boundary.size());
680
681 for (size_t i = offset; i < pressure_boundary_ids_.size(); ++i)
682 {
683 pressure_boundary_ids_[i] = j_boundary[i - offset]["id"];
684
685 auto ff = j_boundary[i - offset]["value"];
686 pressures_[i].value.init(ff);
687 if (j_boundary[i - offset].contains("time_reference") && j_boundary[i - offset]["time_reference"].size() > 0)
688 pressures_[i].value.set_t(j_boundary[i - offset]["time_reference"]);
689
690 pressures_[i].interpolation = std::make_shared<NoInterpolation>();
691 }
692 }
693
694 if (is_param_valid(params, "pressure_cavity"))
695 {
696 const int offset = pressure_cavity_ids_.size();
697
698 auto j_boundary_tmp = params["pressure_cavity"];
699 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
700
701 pressure_cavity_ids_.resize(offset + j_boundary.size());
702
703 for (size_t i = offset; i < pressure_cavity_ids_.size(); ++i)
704 {
705 int boundary_id = j_boundary[i - offset]["id"];
706 pressure_cavity_ids_[i] = boundary_id;
707
708 if (cavity_pressures_.find(boundary_id) == cavity_pressures_.end())
709 {
710 cavity_pressures_[boundary_id] = ScalarBCValue();
711
712 auto ff = j_boundary[i - offset]["value"];
713 cavity_pressures_[boundary_id].value.init(ff);
714
715 cavity_pressures_[boundary_id].interpolation = std::make_shared<NoInterpolation>();
716 }
717 }
718 }
719
720 if (is_param_valid(params, "solution"))
721 {
722 auto rr = params["solution"];
723 initial_position_.resize(rr.size());
724 assert(rr.is_array());
725
726 for (size_t k = 0; k < rr.size(); ++k)
727 {
728 initial_position_[k].first = rr[k]["id"];
729 const auto v = rr[k]["value"];
730 for (size_t d = 0; d < v.size(); ++d)
731 initial_position_[k].second[d].init(v[d]);
732 }
733 }
734
735 if (is_param_valid(params, "velocity"))
736 {
737 auto rr = params["velocity"];
738 initial_velocity_.resize(rr.size());
739 assert(rr.is_array());
740
741 for (size_t k = 0; k < rr.size(); ++k)
742 {
743 initial_velocity_[k].first = rr[k]["id"];
744 const auto v = rr[k]["value"];
745 for (size_t d = 0; d < v.size(); ++d)
746 initial_velocity_[k].second[d].init(v[d]);
747 }
748 }
749
750 if (is_param_valid(params, "acceleration"))
751 {
752 auto rr = params["acceleration"];
753 initial_acceleration_.resize(rr.size());
754 assert(rr.is_array());
755
756 for (size_t k = 0; k < rr.size(); ++k)
757 {
758 initial_acceleration_[k].first = rr[k]["id"];
759 const auto v = rr[k]["value"];
760 for (size_t d = 0; d < v.size(); ++d)
761 initial_acceleration_[k].second[d].init(v[d]);
762 }
763 }
764 }
765
766 void GenericTensorProblem::initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const
767 {
768 val.resize(pts.rows(), pts.cols());
769 if (initial_position_.empty())
770 {
771 val.setZero();
772 return;
773 }
774
775 const bool planar = pts.cols() == 2;
776 for (int i = 0; i < pts.rows(); ++i)
777 {
778 const int id = mesh.get_body_id(global_ids(i));
779 int index = -1;
780 for (int j = 0; j < initial_position_.size(); ++j)
781 {
782 if (initial_position_[j].first == id)
783 {
784 index = j;
785 break;
786 }
787 }
788 if (index < 0)
789 {
790 val.row(i).setZero();
791 continue;
792 }
793
794 for (int j = 0; j < pts.cols(); ++j)
795 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));
796 }
797 }
798
799 void GenericTensorProblem::initial_velocity(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const
800 {
801 val.resize(pts.rows(), pts.cols());
802 if (initial_velocity_.empty())
803 {
804 val.setZero();
805 return;
806 }
807
808 const bool planar = pts.cols() == 2;
809 for (int i = 0; i < pts.rows(); ++i)
810 {
811 const int id = mesh.get_body_id(global_ids(i));
812 int index = -1;
813 for (int j = 0; j < initial_velocity_.size(); ++j)
814 {
815 if (initial_velocity_[j].first == id)
816 {
817 index = j;
818 break;
819 }
820 }
821 if (index < 0)
822 {
823 val.row(i).setZero();
824 continue;
825 }
826
827 for (int j = 0; j < pts.cols(); ++j)
828 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));
829 }
830 }
831
832 void GenericTensorProblem::initial_acceleration(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const
833 {
834 val.resize(pts.rows(), pts.cols());
835 if (initial_acceleration_.empty())
836 {
837 val.setZero();
838 return;
839 }
840
841 const bool planar = pts.cols() == 2;
842 for (int i = 0; i < pts.rows(); ++i)
843 {
844 const int id = mesh.get_body_id(global_ids(i));
845 int index = -1;
846 for (int j = 0; j < initial_acceleration_.size(); ++j)
847 {
848 if (initial_acceleration_[j].first == id)
849 {
850 index = j;
851 break;
852 }
853 }
854 if (index < 0)
855 {
856 val.row(i).setZero();
857 continue;
858 }
859
860 for (int j = 0; j < pts.cols(); ++j)
861 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));
862 }
863 }
864
866 {
868 has_exact_ = false;
869 has_exact_grad_ = false;
870 is_time_dept_ = false;
871
872 forces_.clear();
873 displacements_.clear();
875 pressures_.clear();
876 cavity_pressures_.clear();
877
878 nodal_dirichlet_.clear();
879 nodal_neumann_.clear();
880 nodal_dirichlet_mat_.clear();
881
882 initial_position_.clear();
883 initial_velocity_.clear();
884 initial_acceleration_.clear();
885
886 for (int i = 0; i < rhs_.size(); ++i)
887 rhs_[i].clear();
888 for (int i = 0; i < exact_.size(); ++i)
889 exact_[i].clear();
890 for (int i = 0; i < exact_grad_.size(); ++i)
891 exact_grad_[i].clear();
892 is_all_ = false;
894 }
895
897 : Problem(name), is_all_(false)
898 {
899 }
900
901 void GenericScalarProblem::set_units(const assembler::Assembler &assembler, const Units &units)
902 {
903 // TODO?
904
905 for (auto &v : neumann_)
906 v.set_unit_type("");
907
908 for (auto &v : dirichlet_)
909 v.set_unit_type("");
910
911 for (auto &v : initial_solution_)
912 v.second.set_unit_type("");
913
916
917 for (int i = 0; i < 3; ++i)
918 exact_grad_[i].set_unit_type("");
919 }
920
921 void GenericScalarProblem::rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
922 {
923 val.resize(pts.rows(), 1);
924 if (is_rhs_zero())
925 {
926 val.setZero();
927 return;
928 }
929 const bool planar = pts.cols() == 2;
930 for (int i = 0; i < pts.rows(); ++i)
931 {
932 double x = pts(i, 0), y = pts(i, 1), z = planar ? 0 : pts(i, 2);
933 val(i) = rhs_(x, y, z, t);
934 }
935 }
936
937 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
938 {
939 val = Eigen::MatrixXd::Zero(pts.rows(), 1);
940
941 for (long i = 0; i < pts.rows(); ++i)
942 {
943 const int id = mesh.get_boundary_id(global_ids(i));
944 if (is_all_)
945 {
946 assert(dirichlet_.size() == 1);
947 val(i) = dirichlet_[0].eval(pts.row(i), t);
948 }
949 else
950 {
951 for (size_t b = 0; b < boundary_ids_.size(); ++b)
952 {
953 if (id == boundary_ids_[b])
954 {
955 val(i) = dirichlet_[b].eval(pts.row(i), t);
956 break;
957 }
958 }
959 }
960 }
961 }
962
963 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
964 {
965 val = Eigen::MatrixXd::Zero(pts.rows(), 1);
966
967 for (long i = 0; i < pts.rows(); ++i)
968 {
969 const int id = mesh.get_boundary_id(global_ids(i));
970
971 for (size_t b = 0; b < neumann_boundary_ids_.size(); ++b)
972 {
973 if (id == neumann_boundary_ids_[b])
974 {
975 double x = pts(i, 0), y = pts(i, 1), z = pts.cols() == 2 ? 0 : pts(i, 2);
976 val(i) = neumann_[b].eval(pts.row(i), t);
977 break;
978 }
979 }
980 }
981 }
982
983 void GenericScalarProblem::initial_solution(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &pts, Eigen::MatrixXd &val) const
984 {
985 val.resize(pts.rows(), 1);
986 if (initial_solution_.empty())
987 {
988 val.setZero();
989 return;
990 }
991
992 const bool planar = pts.cols() == 2;
993 for (int i = 0; i < pts.rows(); ++i)
994 {
995 const int id = mesh.get_body_id(global_ids(i));
996 int index = -1;
997 for (int j = 0; j < initial_solution_.size(); ++j)
998 {
999 if (initial_solution_[j].first == id)
1000 {
1001 index = j;
1002 break;
1003 }
1004 }
1005 if (index < 0)
1006 {
1007 val(i) = 0;
1008 continue;
1009 }
1010
1011 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));
1012 }
1013 }
1014
1015 void GenericScalarProblem::exact(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
1016 {
1017 assert(has_exact_sol());
1018 const bool planar = pts.cols() == 2;
1019 val.resize(pts.rows(), 1);
1020
1021 for (int i = 0; i < pts.rows(); ++i)
1022 {
1023 double x = pts(i, 0), y = pts(i, 1), z = pts.cols() == 2 ? 0 : pts(i, 2);
1024 val(i) = exact_(x, y, z, t);
1025 }
1026 }
1027
1028 void GenericScalarProblem::exact_grad(const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const
1029 {
1030 val.resize(pts.rows(), pts.cols());
1031 if (!has_exact_grad_)
1032 return;
1033
1034 const bool planar = pts.cols() == 2;
1035 for (int i = 0; i < pts.rows(); ++i)
1036 {
1037 for (int j = 0; j < pts.cols(); ++j)
1038 {
1039 double x = pts(i, 0), y = pts(i, 1), z = pts.cols() == 2 ? 0 : pts(i, 2);
1040 val(i, j) = exact_grad_[j](x, y, z, t);
1041 }
1042 }
1043 }
1044
1045 void GenericTensorProblem::update_pressure_boundary(const int id, const int time_step, const double val)
1046 {
1047 int index = -1;
1048 for (int i = 0; i < pressure_boundary_ids_.size(); ++i)
1049 {
1050 if (pressure_boundary_ids_[i] == id)
1051 {
1052 index = i;
1053 break;
1054 }
1055 }
1056 if (index == -1)
1057 {
1058 throw "Invalid boundary id for pressure update";
1059 }
1060
1061 if (pressures_[index].value.is_mat())
1062 {
1063 Eigen::MatrixXd curr_val = pressures_[index].value.get_mat();
1064 assert(time_step <= curr_val.size());
1065 assert(curr_val.cols() == 1);
1066 curr_val(time_step) = val;
1067 pressures_[index].value.set_mat(curr_val);
1068 }
1069 else
1070 {
1071 pressures_[index].value.init(val);
1072 }
1073 }
1074
1075 void GenericTensorProblem::update_dirichlet_boundary(const int id, const int time_step, const Eigen::VectorXd &val)
1076 {
1077 int index = -1;
1078 for (int i = 0; i < boundary_ids_.size(); ++i)
1079 {
1080 if (boundary_ids_[i] == id)
1081 {
1082 index = i;
1083 break;
1084 }
1085 }
1086 if (index == -1)
1087 {
1088 throw "Invalid boundary id for dirichlet update";
1089 }
1090
1091 for (int i = 0; i < val.size(); ++i)
1092 {
1093 if (displacements_[index].value[i].is_mat())
1094 {
1095 Eigen::MatrixXd curr_val = displacements_[index].value[i].get_mat();
1096 assert(time_step <= curr_val.size());
1097 assert(curr_val.cols() == 1);
1098 curr_val(time_step) = val(i);
1099 displacements_[index].value[i].set_mat(curr_val);
1100 }
1101 else
1102 {
1103 displacements_[index].value[i].init(val(i));
1104 }
1105 }
1106 }
1107
1108 void GenericTensorProblem::update_dirichlet_nodes(const Eigen::VectorXi &in_node_to_node, const Eigen::VectorXi &node_ids, const Eigen::MatrixXd &nodal_dirichlet)
1109 {
1110 assert(node_ids.size() == nodal_dirichlet.rows());
1111 // NOTE!!! Update nodes called in `State::build_basis()` so the first row of this mat
1112 // must be set to input node ordering if `build_basis` will be called, like in optimization
1113 for (auto &n_dirichlet : nodal_dirichlet_mat_)
1114 for (int i = 0; i < node_ids.size(); ++i)
1115 {
1116 int mapped_node_id = in_node_to_node(node_ids(i));
1117 for (int j = 0; j < n_dirichlet.rows(); ++j)
1118 if (mapped_node_id == n_dirichlet(j, 0))
1119 for (int k = 0; k < n_dirichlet.cols() - 1; ++k)
1120 n_dirichlet(j, k + 1) = nodal_dirichlet(i, k);
1121 }
1122 }
1123
1124 void GenericScalarProblem::dirichlet_nodal_value(const mesh::Mesh &mesh, const int node_id, const RowVectorNd &pt, const double t, Eigen::MatrixXd &val) const
1125 {
1126 val = Eigen::MatrixXd::Zero(1, 1);
1127 const int tag = mesh.get_node_id(node_id);
1128
1129 if (is_all_)
1130 {
1131 assert(nodal_dirichlet_.size() == 1);
1132 const auto &tmp = nodal_dirichlet_.begin()->second;
1133
1134 val(0) = tmp.eval(pt, t);
1135 return;
1136 }
1137
1138 const auto it = nodal_dirichlet_.find(tag);
1139 if (it != nodal_dirichlet_.end())
1140 {
1141 val(0) = it->second.eval(pt, t);
1142 return;
1143 }
1144
1145 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
1146 {
1147 for (int i = 0; i < n_dirichlet.rows(); ++i)
1148 {
1149 if (n_dirichlet(i, 0) == node_id)
1150 {
1151 val(0) = n_dirichlet(i, 1);
1152 return;
1153 }
1154 }
1155 }
1156
1157 assert(false);
1158 }
1159
1160 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
1161 {
1162 // TODO implement me;
1163 log_and_throw_error("Nodal neumann not implemented");
1164 }
1165
1166 bool GenericScalarProblem::is_nodal_dirichlet_boundary(const int n_id, const int tag)
1167 {
1168 if (nodal_dirichlet_.find(tag) != nodal_dirichlet_.end())
1169 return true;
1170
1171 for (const auto &n_dirichlet : nodal_dirichlet_mat_)
1172 {
1173 for (int i = 0; i < n_dirichlet.rows(); ++i)
1174 {
1175 if (n_dirichlet(i, 0) == n_id)
1176 return true;
1177 }
1178 }
1179
1180 return false;
1181 }
1182
1183 bool GenericScalarProblem::is_nodal_neumann_boundary(const int n_id, const int tag)
1184 {
1185 return nodal_neumann_.find(tag) != nodal_neumann_.end();
1186 }
1187
1189 {
1190 return nodal_dirichlet_mat_.size() > 0;
1191 }
1192
1194 {
1195 return false; // nodal_neumann_.size() > 0;
1196 }
1197
1198 void GenericScalarProblem::update_nodes(const Eigen::VectorXi &in_node_to_node)
1199 {
1201 return;
1202 for (auto &n_dirichlet : nodal_dirichlet_mat_)
1203 {
1204 for (int n = 0; n < n_dirichlet.rows(); ++n)
1205 {
1206 const int node_id = in_node_to_node[n_dirichlet(n, 0)];
1207 n_dirichlet(n, 0) = node_id;
1208 }
1209 }
1211 }
1212
1214 {
1215 if (is_param_valid(params, "is_time_dependent"))
1216 {
1217 is_time_dept_ = params["is_time_dependent"];
1218 }
1219
1220 if (is_param_valid(params, "rhs"))
1221 {
1222 rhs_.init(params["rhs"]);
1223 }
1224
1225 if (is_param_valid(params, "reference") && is_param_valid(params["reference"], "solution"))
1226 {
1227 has_exact_ = !params["reference"]["solution"].empty();
1228 exact_.init(params["reference"]["solution"]);
1229 }
1230
1231 if (is_param_valid(params, "reference") && is_param_valid(params["reference"], "gradient"))
1232 {
1233 auto ex = params["reference"]["gradient"];
1234 has_exact_grad_ = ex.size() > 0;
1235 if (ex.is_array())
1236 {
1237 for (size_t k = 0; k < ex.size(); ++k)
1238 exact_grad_[k].init(ex[k]);
1239 }
1240 else
1241 {
1242 assert(false);
1243 }
1244 }
1245
1246 if (is_param_valid(params, "dirichlet_boundary"))
1247 {
1248 // boundary_ids_.clear();
1249 const int offset = boundary_ids_.size();
1250 std::vector<json> j_boundary = flatten_ids(params["dirichlet_boundary"]);
1251
1252 boundary_ids_.resize(offset + j_boundary.size());
1253 dirichlet_.resize(offset + j_boundary.size());
1254
1255 for (size_t i = offset; i < boundary_ids_.size(); ++i)
1256 {
1257 if (j_boundary[i - offset].is_string())
1258 {
1259 const std::string path = resolve_path(j_boundary[i - offset], params["root_path"]);
1260 if (!std::filesystem::is_regular_file(path))
1261 log_and_throw_error(fmt::format("unable to open {} file", path));
1262
1263 Eigen::MatrixXd tmp;
1264 io::read_matrix(path, tmp);
1265 nodal_dirichlet_mat_.emplace_back(tmp);
1266
1267 continue;
1268 }
1269
1270 int current_id = -1;
1271
1272 if (j_boundary[i - offset]["id"] == "all")
1273 {
1274 assert(boundary_ids_.size() == 1);
1275
1276 is_all_ = true;
1277 boundary_ids_.clear();
1278 nodal_dirichlet_[current_id] = ScalarBCValue();
1279 }
1280 else
1281 {
1282 boundary_ids_[i] = j_boundary[i - offset]["id"];
1283 current_id = boundary_ids_[i];
1284 nodal_dirichlet_[current_id] = ScalarBCValue();
1285 }
1286
1287 auto ff = j_boundary[i - offset]["value"];
1288 dirichlet_[i].value.init(ff);
1289 nodal_dirichlet_[current_id].value.init(ff);
1290
1291 if (j_boundary[i - offset]["interpolation"].is_array())
1292 {
1293 if (j_boundary[i - offset]["interpolation"].size() == 0)
1294 dirichlet_[i].interpolation = std::make_shared<NoInterpolation>();
1295 else if (j_boundary[i - offset]["interpolation"].size() == 1)
1296 dirichlet_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"][0]);
1297 else
1298 log_and_throw_error("Only one Dirichlet interpolation supported");
1299 }
1300 else
1301 dirichlet_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"]);
1302
1303 nodal_dirichlet_[current_id].interpolation = dirichlet_[i].interpolation;
1304 }
1305 }
1306
1307 if (is_param_valid(params, "neumann_boundary"))
1308 {
1309 // neumann_boundary_ids_.clear();
1310 const int offset = neumann_boundary_ids_.size();
1311 auto j_boundary_tmp = params["neumann_boundary"];
1312 std::vector<json> j_boundary = flatten_ids(j_boundary_tmp);
1313
1314 neumann_boundary_ids_.resize(offset + j_boundary.size());
1315 neumann_.resize(offset + j_boundary.size());
1316
1317 for (size_t i = offset; i < neumann_boundary_ids_.size(); ++i)
1318 {
1319 neumann_boundary_ids_[i] = j_boundary[i - offset]["id"];
1320
1321 auto ff = j_boundary[i - offset]["value"];
1322 neumann_[i].value.init(ff);
1323
1324 if (j_boundary[i - offset]["interpolation"].is_array())
1325 {
1326 if (j_boundary[i - offset]["interpolation"].size() == 0)
1327 neumann_[i].interpolation = std::make_shared<NoInterpolation>();
1328 else if (j_boundary[i - offset]["interpolation"].size() == 1)
1329 neumann_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"][0]);
1330 else
1331 log_and_throw_error("Only one Neumann interpolation supported");
1332 }
1333 else
1334 neumann_[i].interpolation = Interpolation::build(j_boundary[i - offset]["interpolation"]);
1335 }
1336 }
1337
1338 if (is_param_valid(params, "solution"))
1339 {
1340 auto rr = params["solution"];
1341 initial_solution_.resize(rr.size());
1342 assert(rr.is_array());
1343
1344 for (size_t k = 0; k < rr.size(); ++k)
1345 {
1346 initial_solution_[k].first = rr[k]["id"];
1347 initial_solution_[k].second.init(rr[k]["value"]);
1348 }
1349 }
1350 }
1351
1353 {
1354 neumann_.clear();
1355 dirichlet_.clear();
1356
1357 nodal_dirichlet_.clear();
1358 nodal_neumann_.clear();
1359 nodal_dirichlet_mat_.clear();
1360
1361 rhs_.clear();
1362 exact_.clear();
1363 for (int i = 0; i < exact_grad_.size(); ++i)
1364 exact_grad_[i].clear();
1365 is_all_ = false;
1366 has_exact_ = false;
1367 has_exact_grad_ = false;
1368 is_time_dept_ = false;
1370 }
1371 } // namespace assembler
1372} // namespace polyfem
double val
Definition Assembler.cpp:86
int y
int z
int x
std::string force() const
Definition Units.hpp:31
std::string pressure() const
Definition Units.hpp:32
std::string velocity() const
Definition Units.hpp:29
std::string acceleration() const
Definition Units.hpp:30
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 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
GenericScalarProblem(const std::string &name)
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
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_
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_
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:501
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
virtual int get_node_id(const int node_id) const
Get the boundary selection of a node.
Definition Mesh.hpp:484
void set_unit_type(const std::string &unit_type)
static std::shared_ptr< Interpolation > build(const json &params)
list tmp
Definition p_bases.py:365
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:13
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