PolyFEM
Loading...
Searching...
No Matches
OutData.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
4
6
8
10
12
13#include <paraviewo/ParaviewWriter.hpp>
14#include <paraviewo/VTUWriter.hpp>
15#include <paraviewo/HDF5VTUWriter.hpp>
16
18
19#include <Eigen/Dense>
20
21#include <ipc/potentials/normal_adhesion_potential.hpp>
22#include <ipc/potentials/tangential_adhesion_potential.hpp>
23
24namespace polyfem
25{
26 class State;
27}
28
29namespace polyfem::io
30{
33 {
34 public:
37 {
38 std::vector<std::string> fields; // fields to export, empty means all
39
40 bool volume;
41 bool surface;
42 bool wire;
43 bool points;
48 bool forces;
50
61 bool nodes;
62
65
67
73 ExportOptions(const json &args,
74 const bool is_mesh_linear,
75 const bool mesh_has_prisms,
76 const bool is_problem_scalar);
77
80 inline std::string file_extension() const
81 {
82 return use_hdf5 ? ".hdf" : ".vtu";
83 }
84
85 bool export_field(const std::string &field) const;
86 };
87
97 static void extract_boundary_mesh(
98 const mesh::Mesh &mesh,
99 const int n_bases,
100 const std::vector<basis::ElementBases> &bases,
101 const std::vector<mesh::LocalBoundary> &total_local_boundary,
102 Eigen::MatrixXd &node_positions,
103 Eigen::MatrixXi &boundary_edges,
104 Eigen::MatrixXi &boundary_triangles,
105 std::vector<Eigen::Triplet<double>> &displacement_map_entries);
106
110 void init_sampler(const polyfem::mesh::Mesh &mesh, const double vismesh_rel_area);
111
115 void build_grid(const polyfem::mesh::Mesh &mesh, const double spacing);
116
131 void export_data(
132 const State &state,
133 const Eigen::MatrixXd &sol,
134 const Eigen::MatrixXd &pressure,
135 const bool is_time_dependent,
136 const double tend_in,
137 const double dt,
138 const ExportOptions &opts,
139 const std::string &vis_mesh_path,
140 const std::string &nodes_path,
141 const std::string &solution_path,
142 const std::string &stress_path,
143 const std::string &mises_path,
144 const bool is_contact_enabled) const;
145
155 void save_vtu(const std::string &path,
156 const State &state,
157 const Eigen::MatrixXd &sol,
158 const Eigen::MatrixXd &pressure,
159 const double t,
160 const double dt,
161 const ExportOptions &opts,
162 const bool is_contact_enabled) const;
163
172 void save_volume(const std::string &path,
173 const State &state,
174 const Eigen::MatrixXd &sol,
175 const Eigen::MatrixXd &pressure,
176 const double t,
177 const double dt,
178 const ExportOptions &opts) const;
179
189 void save_surface(const std::string &export_surface,
190 const State &state,
191 const Eigen::MatrixXd &sol,
192 const Eigen::MatrixXd &pressure,
193 const double t,
194 const double dt_in,
195 const ExportOptions &opts,
196 const bool is_contact_enabled) const;
197
208 const std::string &export_surface,
209 const State &state,
210 const Eigen::MatrixXd &sol,
211 const Eigen::MatrixXd &pressure,
212 const double t,
213 const double dt_in,
214 const ExportOptions &opts,
215 const bool is_contact_enabled) const;
216
223 void save_wire(const std::string &name,
224 const State &state,
225 const Eigen::MatrixXd &sol,
226 const double t,
227 const ExportOptions &opts) const;
228
234 void save_points(
235 const std::string &path,
236 const State &state,
237 const Eigen::MatrixXd &sol,
238 const ExportOptions &opts) const;
239
247 void save_pvd(const std::string &name, const std::function<std::string(int)> &vtu_names,
248 int time_steps, double t0, double dt, int skip_frame = 1) const;
249
250 private:
253
255 Eigen::MatrixXd grid_points;
257 Eigen::MatrixXi grid_points_to_elements;
259 Eigen::MatrixXd grid_points_bc;
260
276 const mesh::Mesh &mesh,
277 const std::vector<basis::ElementBases> &bases,
278 const std::vector<basis::ElementBases> &gbases,
279 const std::vector<mesh::LocalBoundary> &total_local_boundary,
280 const Eigen::MatrixXd &solution,
281 const int problem_dim,
282 Eigen::MatrixXd &boundary_vis_vertices,
283 Eigen::MatrixXd &boundary_vis_local_vertices,
284 Eigen::MatrixXi &boundary_vis_elements,
285 Eigen::MatrixXi &boundary_vis_elements_ids,
286 Eigen::MatrixXi &boundary_vis_primitive_ids,
287 Eigen::MatrixXd &boundary_vis_normals,
288 Eigen::MatrixXd &displaced_boundary_vis_normals) const;
289
304 void build_vis_mesh(
305 const mesh::Mesh &mesh,
306 const Eigen::VectorXi &disc_orders,
307 const std::vector<basis::ElementBases> &gbases,
308 const std::map<int, Eigen::MatrixXd> &polys,
309 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
310 const bool boundary_only,
311 Eigen::MatrixXd &points,
312 Eigen::MatrixXi &tets,
313 Eigen::MatrixXi &el_id,
314 Eigen::MatrixXd &discr) const;
315
328 const mesh::Mesh &mesh,
329 const Eigen::VectorXi &disc_orders,
330 const Eigen::VectorXi &disc_ordersq,
331 const std::vector<basis::ElementBases> &bases,
332 Eigen::MatrixXd &points,
333 std::vector<std::vector<int>> &elements,
334 Eigen::MatrixXi &el_id,
335 Eigen::MatrixXd &discr) const;
336
338 const State &state,
339 const Eigen::MatrixXd &points,
340 const ExportOptions &opts,
341 const std::string &name,
342 const Eigen::VectorXd &field,
343 paraviewo::ParaviewWriter &writer) const;
344 };
345
372
375 {
376 public:
378 Eigen::Vector4d spectrum;
379
383
385 double mesh_size;
390
393
397
399 double max_angle;
402
405
428
437 void compute_errors(const int n_bases,
438 const std::vector<polyfem::basis::ElementBases> &bases,
439 const std::vector<polyfem::basis::ElementBases> &gbases,
440 const polyfem::mesh::Mesh &mesh,
441 const assembler::Problem &problem,
442 const double tend,
443 const Eigen::MatrixXd &sol);
444
447 void compute_mesh_stats(const polyfem::mesh::Mesh &mesh);
448
456 void compute_mesh_size(const polyfem::mesh::Mesh &mesh_in, const std::vector<polyfem::basis::ElementBases> &bases_in, const int n_samples, const bool use_curved_mesh_size);
457
459 void reset();
460
464 void count_flipped_elements(const polyfem::mesh::Mesh &mesh, const std::vector<polyfem::basis::ElementBases> &gbases);
465
468
483 void save_json(const nlohmann::json &args,
484 const int n_bases, const int n_pressure_bases,
485 const Eigen::MatrixXd &sol,
486 const mesh::Mesh &mesh,
487 const Eigen::VectorXi &disc_orders,
488 const Eigen::VectorXi &disc_ordersq,
489 const assembler::Problem &problem,
490 const OutRuntimeData &runtime,
491 const std::string &formulation,
492 const bool isoparametric,
493 const int sol_at_node_id,
494 nlohmann::json &j);
495 };
496
498 {
499 public:
500 EnergyCSVWriter(const std::string &path, const solver::SolveData &solve_data);
502
503 void write(const int i, const Eigen::MatrixXd &sol);
504
505 protected:
506 std::ofstream file;
508 };
509
511 {
512 public:
513 RuntimeStatsCSVWriter(const std::string &path, const State &state, const double t0, const double dt);
515
516 void write(const int t, const double forward, const double remeshing, const double global_relaxation, const Eigen::MatrixXd &sol);
517
518 protected:
519 std::ofstream file;
520 const State &state;
521 const double t0;
522 const double dt;
526 };
527} // namespace polyfem::io
main class that contains the polyfem solver and all its state
Definition State.hpp:79
const solver::SolveData & solve_data
Definition OutData.hpp:507
void write(const int i, const Eigen::MatrixXd &sol)
Definition OutData.cpp:3197
Utilies related to export of geometry.
Definition OutData.hpp:33
void build_vis_mesh(const mesh::Mesh &mesh, const Eigen::VectorXi &disc_orders, const std::vector< basis::ElementBases > &gbases, const std::map< int, Eigen::MatrixXd > &polys, const std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > &polys_3d, const bool boundary_only, Eigen::MatrixXd &points, Eigen::MatrixXi &tets, Eigen::MatrixXi &el_id, Eigen::MatrixXd &discr) const
builds visualzation mesh, upsampled mesh used for visualization the visualization mesh is a dense mes...
Definition OutData.cpp:753
void save_volume_vector_field(const State &state, const Eigen::MatrixXd &points, const ExportOptions &opts, const std::string &name, const Eigen::VectorXd &field, paraviewo::ParaviewWriter &writer) const
Definition OutData.cpp:1901
Eigen::MatrixXd grid_points_bc
grid mesh boundaries
Definition OutData.hpp:259
Eigen::MatrixXd grid_points
grid mesh points to export solution sampled on a grid
Definition OutData.hpp:255
void build_vis_boundary_mesh(const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const std::vector< mesh::LocalBoundary > &total_local_boundary, const Eigen::MatrixXd &solution, const int problem_dim, Eigen::MatrixXd &boundary_vis_vertices, Eigen::MatrixXd &boundary_vis_local_vertices, Eigen::MatrixXi &boundary_vis_elements, Eigen::MatrixXi &boundary_vis_elements_ids, Eigen::MatrixXi &boundary_vis_primitive_ids, Eigen::MatrixXd &boundary_vis_normals, Eigen::MatrixXd &displaced_boundary_vis_normals) const
builds the boundary mesh for visualization
Definition OutData.cpp:498
void save_points(const std::string &path, const State &state, const Eigen::MatrixXd &sol, const ExportOptions &opts) const
saves the nodal values
Definition OutData.cpp:2509
void build_grid(const polyfem::mesh::Mesh &mesh, const double spacing)
builds the grid to export the solution
Definition OutData.cpp:2575
void save_contact_surface(const std::string &export_surface, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt_in, const ExportOptions &opts, const bool is_contact_enabled) const
saves the surface vtu file for for constact quantites, eg contact or friction forces
Definition OutData.cpp:2117
static void extract_boundary_mesh(const mesh::Mesh &mesh, const int n_bases, const std::vector< basis::ElementBases > &bases, const std::vector< mesh::LocalBoundary > &total_local_boundary, Eigen::MatrixXd &node_positions, Eigen::MatrixXi &boundary_edges, Eigen::MatrixXi &boundary_triangles, std::vector< Eigen::Triplet< double > > &displacement_map_entries)
extracts the boundary mesh
Definition OutData.cpp:146
void save_volume(const std::string &path, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt, const ExportOptions &opts) const
saves the volume vtu file
Definition OutData.cpp:1331
void save_pvd(const std::string &name, const std::function< std::string(int)> &vtu_names, int time_steps, double t0, double dt, int skip_frame=1) const
save a PVD of a time dependent simulation
Definition OutData.cpp:2562
void build_high_order_vis_mesh(const mesh::Mesh &mesh, const Eigen::VectorXi &disc_orders, const Eigen::VectorXi &disc_ordersq, const std::vector< basis::ElementBases > &bases, Eigen::MatrixXd &points, std::vector< std::vector< int > > &elements, Eigen::MatrixXi &el_id, Eigen::MatrixXd &discr) const
builds high-der visualzation mesh per element all disconnected it also retuns the mapping to element ...
Definition OutData.cpp:902
void save_wire(const std::string &name, const State &state, const Eigen::MatrixXd &sol, const double t, const ExportOptions &opts) const
saves the wireframe
Definition OutData.cpp:2301
void save_vtu(const std::string &path, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt, const ExportOptions &opts, const bool is_contact_enabled) const
saves the vtu file for time t
Definition OutData.cpp:1247
void init_sampler(const polyfem::mesh::Mesh &mesh, const double vismesh_rel_area)
unitalize the ref element sampler
Definition OutData.cpp:2570
void export_data(const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const bool is_time_dependent, const double tend_in, const double dt, const ExportOptions &opts, const std::string &vis_mesh_path, const std::string &nodes_path, const std::string &solution_path, const std::string &stress_path, const std::string &mises_path, const bool is_contact_enabled) const
exports everytihng, txt, vtu, etc
Definition OutData.cpp:1059
void save_surface(const std::string &export_surface, const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure, const double t, const double dt_in, const ExportOptions &opts, const bool is_contact_enabled) const
saves the surface vtu file for for surface quantites, eg traction forces
Definition OutData.cpp:1926
Eigen::MatrixXi grid_points_to_elements
grid mesh mapping to fe elements
Definition OutData.hpp:257
utils::RefElementSampler ref_element_sampler
used to sample the solution
Definition OutData.hpp:252
stores all runtime data
Definition OutData.hpp:348
double loading_mesh_time
time to load the mesh
Definition OutData.hpp:353
double total_time()
computes total time
Definition OutData.hpp:367
double assembling_stiffness_mat_time
time to assembly
Definition OutData.hpp:357
double assigning_rhs_time
time to computing the rhs
Definition OutData.hpp:361
double assembling_mass_mat_time
time to assembly mass
Definition OutData.hpp:359
double building_basis_time
time to construct the basis
Definition OutData.hpp:351
double solving_time
time to solve
Definition OutData.hpp:363
double computing_poly_basis_time
time to build the polygonal/polyhedral bases
Definition OutData.hpp:355
all stats from polyfem
Definition OutData.hpp:375
int n_flipped
number of flipped elements, compute only when using count_flipped_els (false by default)
Definition OutData.hpp:404
int simplex_count
statiscs on the mesh (simplices)
Definition OutData.hpp:407
json solver_info
information of the solver, eg num iteration, time, errors, etc the informations varies depending on t...
Definition OutData.hpp:382
int simple_singular_count
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition OutData.hpp:415
Eigen::Vector4d spectrum
spectrum of the stiffness matrix, enable only if POLYSOLVE_WITH_SPECTRA is ON (off by default)
Definition OutData.hpp:378
int prism_count
statiscs on the mesh (simplices)
Definition OutData.hpp:409
void count_flipped_elements(const polyfem::mesh::Mesh &mesh, const std::vector< polyfem::basis::ElementBases > &gbases)
counts the number of flipped elements
Definition OutData.cpp:2760
double l2_err
errors, lp_err is in fact an L8 error
Definition OutData.hpp:392
void save_json(const nlohmann::json &args, const int n_bases, const int n_pressure_bases, const Eigen::MatrixXd &sol, const mesh::Mesh &mesh, const Eigen::VectorXi &disc_orders, const Eigen::VectorXi &disc_ordersq, const assembler::Problem &problem, const OutRuntimeData &runtime, const std::string &formulation, const bool isoparametric, const int sol_at_node_id, nlohmann::json &j)
saves the output statistic to a json object
Definition OutData.cpp:3043
int undefined_count
statiscs on the mesh (not quad/hex simplex), see Polyspline paper for desciption
Definition OutData.hpp:425
void compute_errors(const int n_bases, const std::vector< polyfem::basis::ElementBases > &bases, const std::vector< polyfem::basis::ElementBases > &gbases, const polyfem::mesh::Mesh &mesh, const assembler::Problem &problem, const double tend, const Eigen::MatrixXd &sol)
compute errors
Definition OutData.cpp:2805
void compute_mesh_size(const polyfem::mesh::Mesh &mesh_in, const std::vector< polyfem::basis::ElementBases > &bases_in, const int n_samples, const bool use_curved_mesh_size)
computes the mesh size, it samples every edges n_samples times uses curved_mesh_size (false by defaul...
Definition OutData.cpp:2673
double average_edge_length
avg edge lenght
Definition OutData.hpp:389
double max_angle
statiscs on angle, compute only when using p_ref (false by default)
Definition OutData.hpp:399
int non_regular_boundary_count
statiscs on the mesh (irregular boundary quad/hex part of the mesh), see Polyspline paper for descipt...
Definition OutData.hpp:421
int non_regular_count
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition OutData.hpp:423
int regular_boundary_count
statiscs on the mesh (regular quad/hex boundary part of the mesh), see Polyspline paper for desciptio...
Definition OutData.hpp:413
long long nn_zero
non zeros and sytem matrix size num dof is the total dof in the system
Definition OutData.hpp:396
double mesh_size
max edge lenght
Definition OutData.hpp:385
void reset()
clears all stats
Definition OutData.cpp:2751
int multi_singular_boundary_count
statiscs on the mesh (irregular boundary quad/hex part of the mesh), see Polyspline paper for descipt...
Definition OutData.hpp:427
void compute_mesh_stats(const polyfem::mesh::Mesh &mesh)
compute stats (counts els type, mesh lenght, etc), step 1 of solve
Definition OutData.cpp:2966
double sigma_max
statiscs on tri/tet quality, compute only when using p_ref (false by default)
Definition OutData.hpp:401
int boundary_count
statiscs on the mesh (boundary quads/hexs), see Polyspline paper for desciption
Definition OutData.hpp:419
int regular_count
statiscs on the mesh (regular quad/hex part of the mesh), see Polyspline paper for desciption
Definition OutData.hpp:411
double min_edge_length
min edge lenght
Definition OutData.hpp:387
int multi_singular_count
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition OutData.hpp:417
void write(const int t, const double forward, const double remeshing, const double global_relaxation, const Eigen::MatrixXd &sol)
Definition OutData.cpp:3223
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:40
class to store time stepping data
Definition SolveData.hpp:59
nlohmann::json json
Definition Common.hpp:9
std::string file_extension() const
return the extension of the output paraview files depending on use_hdf5
Definition OutData.hpp:80
std::vector< std::string > fields
Definition OutData.hpp:38
bool export_field(const std::string &field) const
Definition OutData.cpp:1201