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
82 static void interpolate_function(
83 const mesh::Mesh &mesh,
84 const int actual_dim,
85 const std::vector<basis::ElementBases> &basis,
86 const Eigen::VectorXi &disc_orders,
87 const std::map<int, Eigen::MatrixXd> &polys,
88 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
89 const utils::RefElementSampler &sampler,
90 const int n_points,
91 const Eigen::MatrixXd &fun,
92 Eigen::MatrixXd &result,
93 const bool use_sampler,
94 const bool boundary_only);
95
106 static void interpolate_at_local_vals(
107 const mesh::Mesh &mesh,
108 const bool is_problem_scalar,
109 const std::vector<basis::ElementBases> &bases,
110 const std::vector<basis::ElementBases> &gbases,
111 const int el_index,
112 const Eigen::MatrixXd &local_pts,
113 const Eigen::MatrixXd &fun,
114 Eigen::MatrixXd &result,
115 Eigen::MatrixXd &result_grad);
116
128 static void interpolate_at_local_vals(
129 const mesh::Mesh &mesh,
130 const int actual_dim,
131 const std::vector<basis::ElementBases> &bases,
132 const std::vector<basis::ElementBases> &gbases,
133 const int el_index,
134 const Eigen::MatrixXd &local_pts,
135 const Eigen::MatrixXd &fun,
136 Eigen::MatrixXd &result,
137 Eigen::MatrixXd &result_grad);
138
139 static void interpolate_at_local_vals(
140 const int el_index,
141 const int dim,
142 const int actual_dim,
144 const Eigen::MatrixXd &fun,
145 Eigen::MatrixXd &result,
146 Eigen::MatrixXd &result_grad);
147
163 const mesh::Mesh &mesh,
164 const bool is_problem_scalar,
165 const std::vector<basis::ElementBases> &bases,
166 const std::vector<basis::ElementBases> &gbases,
167 const Eigen::VectorXi &disc_orders,
168 const std::map<int, Eigen::MatrixXd> &polys,
169 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
170 const assembler::Assembler &assembler,
171 const utils::RefElementSampler &sampler,
172 const Eigen::MatrixXd &fun,
173 const double t,
174 const bool use_sampler,
175 const bool boundary_only);
176
192 static void compute_scalar_value(
193 const mesh::Mesh &mesh,
194 const bool is_problem_scalar,
195 const std::vector<basis::ElementBases> &bases,
196 const std::vector<basis::ElementBases> &gbases,
197 const Eigen::VectorXi &disc_orders,
198 const std::map<int, Eigen::MatrixXd> &polys,
199 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
200 const assembler::Assembler &assembler,
201 const utils::RefElementSampler &sampler,
202 const int n_points,
203 const Eigen::MatrixXd &fun,
204 const double t,
205 std::vector<assembler::Assembler::NamedMatrix> &result,
206 const bool use_sampler,
207 const bool boundary_only);
208
227 static void average_grad_based_function(
228 const mesh::Mesh &mesh,
229 const bool is_problem_scalar,
230 const int n_bases,
231 const std::vector<basis::ElementBases> &bases,
232 const std::vector<basis::ElementBases> &gbases,
233 const Eigen::VectorXi &disc_orders,
234 const std::map<int, Eigen::MatrixXd> &polys,
235 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
236 const assembler::Assembler &assembler,
237 const utils::RefElementSampler &sampler,
238 const double t,
239 const int n_points,
240 const Eigen::MatrixXd &fun,
241 std::vector<assembler::Assembler::NamedMatrix> &result_scalar,
242 std::vector<assembler::Assembler::NamedMatrix> &result_tensor,
243 const bool use_sampler,
244 const bool boundary_only);
245
261 static void compute_tensor_value(
262 const mesh::Mesh &mesh,
263 const bool is_problem_scalar,
264 const std::vector<basis::ElementBases> &bases,
265 const std::vector<basis::ElementBases> &gbases,
266 const Eigen::VectorXi &disc_orders,
267 const std::map<int, Eigen::MatrixXd> &polys,
268 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
269 const assembler::Assembler &assembler,
270 const utils::RefElementSampler &sampler,
271 const int n_points,
272 const Eigen::MatrixXd &fun,
273 const double t,
274 std::vector<assembler::Assembler::NamedMatrix> &result,
275 const bool use_sampler,
276 const bool boundary_only);
277
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::MatrixXd &pts,
294 const Eigen::MatrixXi &faces,
295 const Eigen::MatrixXd &fun,
296 const bool compute_avg,
297 Eigen::MatrixXd &result);
298
299 static Eigen::MatrixXd generate_linear_field(
300 const int n_bases,
301 const std::shared_ptr<mesh::MeshNodes> mesh_nodes,
302 const Eigen::MatrixXd &grad);
303
304 static Eigen::MatrixXd get_bases_position(
305 const int n_bases,
306 const std::shared_ptr<mesh::MeshNodes> mesh_nodes);
307
308 static Eigen::VectorXd integrate_function(
309 const std::vector<basis::ElementBases> &bases,
310 const std::vector<basis::ElementBases> &gbases,
311 const Eigen::MatrixXd &fun,
312 const int dim,
313 const int actual_dim);
314 };
315} // 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:56
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