Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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:
32 const mesh::Mesh &mesh,
33 const bool is_problem_scalar,
34 const std::vector<basis::ElementBases> &bases,
35 const std::vector<basis::ElementBases> &gbases,
36 const Eigen::VectorXi &disc_orders,
37 const Eigen::VectorXi &disc_ordersq,
38 const assembler::Assembler &assembler,
39 const Eigen::MatrixXd &fun,
40 const double t,
41 Eigen::MatrixXd &result,
42 Eigen::VectorXd &von_mises);
43
58 static void interpolate_function(
59 const mesh::Mesh &mesh,
60 const bool is_problem_scalar,
61 const std::vector<basis::ElementBases> &bases,
62 const Eigen::VectorXi &disc_orders,
63 const Eigen::VectorXi &disc_ordersq,
64 const std::map<int, Eigen::MatrixXd> &polys,
65 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
66 const utils::RefElementSampler &sampler,
67 const int n_points,
68 const Eigen::MatrixXd &fun,
69 Eigen::MatrixXd &result,
70 const bool use_sampler,
71 const bool boundary_only);
72
73 static void mark_flipped_cells(
74 const mesh::Mesh &mesh,
75 const std::vector<basis::ElementBases> &gbasis,
76 const std::vector<basis::ElementBases> &basis,
77 const Eigen::VectorXi &disc_orders,
78 const std::map<int, Eigen::MatrixXd> &polys,
79 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
80 const utils::RefElementSampler &sampler,
81 const int n_points,
82 const Eigen::MatrixXd &fun,
83 Eigen::Vector<bool, -1> &result,
84 const bool use_sampler,
85 const bool boundary_only);
86
101 static void interpolate_function(
102 const mesh::Mesh &mesh,
103 const int actual_dim,
104 const std::vector<basis::ElementBases> &basis,
105 const Eigen::VectorXi &disc_orders,
106 const Eigen::VectorXi &disc_ordersq,
107 const std::map<int, Eigen::MatrixXd> &polys,
108 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
109 const utils::RefElementSampler &sampler,
110 const int n_points,
111 const Eigen::MatrixXd &fun,
112 Eigen::MatrixXd &result,
113 const bool use_sampler,
114 const bool boundary_only);
115
126 static void interpolate_at_local_vals(
127 const mesh::Mesh &mesh,
128 const bool is_problem_scalar,
129 const std::vector<basis::ElementBases> &bases,
130 const std::vector<basis::ElementBases> &gbases,
131 const int el_index,
132 const Eigen::MatrixXd &local_pts,
133 const Eigen::MatrixXd &fun,
134 Eigen::MatrixXd &result,
135 Eigen::MatrixXd &result_grad);
136
148 static void interpolate_at_local_vals(
149 const mesh::Mesh &mesh,
150 const int actual_dim,
151 const std::vector<basis::ElementBases> &bases,
152 const std::vector<basis::ElementBases> &gbases,
153 const int el_index,
154 const Eigen::MatrixXd &local_pts,
155 const Eigen::MatrixXd &fun,
156 Eigen::MatrixXd &result,
157 Eigen::MatrixXd &result_grad);
158
159 static void interpolate_at_local_vals(
160 const int el_index,
161 const int dim,
162 const int actual_dim,
164 const Eigen::MatrixXd &fun,
165 Eigen::MatrixXd &result,
166 Eigen::MatrixXd &result_grad);
167
184 const mesh::Mesh &mesh,
185 const bool is_problem_scalar,
186 const std::vector<basis::ElementBases> &bases,
187 const std::vector<basis::ElementBases> &gbases,
188 const Eigen::VectorXi &disc_orders,
189 const Eigen::VectorXi &disc_ordersq,
190 const std::map<int, Eigen::MatrixXd> &polys,
191 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
192 const assembler::Assembler &assembler,
193 const utils::RefElementSampler &sampler,
194 const Eigen::MatrixXd &fun,
195 const double t,
196 const bool use_sampler,
197 const bool boundary_only);
198
215 static void compute_scalar_value(
216 const mesh::Mesh &mesh,
217 const bool is_problem_scalar,
218 const std::vector<basis::ElementBases> &bases,
219 const std::vector<basis::ElementBases> &gbases,
220 const Eigen::VectorXi &disc_orders,
221 const Eigen::VectorXi &disc_ordersq,
222 const std::map<int, Eigen::MatrixXd> &polys,
223 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
224 const assembler::Assembler &assembler,
225 const utils::RefElementSampler &sampler,
226 const int n_points,
227 const Eigen::MatrixXd &fun,
228 const double t,
229 std::vector<assembler::Assembler::NamedMatrix> &result,
230 const bool use_sampler,
231 const bool boundary_only);
232
252 static void average_grad_based_function(
253 const mesh::Mesh &mesh,
254 const bool is_problem_scalar,
255 const int n_bases,
256 const std::vector<basis::ElementBases> &bases,
257 const std::vector<basis::ElementBases> &gbases,
258 const Eigen::VectorXi &disc_orders,
259 const Eigen::VectorXi &disc_ordersq,
260 const std::map<int, Eigen::MatrixXd> &polys,
261 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
262 const assembler::Assembler &assembler,
263 const utils::RefElementSampler &sampler,
264 const double t,
265 const int n_points,
266 const Eigen::MatrixXd &fun,
267 std::vector<assembler::Assembler::NamedMatrix> &result_scalar,
268 std::vector<assembler::Assembler::NamedMatrix> &result_tensor,
269 const bool use_sampler,
270 const bool boundary_only);
271
288 static void compute_tensor_value(
289 const mesh::Mesh &mesh,
290 const bool is_problem_scalar,
291 const std::vector<basis::ElementBases> &bases,
292 const std::vector<basis::ElementBases> &gbases,
293 const Eigen::VectorXi &disc_orders,
294 const Eigen::VectorXi &disc_ordersq,
295 const std::map<int, Eigen::MatrixXd> &polys,
296 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
297 const assembler::Assembler &assembler,
298 const utils::RefElementSampler &sampler,
299 const int n_points,
300 const Eigen::MatrixXd &fun,
301 const double t,
302 std::vector<assembler::Assembler::NamedMatrix> &result,
303 const bool use_sampler,
304 const bool boundary_only);
305
317 const mesh::Mesh &mesh,
318 const bool is_problem_scalar,
319 const std::vector<basis::ElementBases> &bases,
320 const std::vector<basis::ElementBases> &gbases,
321 const Eigen::MatrixXd &pts,
322 const Eigen::MatrixXi &faces,
323 const Eigen::MatrixXd &fun,
324 const bool compute_avg,
325 Eigen::MatrixXd &result);
326
327 static Eigen::MatrixXd generate_linear_field(
328 const int n_bases,
329 const std::shared_ptr<mesh::MeshNodes> mesh_nodes,
330 const Eigen::MatrixXd &grad);
331
332 static Eigen::MatrixXd get_bases_position(
333 const int n_bases,
334 const std::shared_ptr<mesh::MeshNodes> mesh_nodes);
335
336 static Eigen::VectorXd integrate_function(
337 const std::vector<basis::ElementBases> &bases,
338 const std::vector<basis::ElementBases> &gbases,
339 const Eigen::MatrixXd &fun,
340 const int dim,
341 const int actual_dim);
342 };
343} // 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)
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 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 Eigen::VectorXi &disc_ordersq, 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::MatrixXd get_bases_position(const int n_bases, const std::shared_ptr< mesh::MeshNodes > mesh_nodes)
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 Eigen::VectorXi &disc_ordersq, 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 Eigen::VectorXi &disc_ordersq, 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)
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 Eigen::VectorXi &disc_ordersq, 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_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:59
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 Eigen::VectorXi &disc_ordersq, 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 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 Eigen::VectorXi &disc_ordersq, 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::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)
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:40