PolyFEM
Loading...
Searching...
No Matches
Evaluator.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4
9
11
12namespace polyfem::io
13{
15 {
16 private:
18
19 public:
31 const mesh::Mesh &mesh,
32 const bool is_problem_scalar,
33 const std::vector<basis::ElementBases> &bases,
34 const std::vector<basis::ElementBases> &gbases,
35 const Eigen::VectorXi &disc_orders,
36 const assembler::Assembler &assembler,
37 const Eigen::MatrixXd &fun,
38 const double t,
39 Eigen::MatrixXd &result,
40 Eigen::VectorXd &von_mises);
41
55 static void interpolate_function(
56 const mesh::Mesh &mesh,
57 const bool is_problem_scalar,
58 const std::vector<basis::ElementBases> &bases,
59 const Eigen::VectorXi &disc_orders,
60 const std::map<int, Eigen::MatrixXd> &polys,
61 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
62 const utils::RefElementSampler &sampler,
63 const int n_points,
64 const Eigen::MatrixXd &fun,
65 Eigen::MatrixXd &result,
66 const bool use_sampler,
67 const bool boundary_only);
68
69 static void mark_flipped_cells(
70 const mesh::Mesh &mesh,
71 const std::vector<basis::ElementBases> &gbasis,
72 const std::vector<basis::ElementBases> &basis,
73 const Eigen::VectorXi &disc_orders,
74 const std::map<int, Eigen::MatrixXd> &polys,
75 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
76 const utils::RefElementSampler &sampler,
77 const int n_points,
78 const Eigen::MatrixXd &fun,
79 Eigen::Vector<bool, -1> &result,
80 const bool use_sampler,
81 const bool boundary_only);
82
96 static void interpolate_function(
97 const mesh::Mesh &mesh,
98 const int actual_dim,
99 const std::vector<basis::ElementBases> &basis,
100 const Eigen::VectorXi &disc_orders,
101 const std::map<int, Eigen::MatrixXd> &polys,
102 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
103 const utils::RefElementSampler &sampler,
104 const int n_points,
105 const Eigen::MatrixXd &fun,
106 Eigen::MatrixXd &result,
107 const bool use_sampler,
108 const bool boundary_only);
109
120 static void interpolate_at_local_vals(
121 const mesh::Mesh &mesh,
122 const bool is_problem_scalar,
123 const std::vector<basis::ElementBases> &bases,
124 const std::vector<basis::ElementBases> &gbases,
125 const int el_index,
126 const Eigen::MatrixXd &local_pts,
127 const Eigen::MatrixXd &fun,
128 Eigen::MatrixXd &result,
129 Eigen::MatrixXd &result_grad);
130
142 static void interpolate_at_local_vals(
143 const mesh::Mesh &mesh,
144 const int actual_dim,
145 const std::vector<basis::ElementBases> &bases,
146 const std::vector<basis::ElementBases> &gbases,
147 const int el_index,
148 const Eigen::MatrixXd &local_pts,
149 const Eigen::MatrixXd &fun,
150 Eigen::MatrixXd &result,
151 Eigen::MatrixXd &result_grad);
152
153 static void interpolate_at_local_vals(
154 const int el_index,
155 const int dim,
156 const int actual_dim,
158 const Eigen::MatrixXd &fun,
159 Eigen::MatrixXd &result,
160 Eigen::MatrixXd &result_grad);
161
177 const mesh::Mesh &mesh,
178 const bool is_problem_scalar,
179 const std::vector<basis::ElementBases> &bases,
180 const std::vector<basis::ElementBases> &gbases,
181 const Eigen::VectorXi &disc_orders,
182 const std::map<int, Eigen::MatrixXd> &polys,
183 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
184 const assembler::Assembler &assembler,
185 const utils::RefElementSampler &sampler,
186 const Eigen::MatrixXd &fun,
187 const double t,
188 const bool use_sampler,
189 const bool boundary_only);
190
206 static void compute_scalar_value(
207 const mesh::Mesh &mesh,
208 const bool is_problem_scalar,
209 const std::vector<basis::ElementBases> &bases,
210 const std::vector<basis::ElementBases> &gbases,
211 const Eigen::VectorXi &disc_orders,
212 const std::map<int, Eigen::MatrixXd> &polys,
213 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
214 const assembler::Assembler &assembler,
215 const utils::RefElementSampler &sampler,
216 const int n_points,
217 const Eigen::MatrixXd &fun,
218 const double t,
219 std::vector<assembler::Assembler::NamedMatrix> &result,
220 const bool use_sampler,
221 const bool boundary_only);
222
241 static void average_grad_based_function(
242 const mesh::Mesh &mesh,
243 const bool is_problem_scalar,
244 const int n_bases,
245 const std::vector<basis::ElementBases> &bases,
246 const std::vector<basis::ElementBases> &gbases,
247 const Eigen::VectorXi &disc_orders,
248 const std::map<int, Eigen::MatrixXd> &polys,
249 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
250 const assembler::Assembler &assembler,
251 const utils::RefElementSampler &sampler,
252 const double t,
253 const int n_points,
254 const Eigen::MatrixXd &fun,
255 std::vector<assembler::Assembler::NamedMatrix> &result_scalar,
256 std::vector<assembler::Assembler::NamedMatrix> &result_tensor,
257 const bool use_sampler,
258 const bool boundary_only);
259
275 static void compute_tensor_value(
276 const mesh::Mesh &mesh,
277 const bool is_problem_scalar,
278 const std::vector<basis::ElementBases> &bases,
279 const std::vector<basis::ElementBases> &gbases,
280 const Eigen::VectorXi &disc_orders,
281 const std::map<int, Eigen::MatrixXd> &polys,
282 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
283 const assembler::Assembler &assembler,
284 const utils::RefElementSampler &sampler,
285 const int n_points,
286 const Eigen::MatrixXd &fun,
287 const double t,
288 std::vector<assembler::Assembler::NamedMatrix> &result,
289 const bool use_sampler,
290 const bool boundary_only);
291
303 const mesh::Mesh &mesh,
304 const bool is_problem_scalar,
305 const std::vector<basis::ElementBases> &bases,
306 const std::vector<basis::ElementBases> &gbases,
307 const Eigen::MatrixXd &pts,
308 const Eigen::MatrixXi &faces,
309 const Eigen::MatrixXd &fun,
310 const bool compute_avg,
311 Eigen::MatrixXd &result);
312
313 static Eigen::MatrixXd generate_linear_field(
314 const int n_bases,
315 const std::shared_ptr<mesh::MeshNodes> mesh_nodes,
316 const Eigen::MatrixXd &grad);
317
318 static Eigen::MatrixXd get_bases_position(
319 const int n_bases,
320 const std::shared_ptr<mesh::MeshNodes> mesh_nodes);
321
322 static Eigen::VectorXd integrate_function(
323 const std::vector<basis::ElementBases> &bases,
324 const std::vector<basis::ElementBases> &gbases,
325 const Eigen::MatrixXd &fun,
326 const int dim,
327 const int actual_dim);
328 };
329} // namespace polyfem::io
ElementAssemblyValues vals
Definition Assembler.cpp:22
std::vector< Eigen::VectorXi > faces
stores per element basis values at given quadrature points and geometric mapping
static Eigen::MatrixXd generate_linear_field(const int n_bases, const std::shared_ptr< mesh::MeshNodes > mesh_nodes, const Eigen::MatrixXd &grad)
bool check_scalar_value(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const assembler::Assembler &assembler, const utils::RefElementSampler &sampler, const Eigen::MatrixXd &fun, const double t, const bool use_sampler, const bool boundary_only)
checks if mises are not nan
static void interpolate_at_local_vals(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const int el_index, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &fun, Eigen::MatrixXd &result, Eigen::MatrixXd &result_grad)
interpolate solution and gradient at element (calls interpolate_at_local_vals with sol)
static void compute_stress_at_quadrature_points(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const assembler::Assembler &assembler, const Eigen::MatrixXd &fun, const double t, Eigen::MatrixXd &result, Eigen::VectorXd &von_mises)
compute von mises stress at quadrature points for the function fun, also compute the interpolated fun...
static Eigen::MatrixXd get_bases_position(const int n_bases, const std::shared_ptr< mesh::MeshNodes > mesh_nodes)
static void interpolate_boundary_function(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::MatrixXd &pts, const Eigen::MatrixXi &faces, const Eigen::MatrixXd &fun, const bool compute_avg, Eigen::MatrixXd &result)
computes integrated solution (fun) per surface face.
Definition Evaluator.cpp:57
static void mark_flipped_cells(const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &gbasis, const std::vector< basis::ElementBases > &basis, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const utils::RefElementSampler &sampler, const int n_points, const Eigen::MatrixXd &fun, Eigen::Vector< bool, -1 > &result, const bool use_sampler, const bool boundary_only)
static void interpolate_function(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const utils::RefElementSampler &sampler, const int n_points, const Eigen::MatrixXd &fun, Eigen::MatrixXd &result, const bool use_sampler, const bool boundary_only)
interpolate the function fun.
static void average_grad_based_function(const mesh::Mesh &mesh, const bool is_problem_scalar, const int n_bases, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const assembler::Assembler &assembler, const utils::RefElementSampler &sampler, const double t, const int n_points, const Eigen::MatrixXd &fun, std::vector< assembler::Assembler::NamedMatrix > &result_scalar, std::vector< assembler::Assembler::NamedMatrix > &result_tensor, const bool use_sampler, const bool boundary_only)
computes scalar quantity of funtion (ie von mises for elasticity and norm of velocity for fluid) the ...
static Eigen::VectorXd integrate_function(const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::MatrixXd &fun, const int dim, const int actual_dim)
static void compute_scalar_value(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const assembler::Assembler &assembler, const utils::RefElementSampler &sampler, const int n_points, const Eigen::MatrixXd &fun, const double t, std::vector< assembler::Assembler::NamedMatrix > &result, const bool use_sampler, const bool boundary_only)
computes scalar quantity of funtion (ie von mises for elasticity and norm of velocity for fluid)
static void compute_tensor_value(const mesh::Mesh &mesh, const bool is_problem_scalar, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXi &disc_orders, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const assembler::Assembler &assembler, const utils::RefElementSampler &sampler, const int n_points, const Eigen::MatrixXd &fun, const double t, std::vector< assembler::Assembler::NamedMatrix > &result, const bool use_sampler, const bool boundary_only)
compute tensor quantity (ie stress tensor or velocy)
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39