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
21namespace polyfem
22{
23 class State;
24}
25
26namespace polyfem::io
27{
30 {
31 public:
32 std::string name;
33 Eigen::MatrixXd points;
34 Eigen::MatrixXi connectivity;
35 Eigen::MatrixXd solution;
36 Eigen::MatrixXd pressure;
37 Eigen::MatrixXd exact;
38 Eigen::MatrixXd error;
39 Eigen::MatrixXd scalar_value;
40 Eigen::MatrixXd scalar_value_avg;
41 };
42
45 {
46 public:
49 {
50 bool volume;
51 bool surface;
52 bool wire;
53 bool points;
56 bool forces;
57
68 bool nodes;
69
72
74
76
82 ExportOptions(const json &args,
83 const bool is_mesh_linear,
84 const bool is_problem_scalar,
85 const bool solve_export_to_file);
86
89 inline std::string file_extension() const { return use_hdf5 ? ".hdf" : ".vtu"; }
90 };
91
101 static void extract_boundary_mesh(
102 const mesh::Mesh &mesh,
103 const int n_bases,
104 const std::vector<basis::ElementBases> &bases,
105 const std::vector<mesh::LocalBoundary> &total_local_boundary,
106 Eigen::MatrixXd &node_positions,
107 Eigen::MatrixXi &boundary_edges,
108 Eigen::MatrixXi &boundary_triangles,
109 std::vector<Eigen::Triplet<double>> &displacement_map_entries);
110
114 void init_sampler(const polyfem::mesh::Mesh &mesh, const double vismesh_rel_area);
115
119 void build_grid(const polyfem::mesh::Mesh &mesh, const double spacing);
120
136 void export_data(
137 const State &state,
138 const Eigen::MatrixXd &sol,
139 const Eigen::MatrixXd &pressure,
140 const bool is_time_dependent,
141 const double tend_in,
142 const double dt,
143 const ExportOptions &opts,
144 const std::string &vis_mesh_path,
145 const std::string &nodes_path,
146 const std::string &solution_path,
147 const std::string &stress_path,
148 const std::string &mises_path,
149 const bool is_contact_enabled,
150 std::vector<SolutionFrame> &solution_frames) const;
151
162 void save_vtu(const std::string &path,
163 const State &state,
164 const Eigen::MatrixXd &sol,
165 const Eigen::MatrixXd &pressure,
166 const double t,
167 const double dt,
168 const ExportOptions &opts,
169 const bool is_contact_enabled,
170 std::vector<SolutionFrame> &solution_frames) const;
171
181 void save_volume(const std::string &path,
182 const State &state,
183 const Eigen::MatrixXd &sol,
184 const Eigen::MatrixXd &pressure,
185 const double t,
186 const double dt,
187 const ExportOptions &opts,
188 std::vector<SolutionFrame> &solution_frames) const;
189
200 void save_surface(const std::string &export_surface,
201 const State &state,
202 const Eigen::MatrixXd &sol,
203 const Eigen::MatrixXd &pressure,
204 const double t,
205 const double dt_in,
206 const ExportOptions &opts,
207 const bool is_contact_enabled,
208 std::vector<SolutionFrame> &solution_frames) const;
209
221 const std::string &export_surface,
222 const State &state,
223 const Eigen::MatrixXd &sol,
224 const Eigen::MatrixXd &pressure,
225 const double t,
226 const double dt_in,
227 const ExportOptions &opts,
228 const bool is_contact_enabled,
229 std::vector<SolutionFrame> &solution_frames) const;
230
238 void save_wire(const std::string &name,
239 const State &state,
240 const Eigen::MatrixXd &sol,
241 const double t,
242 const ExportOptions &opts,
243 std::vector<SolutionFrame> &solution_frames) const;
244
251 void save_points(
252 const std::string &path,
253 const State &state,
254 const Eigen::MatrixXd &sol,
255 const ExportOptions &opts,
256 std::vector<SolutionFrame> &solution_frames) const;
257
265 void save_pvd(const std::string &name, const std::function<std::string(int)> &vtu_names,
266 int time_steps, double t0, double dt, int skip_frame = 1) const;
267
268 private:
271
273 Eigen::MatrixXd grid_points;
275 Eigen::MatrixXi grid_points_to_elements;
277 Eigen::MatrixXd grid_points_bc;
278
294 const mesh::Mesh &mesh,
295 const std::vector<basis::ElementBases> &bases,
296 const std::vector<basis::ElementBases> &gbases,
297 const std::vector<mesh::LocalBoundary> &total_local_boundary,
298 const Eigen::MatrixXd &solution,
299 const int problem_dim,
300 Eigen::MatrixXd &boundary_vis_vertices,
301 Eigen::MatrixXd &boundary_vis_local_vertices,
302 Eigen::MatrixXi &boundary_vis_elements,
303 Eigen::MatrixXi &boundary_vis_elements_ids,
304 Eigen::MatrixXi &boundary_vis_primitive_ids,
305 Eigen::MatrixXd &boundary_vis_normals,
306 Eigen::MatrixXd &displaced_boundary_vis_normals) const;
307
322 void build_vis_mesh(
323 const mesh::Mesh &mesh,
324 const Eigen::VectorXi &disc_orders,
325 const std::vector<basis::ElementBases> &gbases,
326 const std::map<int, Eigen::MatrixXd> &polys,
327 const std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> &polys_3d,
328 const bool boundary_only,
329 Eigen::MatrixXd &points,
330 Eigen::MatrixXi &tets,
331 Eigen::MatrixXi &el_id,
332 Eigen::MatrixXd &discr) const;
333
345 const mesh::Mesh &mesh,
346 const Eigen::VectorXi &disc_orders,
347 const std::vector<basis::ElementBases> &bases,
348 Eigen::MatrixXd &points,
349 std::vector<std::vector<int>> &elements,
350 Eigen::MatrixXi &el_id,
351 Eigen::MatrixXd &discr) const;
352
354 const State &state,
355 const Eigen::MatrixXd &points,
356 const ExportOptions &opts,
357 const std::string &name,
358 const Eigen::VectorXd &field,
359 paraviewo::ParaviewWriter &writer) const;
360 };
361
388
391 {
392 public:
394 Eigen::Vector4d spectrum;
395
399
401 double mesh_size;
406
409
413
415 double max_angle;
418
421
442
451 void compute_errors(const int n_bases,
452 const std::vector<polyfem::basis::ElementBases> &bases,
453 const std::vector<polyfem::basis::ElementBases> &gbases,
454 const polyfem::mesh::Mesh &mesh,
455 const assembler::Problem &problem,
456 const double tend,
457 const Eigen::MatrixXd &sol);
458
461 void compute_mesh_stats(const polyfem::mesh::Mesh &mesh);
462
470 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);
471
473 void reset();
474
478 void count_flipped_elements(const polyfem::mesh::Mesh &mesh, const std::vector<polyfem::basis::ElementBases> &gbases);
479
482
496 void save_json(const nlohmann::json &args,
497 const int n_bases, const int n_pressure_bases,
498 const Eigen::MatrixXd &sol,
499 const mesh::Mesh &mesh,
500 const Eigen::VectorXi &disc_orders,
501 const assembler::Problem &problem,
502 const OutRuntimeData &runtime,
503 const std::string &formulation,
504 const bool isoparametric,
505 const int sol_at_node_id,
506 nlohmann::json &j);
507 };
508
510 {
511 public:
512 EnergyCSVWriter(const std::string &path, const solver::SolveData &solve_data);
514
515 void write(const int i, const Eigen::MatrixXd &sol);
516
517 protected:
518 std::ofstream file;
520 };
521
523 {
524 public:
525 RuntimeStatsCSVWriter(const std::string &path, const State &state, const double t0, const double dt);
527
528 void write(const int t, const double forward, const double remeshing, const double global_relaxation, const Eigen::MatrixXd &sol);
529
530 protected:
531 std::ofstream file;
532 const State &state;
533 const double t0;
534 const double dt;
538 };
539} // 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:519
void write(const int i, const Eigen::MatrixXd &sol)
Definition OutData.cpp:2927
Utilies related to export of geometry.
Definition OutData.hpp:45
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:640
void build_high_order_vis_mesh(const mesh::Mesh &mesh, const Eigen::VectorXi &disc_orders, 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:783
void save_points(const std::string &path, const State &state, const Eigen::MatrixXd &sol, const ExportOptions &opts, std::vector< SolutionFrame > &solution_frames) const
saves the nodal values
Definition OutData.cpp:2245
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:1735
Eigen::MatrixXd grid_points_bc
grid mesh boundaries
Definition OutData.hpp:277
Eigen::MatrixXd grid_points
grid mesh points to export solution sampled on a grid
Definition OutData.hpp:273
void save_wire(const std::string &name, const State &state, const Eigen::MatrixXd &sol, const double t, const ExportOptions &opts, std::vector< SolutionFrame > &solution_frames) const
saves the wireframe
Definition OutData.cpp:2057
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, std::vector< SolutionFrame > &solution_frames) const
saves the surface vtu file for for surface quantites, eg traction forces
Definition OutData.cpp:1764
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:392
void build_grid(const polyfem::mesh::Mesh &mesh, const double spacing)
builds the grid to export the solution
Definition OutData.cpp:2314
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, std::vector< SolutionFrame > &solution_frames) const
saves the volume vtu file
Definition OutData.cpp:1191
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_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:2301
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, std::vector< SolutionFrame > &solution_frames) const
exports everytihng, txt, vtu, etc
Definition OutData.cpp:935
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, std::vector< SolutionFrame > &solution_frames) const
saves the vtu file for time t
Definition OutData.cpp:1108
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, std::vector< SolutionFrame > &solution_frames) const
saves the surface vtu file for for constact quantites, eg contact or friction forces
Definition OutData.cpp:1963
void init_sampler(const polyfem::mesh::Mesh &mesh, const double vismesh_rel_area)
unitalize the ref element sampler
Definition OutData.cpp:2309
Eigen::MatrixXi grid_points_to_elements
grid mesh mapping to fe elements
Definition OutData.hpp:275
utils::RefElementSampler ref_element_sampler
used to sample the solution
Definition OutData.hpp:270
stores all runtime data
Definition OutData.hpp:364
double loading_mesh_time
time to load the mesh
Definition OutData.hpp:369
double total_time()
computes total time
Definition OutData.hpp:383
double assembling_stiffness_mat_time
time to assembly
Definition OutData.hpp:373
double assigning_rhs_time
time to computing the rhs
Definition OutData.hpp:377
double assembling_mass_mat_time
time to assembly mass
Definition OutData.hpp:375
double building_basis_time
time to construct the basis
Definition OutData.hpp:367
double solving_time
time to solve
Definition OutData.hpp:379
double computing_poly_basis_time
time to build the polygonal/polyhedral bases
Definition OutData.hpp:371
all stats from polyfem
Definition OutData.hpp:391
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 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:2777
int n_flipped
number of flipped elements, compute only when using count_flipped_els (false by default)
Definition OutData.hpp:420
int simplex_count
statiscs on the mesh (simplices)
Definition OutData.hpp:423
json solver_info
information of the solver, eg num iteration, time, errors, etc the informations varies depending on t...
Definition OutData.hpp:398
int simple_singular_count
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition OutData.hpp:429
Eigen::Vector4d spectrum
spectrum of the stiffness matrix, enable only if POLYSOLVE_WITH_SPECTRA is ON (off by default)
Definition OutData.hpp:394
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:2499
double l2_err
errors, lp_err is in fact an L8 error
Definition OutData.hpp:408
int undefined_count
statiscs on the mesh (not quad/hex simplex), see Polyspline paper for desciption
Definition OutData.hpp:439
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:2544
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:2412
double average_edge_length
avg edge lenght
Definition OutData.hpp:405
double max_angle
statiscs on angle, compute only when using p_ref (false by default)
Definition OutData.hpp:415
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:435
int non_regular_count
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition OutData.hpp:437
int regular_boundary_count
statiscs on the mesh (regular quad/hex boundary part of the mesh), see Polyspline paper for desciptio...
Definition OutData.hpp:427
long long nn_zero
non zeros and sytem matrix size num dof is the total dof in the system
Definition OutData.hpp:412
double mesh_size
max edge lenght
Definition OutData.hpp:401
void reset()
clears all stats
Definition OutData.cpp:2490
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:441
void compute_mesh_stats(const polyfem::mesh::Mesh &mesh)
compute stats (counts els type, mesh lenght, etc), step 1 of solve
Definition OutData.cpp:2705
double sigma_max
statiscs on tri/tet quality, compute only when using p_ref (false by default)
Definition OutData.hpp:417
int boundary_count
statiscs on the mesh (boundary quads/hexs), see Polyspline paper for desciption
Definition OutData.hpp:433
int regular_count
statiscs on the mesh (regular quad/hex part of the mesh), see Polyspline paper for desciption
Definition OutData.hpp:425
double min_edge_length
min edge lenght
Definition OutData.hpp:403
int multi_singular_count
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition OutData.hpp:431
void write(const int t, const double forward, const double remeshing, const double global_relaxation, const Eigen::MatrixXd &sol)
Definition OutData.cpp:2953
class used to save the solution of time dependent problems in code instead of saving it to the disc
Definition OutData.hpp:30
Eigen::MatrixXi connectivity
Definition OutData.hpp:34
Eigen::MatrixXd scalar_value_avg
Definition OutData.hpp:40
Eigen::MatrixXd solution
Definition OutData.hpp:35
Eigen::MatrixXd exact
Definition OutData.hpp:37
Eigen::MatrixXd scalar_value
Definition OutData.hpp:39
Eigen::MatrixXd error
Definition OutData.hpp:38
Eigen::MatrixXd points
Definition OutData.hpp:33
Eigen::MatrixXd pressure
Definition OutData.hpp:36
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
class to store time stepping data
Definition SolveData.hpp:51
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:89