Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Jacobian.hpp
Go to the documentation of this file.
1#pragma once
2
4
5namespace polyfem::utils
6{
7 class Tree
8 {
9 public:
10 Tree() {}
11 Tree(const Tree &T)
12 {
13 if (T.has_children())
14 {
15 for (int i = 0; i < T.n_children(); i++)
16 this->children.push_back(std::make_unique<Tree>(T.child(i)));
17 }
18 }
20 {
21 if (this != &T)
22 {
23 this->children.clear();
24 if (T.has_children())
25 {
26 for (int i = 0; i < T.n_children(); i++)
27 this->children.push_back(std::make_unique<Tree>(T.child(i)));
28 }
29 }
30 return *this;
31 }
32 bool merge(const Tree &T, int max_depth = 2)
33 {
34 bool flag = false;
35 if (!T.has_children() || max_depth <= 0)
36 return flag;
37 if (!this->has_children())
38 {
39 this->add_children(T.n_children());
40 flag = true;
41 max_depth--;
42 }
43 for (int i = 0; i < T.n_children(); i++)
44 flag = this->child(i).merge(T.child(i), max_depth) || flag;
45 return flag;
46 }
47
48 // Debug print
49 friend std::ostream &operator<<(
50 std::ostream &ost, const Tree &T)
51 {
52 ost << "(";
53 if (T.has_children())
54 {
55 for (int i = 0; i < T.n_children(); i++)
56 ost << T.child(i) << ", ";
57 }
58 ost << ")";
59 return ost;
60 }
61
62 bool has_children() const { return !children.empty(); }
63 int n_children() const { return children.size(); }
64 int depth() const
65 {
66 if (!has_children())
67 return 0;
68 int d = 0;
69 for (int i = 0; i < n_children(); i++)
70 d = std::max(d, child(i).depth());
71 return d + 1;
72 }
73 int n_leaves() const
74 {
75 if (!has_children())
76 return 1;
77 int n = 0;
78 for (int i = 0; i < n_children(); i++)
79 n += child(i).n_leaves();
80 return n;
81 }
82 Tree &child(int i) { return *children[i]; }
83 const Tree &child(int i) const { return *children[i]; }
84 void add_children(int n)
85 {
86 for (int i = 0; i < n; i++)
87 children.push_back(std::make_unique<Tree>());
88 }
89
90 private:
91 std::vector<std::unique_ptr<Tree>> children;
92 };
93
94 Eigen::VectorXd robust_evaluate_jacobian(
95 const int order,
96 const Eigen::MatrixXd &cp,
97 const Eigen::MatrixXd &uv);
98
99 std::vector<int> count_invalid(
100 const int dim,
101 const std::vector<basis::ElementBases> &bases,
102 const std::vector<basis::ElementBases> &gbases,
103 const Eigen::VectorXd &u);
104
105 std::tuple<bool, int, Tree> is_valid(
106 const int dim,
107 const std::vector<basis::ElementBases> &bases,
108 const std::vector<basis::ElementBases> &gbases,
109 const Eigen::VectorXd &u,
110 const double threshold = 0);
111
112 bool is_valid(
113 const int dim,
114 const std::vector<basis::ElementBases> &bases,
115 const std::vector<basis::ElementBases> &gbases,
116 const Eigen::VectorXd &u1,
117 const Eigen::VectorXd &u2,
118 const double threshold = 0);
119
120 std::tuple<double, int, double, Tree> max_time_step(
121 const int dim,
122 const std::vector<basis::ElementBases> &bases,
123 const std::vector<basis::ElementBases> &gbases,
124 const Eigen::VectorXd &u1,
125 const Eigen::VectorXd &u2,
126 double precision = .25);
127
128 Eigen::MatrixXd extract_nodes(const int dim, const basis::ElementBases &basis, const basis::ElementBases &gbasis, const Eigen::VectorXd &u, int order);
129 Eigen::MatrixXd extract_nodes(const int dim, const std::vector<basis::ElementBases> &bases, const std::vector<basis::ElementBases> &gbases, const Eigen::VectorXd &u, int order, int n_elem = -1);
130} // namespace polyfem::utils
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
std::vector< std::unique_ptr< Tree > > children
Definition Jacobian.hpp:91
Tree(const Tree &T)
Definition Jacobian.hpp:11
bool merge(const Tree &T, int max_depth=2)
Definition Jacobian.hpp:32
const Tree & child(int i) const
Definition Jacobian.hpp:83
int depth() const
Definition Jacobian.hpp:64
int n_children() const
Definition Jacobian.hpp:63
bool has_children() const
Definition Jacobian.hpp:62
void add_children(int n)
Definition Jacobian.hpp:84
int n_leaves() const
Definition Jacobian.hpp:73
Tree operator=(const Tree &T)
Definition Jacobian.hpp:19
Tree & child(int i)
Definition Jacobian.hpp:82
friend std::ostream & operator<<(std::ostream &ost, const Tree &T)
Definition Jacobian.hpp:49
std::tuple< bool, int, Tree > is_valid(const int dim, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXd &u, const double threshold)
Definition Jacobian.cpp:297
Eigen::VectorXd robust_evaluate_jacobian(const int order, const Eigen::MatrixXd &cp, const Eigen::MatrixXd &uv)
Definition Jacobian.cpp:145
std::tuple< double, int, double, Tree > max_time_step(const int dim, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXd &u1, const Eigen::VectorXd &u2, double precision)
Definition Jacobian.cpp:409
std::vector< int > count_invalid(const int dim, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXd &u)
Definition Jacobian.cpp:212
Eigen::MatrixXd extract_nodes(const int dim, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXd &u, int order, int n_elem)
Definition Jacobian.cpp:98