LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Postprocessing such as computing norms and Mesh Functions and plot scripts. More...
Classes | |
class | MeshFunctionWhitneyOne |
Provides Mesh Function for given basis expansion coefficients. More... | |
class | MeshFunctionWhitneyTwo |
Provides Mesh Function for Whitney two basis expansion coefficients. More... | |
class | MeshFunctionWhitneyZero |
Provides Mesh Function for basis expansion coefficients of the zero form. More... | |
struct | ProblemSolution |
stores solutions for several k and one refinement level More... | |
struct | ProblemSolutionWrapper |
stores solutions and setupinformation for a list of refinement levels and k values More... | |
Functions | |
template<typename MF , typename SQ_F > | |
std::pair< double, lf::mesh::utils::CodimMeshDataSet< double > > | L2norm (const std::shared_ptr< const lf::mesh::Mesh > &mesh_p, const MF &f, const SQ_F &sq_f, const lf::quad::QuadRule &quadrule) |
Computes the L2Norm of some MeshFunction (type MF) f. More... | |
template<typename MF , typename SQ_F > | |
std::pair< double, lf::mesh::utils::CodimMeshDataSet< double > > | SupNorm (const std::shared_ptr< const lf::mesh::Mesh > &mesh_p, const MF &f, const SQ_F &sq_f, const lf::quad::QuadRule &quadrule) |
Computes the supremum norm of the squared values of some MeshFunction (type MF) f. More... | |
template<typename... Args> | |
static std::string | concat (Args &&...args) |
Concatenate objects defining an operator<<(std::ostream&) More... | |
template<typename SCALAR > | |
void | compare_results (std::string name, ProblemSolutionWrapper< SCALAR > &results_lap, ProblemSolutionWrapper< SCALAR > &results_dirac) |
Generate vtk files containing meshdata and a csv table with results of the comparison of the two given basis expansion coefficients. More... | |
template<typename U_ZERO , typename U_ONE , typename U_TWO , typename SCALAR > | |
void | process_results (std::string name, ProblemSolutionWrapper< SCALAR > &results, U_ZERO &u_zero, U_ONE &u_one, U_TWO &u_two, double &k) |
Generate vtk files containing meshdata and a csv table with results. More... | |
Postprocessing such as computing norms and Mesh Functions and plot scripts.
void projects::hldo_sphere::post_processing::compare_results | ( | std::string | name, |
ProblemSolutionWrapper< SCALAR > & | results_lap, | ||
ProblemSolutionWrapper< SCALAR > & | results_dirac | ||
) |
Generate vtk files containing meshdata and a csv table with results of the comparison of the two given basis expansion coefficients.
The Table contains for each refinement level, the number of cells, edges and vertices, as well as the global mesh width (the maximum edge length of a triangle). Moreover it contains the l2 norm of the difference of the two functions as well as the approximate orders for every value of k and all three forms.
SCALAR | base type of the result vectors |
name | identifying the experiment (used to name the outputs) |
results_lap | ProblemSolutionWrapper with the results for all levels and ks computed with the Laplace operator |
results_dirac | ProblemSolutionWrapper with the results for all levels and ks computed with the Dirac operator |
Definition at line 96 of file results_processing.h.
References concat(), projects::hldo_sphere::post_processing::ProblemSolutionWrapper< SCALAR >::k, lf::base::RefEl::kTria(), L2norm(), projects::hldo_sphere::post_processing::ProblemSolutionWrapper< SCALAR >::levels, lf::quad::make_QuadRule(), projects::hldo_sphere::post_processing::ProblemSolutionWrapper< SCALAR >::mesh, projects::hldo_sphere::post_processing::ProblemSolutionWrapper< SCALAR >::solutions, lf::geometry::Volume(), and lf::io::VtkWriter::WriteCellData().
|
static |
Concatenate objects defining an operator<<(std::ostream&)
args | A variadic pack of objects implementing operator<<(std::ostream&) |
Definition at line 68 of file results_processing.h.
Referenced by compare_results(), and process_results().
std::pair< double, lf::mesh::utils::CodimMeshDataSet< double > > projects::hldo_sphere::post_processing::L2norm | ( | const std::shared_ptr< const lf::mesh::Mesh > & | mesh_p, |
const MF & | f, | ||
const SQ_F & | sq_f, | ||
const lf::quad::QuadRule & | quadrule | ||
) |
Computes the L2Norm of some MeshFunction (type MF) f.
MF | type of the mesh function to take to norm form is required to be defined on the mesh |
SQ_F | type of the square function sq_f |
mesh_p | pointer to the mesh |
f | meshfunction to evaluate |
sq_f | square function for the output of f \(sq\_f(x) \mapsto \|x\|^2\) |
quadrule | quadrature rule used for each cell in the mesh |
Iterates over all cells in the mesh and compute the squared integral locally, then add up and square
Definition at line 33 of file norms.h.
References lf::quad::QuadRule::NumPoints(), lf::quad::QuadRule::Points(), and lf::quad::QuadRule::Weights().
Referenced by compare_results(), projects::hldo_sphere::debugging::DiracConvergenceTest::Compute(), projects::hldo_sphere::debugging::WhitneyOneBasisExpansionCoeffs::GetL2Error(), and process_results().
void projects::hldo_sphere::post_processing::process_results | ( | std::string | name, |
ProblemSolutionWrapper< SCALAR > & | results, | ||
U_ZERO & | u_zero, | ||
U_ONE & | u_one, | ||
U_TWO & | u_two, | ||
double & | k | ||
) |
Generate vtk files containing meshdata and a csv table with results.
The Table contains for each refinement level, the number of cells, edges and vetices, as well as the global mesh width (the maximum edgelenght of a triangle). Moreover it contains the l2 error and the Supremum error as well as the approximate orders for every value of k and all three forms.
U_ZERO | functor for the analytical solution of the zero form |
U_ONE | functor for the analytical solution of the one form |
U_TWO | functor for the analytical solution of the two form |
SCALAR | base type of the resultvectors |
name | identifying the experiment (used to name the outputs) |
results | ProblemSolutionWrapper with the results for all levels and ks |
u_zero | vector with the analytical solutions for u_zero of every k |
u_one | vector with the analytical solutions for u_one of every k |
u_two | vector with the analytical solutions for u_two of every k |
k | reference to the variable used in the functions u_one, u_two, u_zero |
Definition at line 302 of file results_processing.h.
References concat(), projects::hldo_sphere::post_processing::ProblemSolutionWrapper< SCALAR >::k, lf::base::RefEl::kTria(), L2norm(), projects::hldo_sphere::post_processing::ProblemSolutionWrapper< SCALAR >::levels, lf::quad::make_QuadRule(), projects::hldo_sphere::post_processing::ProblemSolutionWrapper< SCALAR >::mesh, projects::hldo_sphere::post_processing::ProblemSolutionWrapper< SCALAR >::solutions, SupNorm(), lf::geometry::Volume(), and lf::io::VtkWriter::WriteCellData().
Referenced by projects::hldo_sphere::experiments::DiracOperatorExperiment::Compute(), and projects::hldo_sphere::experiments::HodgeLaplacianExperiment< SCALAR >::Compute().
std::pair< double, lf::mesh::utils::CodimMeshDataSet< double > > projects::hldo_sphere::post_processing::SupNorm | ( | const std::shared_ptr< const lf::mesh::Mesh > & | mesh_p, |
const MF & | f, | ||
const SQ_F & | sq_f, | ||
const lf::quad::QuadRule & | quadrule | ||
) |
Computes the supremum norm of the squared values of some MeshFunction (type MF) f.
MF | type of the mesh function to take to norm form |
SQ_F | type of the square function sq_f |
mesh_p | pointer to the mesh |
f | meshfunction to evaluate |
sq_f | square function for the output of f \( sq\_f(x) \mapsto \|x\|^2 \) |
quadrule | quadrature rule only used to determine the evalueation points |
Definition at line 84 of file norms.h.
References lf::quad::QuadRule::NumPoints(), and lf::quad::QuadRule::Points().
Referenced by process_results().