|
PolyFEM
|
all stats from polyfem More...
#include <OutData.hpp>
Public Member Functions | |
| 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 | |
| void | compute_mesh_stats (const polyfem::mesh::Mesh &mesh) |
| compute stats (counts els type, mesh lenght, etc), step 1 of solve | |
| 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 default) to compute the size of the linear mesh | |
| void | reset () |
| clears all stats | |
| void | count_flipped_elements (const polyfem::mesh::Mesh &mesh, const std::vector< polyfem::basis::ElementBases > &gbases) |
| counts the number of flipped elements | |
| 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 | |
Public Attributes | |
| Eigen::Vector4d | spectrum |
| spectrum of the stiffness matrix, enable only if POLYSOLVE_WITH_SPECTRA is ON (off by default) | |
| json | solver_info |
| information of the solver, eg num iteration, time, errors, etc the informations varies depending on the solver | |
| double | mesh_size |
| max edge lenght | |
| double | min_edge_length |
| min edge lenght | |
| double | average_edge_length |
| avg edge lenght | |
| double | l2_err |
| errors, lp_err is in fact an L8 error | |
| double | linf_err |
| double | lp_err |
| double | h1_err |
| double | h1_semi_err |
| double | grad_max_err |
| long long | nn_zero |
| non zeros and sytem matrix size num dof is the total dof in the system | |
| long long | mat_size |
| long long | num_dofs |
| double | max_angle |
| statiscs on angle, compute only when using p_ref (false by default) | |
| double | sigma_max |
| statiscs on tri/tet quality, compute only when using p_ref (false by default) | |
| double | sigma_min |
| double | sigma_avg |
| int | n_flipped |
| number of flipped elements, compute only when using count_flipped_els (false by default) | |
| int | simplex_count |
| statiscs on the mesh (simplices) | |
| int | prism_count |
| statiscs on the mesh (simplices) | |
| int | regular_count |
| statiscs on the mesh (regular quad/hex part of the mesh), see Polyspline paper for desciption | |
| int | regular_boundary_count |
| statiscs on the mesh (regular quad/hex boundary part of the mesh), see Polyspline paper for desciption | |
| int | simple_singular_count |
| statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption | |
| int | multi_singular_count |
| statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption | |
| int | boundary_count |
| statiscs on the mesh (boundary quads/hexs), see Polyspline paper for desciption | |
| int | non_regular_boundary_count |
| statiscs on the mesh (irregular boundary quad/hex part of the mesh), see Polyspline paper for desciption | |
| int | non_regular_count |
| statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption | |
| int | undefined_count |
| statiscs on the mesh (not quad/hex simplex), see Polyspline paper for desciption | |
| int | multi_singular_boundary_count |
| statiscs on the mesh (irregular boundary quad/hex part of the mesh), see Polyspline paper for desciption | |
all stats from polyfem
Definition at line 374 of file OutData.hpp.
| void polyfem::io::OutStatsData::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
| [in] | n_bases | number of base |
| [in] | bases | bases |
| [in] | gbases | geometric bases |
| [in] | mesh | mesh |
| [in] | problem | problem |
| [in] | tend | end time step |
| [in] | sol | solution |
Definition at line 2805 of file OutData.cpp.
References polyfem::assembler::ElementAssemblyValues::compute(), polyfem::mesh::Mesh::dimension(), polyfem::assembler::Problem::exact(), polyfem::assembler::Problem::exact_grad(), polyfem::assembler::Problem::has_exact_sol(), polyfem::assembler::Problem::is_scalar(), polyfem::mesh::Mesh::is_volume(), polyfem::logger(), val, and vals.
Referenced by polyfem::State::compute_errors().
| void polyfem::io::OutStatsData::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 default) to compute the size of the linear mesh
| [in] | mesh | to compute stats |
| [in] | bases | geom bases |
| [in] | n_samples | used for curved meshes |
| [in] | use_curved_mesh_size | use curved edges to compute mesh size |
Definition at line 2673 of file OutData.cpp.
References polyfem::mesh::Mesh::get_edges(), polyfem::mesh::Mesh::is_polytope(), polyfem::mesh::Mesh::is_simplex(), polyfem::mesh::Mesh::is_volume(), polyfem::logger(), polyfem::utils::EdgeSampler::sample_2d_cube(), polyfem::utils::EdgeSampler::sample_2d_simplex(), polyfem::utils::EdgeSampler::sample_3d_cube(), and polyfem::utils::EdgeSampler::sample_3d_simplex().
Referenced by polyfem::State::build_basis().
| void polyfem::io::OutStatsData::compute_mesh_stats | ( | const polyfem::mesh::Mesh & | mesh | ) |
compute stats (counts els type, mesh lenght, etc), step 1 of solve
| mesh | mesh |
Definition at line 2966 of file OutData.cpp.
References polyfem::mesh::Mesh::elements_tag(), polyfem::logger(), and polyfem::mesh::Mesh::n_elements().
Referenced by forward_simulation(), and polyfem::State::solve().
| void polyfem::io::OutStatsData::count_flipped_elements | ( | const polyfem::mesh::Mesh & | mesh, |
| const std::vector< polyfem::basis::ElementBases > & | gbases | ||
| ) |
counts the number of flipped elements
| [in] | mesh | mesh |
| [in] | gbases | geometric bases |
Definition at line 2760 of file OutData.cpp.
References polyfem::mesh::Mesh::elements_tag(), polyfem::mesh::Mesh::is_polytope(), polyfem::mesh::Mesh::is_volume(), polyfem::log_and_throw_error(), polyfem::logger(), and vals.
Referenced by polyfem::State::build_basis().
| void polyfem::io::OutStatsData::reset | ( | ) |
clears all stats
Definition at line 2751 of file OutData.cpp.
Referenced by polyfem::State::build_basis().
| void polyfem::io::OutStatsData::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
| [in] | j | output json |
save json
| [in] | args | input argumeents |
| [in] | n_bases | number of bases |
| [in] | n_pressure_bases | number fo pressure bases |
| [in] | sol | solution |
| [in] | mesh | mesh |
| [in] | disc_orders | discretization order |
| [in] | disc_ordersq | discretization order |
| [in] | problem | problem |
| [in] | runtime | rumtime |
| [in] | formulation | formulation |
| [in] | isoparametric | if isoparametric |
| [in] | sol_at_node_id | export solution at node |
| [out] | j | output json |
Definition at line 3043 of file OutData.cpp.
References polyfem::io::OutRuntimeData::assembling_mass_mat_time, polyfem::io::OutRuntimeData::assembling_stiffness_mat_time, polyfem::io::OutRuntimeData::assigning_rhs_time, polyfem::io::OutRuntimeData::building_basis_time, polyfem::io::OutRuntimeData::computing_poly_basis_time, polyfem::mesh::Mesh::dimension(), polyfem::utils::get_n_threads(), getPeakRSS(), polyfem::assembler::Problem::is_scalar(), polyfem::io::OutRuntimeData::loading_mesh_time, polyfem::logger(), polyfem::mesh::Mesh::n_elements(), polyfem::mesh::Mesh::n_vertices(), polyfem::assembler::Problem::name(), polyfem::mesh::Mesh::orders(), and polyfem::io::OutRuntimeData::solving_time.
Referenced by polyfem::State::save_json().
| double polyfem::io::OutStatsData::average_edge_length |
avg edge lenght
Definition at line 389 of file OutData.hpp.
| int polyfem::io::OutStatsData::boundary_count |
statiscs on the mesh (boundary quads/hexs), see Polyspline paper for desciption
Definition at line 419 of file OutData.hpp.
| double polyfem::io::OutStatsData::grad_max_err |
Definition at line 392 of file OutData.hpp.
| double polyfem::io::OutStatsData::h1_err |
Definition at line 392 of file OutData.hpp.
| double polyfem::io::OutStatsData::h1_semi_err |
Definition at line 392 of file OutData.hpp.
| double polyfem::io::OutStatsData::l2_err |
errors, lp_err is in fact an L8 error
Definition at line 392 of file OutData.hpp.
| double polyfem::io::OutStatsData::linf_err |
Definition at line 392 of file OutData.hpp.
| double polyfem::io::OutStatsData::lp_err |
Definition at line 392 of file OutData.hpp.
| long long polyfem::io::OutStatsData::mat_size |
Definition at line 396 of file OutData.hpp.
Referenced by polyfem::State::assemble_mass_mat(), and polyfem::State::build_stiffness_mat().
| double polyfem::io::OutStatsData::max_angle |
statiscs on angle, compute only when using p_ref (false by default)
Definition at line 399 of file OutData.hpp.
Referenced by polyfem::refinement::APriori::p_refine(), and polyfem::refinement::APriori::p_refine().
| double polyfem::io::OutStatsData::mesh_size |
max edge lenght
Definition at line 385 of file OutData.hpp.
Referenced by polyfem::State::build_basis().
| double polyfem::io::OutStatsData::min_edge_length |
min edge lenght
Definition at line 387 of file OutData.hpp.
Referenced by polyfem::State::build_basis().
| int polyfem::io::OutStatsData::multi_singular_boundary_count |
statiscs on the mesh (irregular boundary quad/hex part of the mesh), see Polyspline paper for desciption
Definition at line 427 of file OutData.hpp.
| int polyfem::io::OutStatsData::multi_singular_count |
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition at line 417 of file OutData.hpp.
| int polyfem::io::OutStatsData::n_flipped |
number of flipped elements, compute only when using count_flipped_els (false by default)
Definition at line 404 of file OutData.hpp.
Referenced by polyfem::State::build_basis().
| long long polyfem::io::OutStatsData::nn_zero |
non zeros and sytem matrix size num dof is the total dof in the system
Definition at line 396 of file OutData.hpp.
Referenced by polyfem::State::assemble_mass_mat(), and polyfem::State::build_stiffness_mat().
| int polyfem::io::OutStatsData::non_regular_boundary_count |
statiscs on the mesh (irregular boundary quad/hex part of the mesh), see Polyspline paper for desciption
Definition at line 421 of file OutData.hpp.
| int polyfem::io::OutStatsData::non_regular_count |
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition at line 423 of file OutData.hpp.
| long long polyfem::io::OutStatsData::num_dofs |
Definition at line 396 of file OutData.hpp.
Referenced by polyfem::State::assemble_mass_mat(), and polyfem::State::build_stiffness_mat().
| int polyfem::io::OutStatsData::prism_count |
statiscs on the mesh (simplices)
Definition at line 409 of file OutData.hpp.
| int polyfem::io::OutStatsData::regular_boundary_count |
statiscs on the mesh (regular quad/hex boundary part of the mesh), see Polyspline paper for desciption
Definition at line 413 of file OutData.hpp.
| int polyfem::io::OutStatsData::regular_count |
statiscs on the mesh (regular quad/hex part of the mesh), see Polyspline paper for desciption
Definition at line 411 of file OutData.hpp.
| double polyfem::io::OutStatsData::sigma_avg |
Definition at line 401 of file OutData.hpp.
Referenced by polyfem::refinement::APriori::p_refine(), and polyfem::refinement::APriori::p_refine().
| double polyfem::io::OutStatsData::sigma_max |
statiscs on tri/tet quality, compute only when using p_ref (false by default)
Definition at line 401 of file OutData.hpp.
Referenced by polyfem::refinement::APriori::p_refine(), and polyfem::refinement::APriori::p_refine().
| double polyfem::io::OutStatsData::sigma_min |
Definition at line 401 of file OutData.hpp.
Referenced by polyfem::refinement::APriori::p_refine(), and polyfem::refinement::APriori::p_refine().
| int polyfem::io::OutStatsData::simple_singular_count |
statiscs on the mesh (irregular quad/hex part of the mesh), see Polyspline paper for desciption
Definition at line 415 of file OutData.hpp.
| int polyfem::io::OutStatsData::simplex_count |
statiscs on the mesh (simplices)
Definition at line 407 of file OutData.hpp.
| json polyfem::io::OutStatsData::solver_info |
information of the solver, eg num iteration, time, errors, etc the informations varies depending on the solver
Definition at line 382 of file OutData.hpp.
Referenced by polyfem::State::init_nonlinear_tensor_solve(), polyfem::State::solve_linear(), and polyfem::State::solve_tensor_nonlinear().
| Eigen::Vector4d polyfem::io::OutStatsData::spectrum |
spectrum of the stiffness matrix, enable only if POLYSOLVE_WITH_SPECTRA is ON (off by default)
Definition at line 378 of file OutData.hpp.
Referenced by polyfem::State::solve_linear(), and polyfem::State::solve_problem().
| int polyfem::io::OutStatsData::undefined_count |
statiscs on the mesh (not quad/hex simplex), see Polyspline paper for desciption
Definition at line 425 of file OutData.hpp.