PolyFEM
Loading...
Searching...
No Matches
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 if (T.has_children())
13 {
14 for (int i = 0; i < T.n_children(); i++)
15 this->children.push_back(std::make_unique<Tree>(T.child(i)));
16 }
17 }
18 Tree operator=(const Tree &T) {
19 if (this != &T)
20 {
21 this->children.clear();
22 if (T.has_children())
23 {
24 for (int i = 0; i < T.n_children(); i++)
25 this->children.push_back(std::make_unique<Tree>(T.child(i)));
26 }
27 }
28 return *this;
29 }
30 bool merge(const Tree &T, int max_depth = 2) {
31 bool flag = false;
32 if (!T.has_children() || max_depth <= 0)
33 return flag;
34 if (!this->has_children())
35 {
36 this->add_children(T.n_children());
37 flag = true;
38 max_depth--;
39 }
40 for (int i = 0; i < T.n_children(); i++)
41 flag = this->child(i).merge(T.child(i), max_depth) || flag;
42 return flag;
43 }
44
45 // Debug print
46 friend std::ostream& operator<<(
47 std::ostream& ost, const Tree& T
48 ) {
49 ost << "(";
50 if (T.has_children())
51 {
52 for (int i = 0; i < T.n_children(); i++)
53 ost << T.child(i) << ", ";
54 }
55 ost << ")";
56 return ost;
57 }
58
59 bool has_children() const { return !children.empty(); }
60 int n_children() const { return children.size(); }
61 int depth() const {
62 if (!has_children())
63 return 0;
64 int d = 0;
65 for (int i = 0; i < n_children(); i++)
66 d = std::max(d, child(i).depth());
67 return d + 1;
68 }
69 int n_leaves() const {
70 if (!has_children())
71 return 1;
72 int n = 0;
73 for (int i = 0; i < n_children(); i++)
74 n += child(i).n_leaves();
75 return n;
76 }
77 Tree& child(int i) { return *children[i]; }
78 const Tree& child(int i) const { return *children[i]; }
79 void add_children(int n) {
80 for (int i = 0; i < n; i++)
81 children.push_back(std::make_unique<Tree>());
82 }
83
84 private:
85 std::vector<std::unique_ptr<Tree>> children;
86 };
87
88 Eigen::VectorXd robust_evaluate_jacobian(
89 const int order,
90 const Eigen::MatrixXd &cp,
91 const Eigen::MatrixXd &uv);
92
93 std::vector<int> count_invalid(
94 const int dim,
95 const std::vector<basis::ElementBases> &bases,
96 const std::vector<basis::ElementBases> &gbases,
97 const Eigen::VectorXd &u);
98
99 std::tuple<bool, int, Tree> is_valid(
100 const int dim,
101 const std::vector<basis::ElementBases> &bases,
102 const std::vector<basis::ElementBases> &gbases,
103 const Eigen::VectorXd &u,
104 const double threshold = 0);
105
106 bool is_valid(
107 const int dim,
108 const std::vector<basis::ElementBases> &bases,
109 const std::vector<basis::ElementBases> &gbases,
110 const Eigen::VectorXd &u1,
111 const Eigen::VectorXd &u2,
112 const double threshold = 0);
113
114 std::tuple<double, int, double, Tree> max_time_step(
115 const int dim,
116 const std::vector<basis::ElementBases> &bases,
117 const std::vector<basis::ElementBases> &gbases,
118 const Eigen::VectorXd &u1,
119 const Eigen::VectorXd &u2,
120 double precision = .25);
121
122 Eigen::MatrixXd extract_nodes(const int dim, const basis::ElementBases &basis, const basis::ElementBases &gbasis, const Eigen::VectorXd &u, int order);
123 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);
124}
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:85
Tree(const Tree &T)
Definition Jacobian.hpp:11
bool merge(const Tree &T, int max_depth=2)
Definition Jacobian.hpp:30
const Tree & child(int i) const
Definition Jacobian.hpp:78
int depth() const
Definition Jacobian.hpp:61
int n_children() const
Definition Jacobian.hpp:60
bool has_children() const
Definition Jacobian.hpp:59
void add_children(int n)
Definition Jacobian.hpp:79
int n_leaves() const
Definition Jacobian.hpp:69
Tree operator=(const Tree &T)
Definition Jacobian.hpp:18
Tree & child(int i)
Definition Jacobian.hpp:77
friend std::ostream & operator<<(std::ostream &ost, const Tree &T)
Definition Jacobian.hpp:46
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:241
Eigen::VectorXd robust_evaluate_jacobian(const int order, const Eigen::MatrixXd &cp, const Eigen::MatrixXd &uv)
Definition Jacobian.cpp:124
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:345
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:173
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:77