PolyFEM
Loading...
Searching...
No Matches
PointBasedProblem.hpp
Go to the documentation of this file.
1#pragma once
2
6
7#include <polyfem/Common.hpp>
8
9#include <vector>
10#include <Eigen/Dense>
11
12namespace polyfem
13{
14 namespace problem
15 {
17 {
18 private:
19 class BCValue
20 {
21 public:
23 {
24 Eigen::Matrix<bool, 3, 1> dd;
25 dd.setConstant(true);
26 init(0, 0, 0, dd);
27 }
28
29 Eigen::RowVector3d operator()(const Eigen::RowVector3d &pt) const;
30
31 bool init(const json &data);
32
33 void init(const double x, const double y, const double z, const Eigen::Matrix<bool, 3, 1> &dd)
34 {
35 this->val << x, y, z;
36 dirichlet_dims = dd;
37 is_val = true;
38 }
39
40 void init(const Eigen::Vector3d &v, const Eigen::Matrix<bool, 3, 1> &dd)
41 {
42 this->val = v;
43 dirichlet_dims = dd;
44 is_val = true;
45 }
46
47 void init(const Eigen::MatrixXd &pts, const Eigen::MatrixXi &tri, const Eigen::MatrixXd &fun, const int coord, const Eigen::Matrix<bool, 3, 1> &dd)
48 {
49 tri_func.init(fun, pts, tri);
50 is_tri = true;
51
52 coordiante_0 = (coord + 1) % 3;
53 coordiante_1 = (coord + 2) % 3;
54
55 dirichlet_dims = dd;
56
57 is_val = false;
58 }
59
60 void init(const Eigen::MatrixXd &pts, const Eigen::MatrixXd &fun, const std::string &rbf, const double eps, const int coord, const Eigen::Matrix<bool, 3, 1> &dd)
61 {
62 rbf_func.init(fun, pts, rbf, eps);
63 is_tri = false;
64
65 if (coord >= 0)
66 {
67 coordiante_0 = (coord + 1) % 3;
68 coordiante_1 = (coord + 2) % 3;
69 }
70 else
71 {
72 coordiante_0 = -1;
73 coordiante_1 = -1;
74 }
75
76 dirichlet_dims = dd;
77
78 is_val = false;
79 }
80
81 bool is_dirichet_dim(const int d) const { return dirichlet_dims(d); }
82
83 private:
84 Eigen::Vector3d val;
87 bool is_val;
88 bool is_tri;
89 int coordiante_0 = 0;
90 int coordiante_1 = 1;
91 Eigen::Matrix<bool, 3, 1> dirichlet_dims;
92 };
93
94 public:
95 PointBasedTensorProblem(const std::string &name);
96
97 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
98 bool is_rhs_zero() const override { return abs(rhs_) < 1e-10; }
99
100 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;
101 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;
102
103 bool has_exact_sol() const override { return false; }
104 bool is_scalar() const override { return false; }
105
106 void set_parameters(const json &params) override;
107
108 void add_constant(const int bc_tag, const Eigen::Vector3d &value)
109 {
110 Eigen::Matrix<bool, 3, 1> dd;
111 dd.setConstant(true);
112 add_constant(bc_tag, value, dd);
113 }
114
115 void add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const Eigen::MatrixXi &tri, const int coord)
116 {
117 Eigen::Matrix<bool, 3, 1> dd;
118 dd.setConstant(true);
119 add_function(bc_tag, func, pts, tri, coord, dd);
120 }
121
122 void add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const std::string &rbf, const double eps, const int coord)
123 {
124 Eigen::Matrix<bool, 3, 1> dd;
125 dd.setConstant(true);
126 add_function(bc_tag, func, pts, rbf, eps, coord, dd);
127 }
128
129 void add_constant(const int bc_tag, const Eigen::Vector3d &value, const Eigen::Matrix<bool, 3, 1> &dd, const bool is_neumann = false);
130 void add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const Eigen::MatrixXi &tri, const int coord, const Eigen::Matrix<bool, 3, 1> &dd, const bool is_neumann = false);
131 void add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const std::string &rbf, const double eps, const int coord, const Eigen::Matrix<bool, 3, 1> &dd, const bool is_neumann = false);
132
133 bool is_dimension_dirichet(const int tag, const int dim) const override;
134 bool all_dimensions_dirichlet() const override { return all_dimensions_dirichlet_; }
135
136 private:
137 bool initialized_ = false;
139 double rhs_;
140 double scaling_;
141 Eigen::Vector3d translation_;
142 std::vector<BCValue> bc_;
143 std::vector<BCValue> neumann_bc_;
144 };
145 } // namespace problem
146} // namespace polyfem
double val
Definition Assembler.cpp:86
int y
int z
int x
const std::string & name() const
Definition Problem.hpp:23
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
void init(const double x, const double y, const double z, const Eigen::Matrix< bool, 3, 1 > &dd)
void init(const Eigen::Vector3d &v, const Eigen::Matrix< bool, 3, 1 > &dd)
void init(const Eigen::MatrixXd &pts, const Eigen::MatrixXd &fun, const std::string &rbf, const double eps, const int coord, const Eigen::Matrix< bool, 3, 1 > &dd)
Eigen::RowVector3d operator()(const Eigen::RowVector3d &pt) const
void init(const Eigen::MatrixXd &pts, const Eigen::MatrixXi &tri, const Eigen::MatrixXd &fun, const int coord, const Eigen::Matrix< bool, 3, 1 > &dd)
void add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const std::string &rbf, const double eps, const int coord)
void set_parameters(const json &params) override
void add_constant(const int bc_tag, const Eigen::Vector3d &value)
void dirichlet_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void add_function(const int bc_tag, const Eigen::MatrixXd &func, const Eigen::MatrixXd &pts, const Eigen::MatrixXi &tri, const int coord)
bool is_dimension_dirichet(const int tag, const int dim) const override
void neumann_bc(const mesh::Mesh &mesh, const Eigen::MatrixXi &global_ids, const Eigen::MatrixXd &uv, const Eigen::MatrixXd &pts, const Eigen::MatrixXd &normals, const double t, Eigen::MatrixXd &val) const override
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override
void init(const Eigen::MatrixXd &fun, const Eigen::MatrixXd &pts, const Eigen::MatrixXi &tris)
void init(const Eigen::MatrixXd &fun, const Eigen::MatrixXd &pts, const std::function< double(double)> &rbf)
nlohmann::json json
Definition Common.hpp:9