PolyFEM
Loading...
Searching...
No Matches
NodeProblem.hpp
Go to the documentation of this file.
1#pragma once
2
6
7#include <vector>
8#include <Eigen/Dense>
9
10namespace polyfem
11{
12 namespace problem
13 {
15 {
16 public:
17 NodeValues();
18
19 void load(const std::string &path);
20 void init(const mesh::Mesh &mesh);
21
22 double dirichlet_interpolate(const int p_id, const Eigen::MatrixXd &uv) const
23 {
24 return interpolate(p_id, uv, true);
25 }
26 double neumann_interpolate(const int p_id, const Eigen::MatrixXd &uv) const
27 {
28 return interpolate(p_id, uv, false);
29 }
30
31 private:
32 double interpolate(const int p_id, const Eigen::MatrixXd &uv, bool is_dirichlet) const;
33
34 std::vector<int> raw_ids_;
35 std::vector<std::vector<double>> raw_data_;
36 std::vector<bool> raw_dirichlet_;
37
38 std::vector<Eigen::VectorXd> data_;
39 std::vector<bool> dirichlet_;
40 };
41
43 {
44 public:
45 NodeProblem(const std::string &name);
46 void init(const mesh::Mesh &mesh) override;
47
48 void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, const double t, Eigen::MatrixXd &val) const override;
49 bool is_rhs_zero() const override { return abs(rhs_) < 1e-10; }
50
51 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;
52 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;
53
54 bool has_exact_sol() const override { return false; }
55 bool is_scalar() const override { return true; }
56
57 void set_parameters(const json &params) override;
58
59 bool is_dimension_dirichet(const int tag, const int dim) const override;
60 bool all_dimensions_dirichlet() const override { return all_dimensions_dirichlet_; }
61
62 private:
64 std::vector<Eigen::Matrix<bool, 1, 3, Eigen::RowMajor>> dirichlet_dimensions_;
65 double rhs_;
67 bool is_all_;
68 };
69 } // namespace problem
70} // namespace polyfem
double val
Definition Assembler.cpp:86
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
bool has_exact_sol() const override
void init(const mesh::Mesh &mesh) override
bool is_dimension_dirichet(const int tag, const int dim) const override
bool is_scalar() const override
void set_parameters(const json &params) 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< Eigen::Matrix< bool, 1, 3, Eigen::RowMajor > > dirichlet_dimensions_
void rhs(const assembler::Assembler &assembler, const Eigen::MatrixXd &pts, 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
bool is_rhs_zero() const override
bool all_dimensions_dirichlet() const override
std::vector< bool > dirichlet_
void load(const std::string &path)
double neumann_interpolate(const int p_id, const Eigen::MatrixXd &uv) const
std::vector< int > raw_ids_
std::vector< bool > raw_dirichlet_
double interpolate(const int p_id, const Eigen::MatrixXd &uv, bool is_dirichlet) const
void init(const mesh::Mesh &mesh)
std::vector< std::vector< double > > raw_data_
double dirichlet_interpolate(const int p_id, const Eigen::MatrixXd &uv) const
std::vector< Eigen::VectorXd > data_
nlohmann::json json
Definition Common.hpp:9