PolyFEM
Loading...
Searching...
No Matches
polyfem::assembler::ElementAssemblyValues Class Reference

stores per element basis values at given quadrature points and geometric mapping More...

#include <ElementAssemblyValues.hpp>

Collaboration diagram for polyfem::assembler::ElementAssemblyValues:
[legend]

Public Member Functions

void compute (const int el_index, const bool is_volume, const Eigen::MatrixXd &pts, const basis::ElementBases &basis, const basis::ElementBases &gbasis)
 computes the per element values at the local (ref el) points (pts) sets basis_values, jac_it, val, and det members
 
void compute (const int el_index, const bool is_volume, const basis::ElementBases &basis, const basis::ElementBases &gbasis)
 computes quadrature points for given element then calls above (overloaded) compute function
 
bool is_geom_mapping_positive (const bool is_volume, const basis::ElementBases &gbasis) const
 check if the element is flipped
 

Public Attributes

std::vector< AssemblyValuesbasis_values
 
std::vector< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > > jac_it
 
quadrature::Quadrature quadrature
 
int element_id
 
Eigen::MatrixXd val
 
Eigen::VectorXd det
 
bool has_parameterization = true
 

Private Member Functions

void finalize_global_element (const Eigen::MatrixXd &v)
 
void finalize2d (const basis::ElementBases &gbasis, const std::vector< AssemblyValues > &gbasis_values)
 void finalize(const Eigen::MatrixXd &v, const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy); void finalize(const Eigen::MatrixXd &v, const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy, const Eigen::MatrixXd &dz);
 
void finalize3d (const basis::ElementBases &gbasis, const std::vector< AssemblyValues > &gbasis_values)
 
bool is_geom_mapping_positive (const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy, const Eigen::MatrixXd &dz) const
 
bool is_geom_mapping_positive (const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy) const
 

Private Attributes

std::vector< AssemblyValuesg_basis_values_cache_
 

Detailed Description

stores per element basis values at given quadrature points and geometric mapping

Definition at line 13 of file ElementAssemblyValues.hpp.

Member Function Documentation

◆ compute() [1/2]

void polyfem::assembler::ElementAssemblyValues::compute ( const int  el_index,
const bool  is_volume,
const basis::ElementBases basis,
const basis::ElementBases gbasis 
)

computes quadrature points for given element then calls above (overloaded) compute function

Definition at line 145 of file ElementAssemblyValues.cpp.

References compute(), polyfem::basis::ElementBases::compute_quadrature(), polyfem::quadrature::Quadrature::points, and quadrature.

Here is the call graph for this function:

◆ compute() [2/2]

void polyfem::assembler::ElementAssemblyValues::compute ( const int  el_index,
const bool  is_volume,
const Eigen::MatrixXd &  pts,
const basis::ElementBases basis,
const basis::ElementBases gbasis 
)

computes the per element values at the local (ref el) points (pts) sets basis_values, jac_it, val, and det members

Definition at line 151 of file ElementAssemblyValues.cpp.

References polyfem::basis::ElementBases::bases, basis_values, element_id, polyfem::basis::ElementBases::evaluate_bases(), polyfem::basis::ElementBases::evaluate_grads(), finalize2d(), finalize3d(), finalize_global_element(), g_basis_values_cache_, polyfem::assembler::AssemblyValues::global, polyfem::basis::Basis::global(), polyfem::assembler::AssemblyValues::grad, polyfem::basis::ElementBases::has_parameterization, polyfem::assembler::AssemblyValues::val, and val.

Referenced by polyfem::solver::OperatorSplittingSolver::advection_FLIP(), polyfem::solver::OperatorSplittingSolver::advection_PIC(), polyfem::assembler::HookeLinearElasticity::assign_stress_tensor(), polyfem::assembler::IncompressibleLinearElasticityDispacement::assign_stress_tensor(), polyfem::assembler::LinearElasticity::assign_stress_tensor(), polyfem::assembler::MooneyRivlin3ParamSymbolic::assign_stress_tensor(), polyfem::assembler::NeoHookeanElasticity::assign_stress_tensor(), polyfem::assembler::SaintVenantElasticity::assign_stress_tensor(), polyfem::State::build_basis(), compute(), polyfem::io::OutStatsData::compute_errors(), polyfem::basis::PolygonalBasis3d::compute_integral_constraints(), polyfem::basis::PolygonalBasis2d::compute_integral_constraints(), polyfem::compute_integral_constraints(), polyfem::solver::AdjointTools::compute_shape_derivative_functional_term(), polyfem::solver::InertiaForm::force_shape_derivative(), polyfem::solver::ElasticForm::force_shape_derivative(), polyfem::io::Evaluator::integrate_function(), polyfem::io::Evaluator::interpolate_at_local_vals(), polyfem::io::Evaluator::interpolate_boundary_function(), polyfem::io::Evaluator::interpolate_boundary_function_at_vertices(), polyfem::io::Evaluator::interpolate_boundary_tensor_function(), polyfem::solver::OperatorSplittingSolver::interpolator(), polyfem::solver::OperatorSplittingSolver::projection(), and polyfem::io::OutGeometryData::save_volume().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalize2d()

void polyfem::assembler::ElementAssemblyValues::finalize2d ( const basis::ElementBases gbasis,
const std::vector< AssemblyValues > &  gbasis_values 
)
private

void finalize(const Eigen::MatrixXd &v, const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy); void finalize(const Eigen::MatrixXd &v, const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy, const Eigen::MatrixXd &dz);

compute Jacobians

Definition at line 104 of file ElementAssemblyValues.cpp.

References polyfem::basis::ElementBases::bases, basis_values, det, polyfem::basis::Basis::global(), polyfem::basis::ElementBases::has_parameterization, jac_it, and val.

Referenced by compute().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalize3d()

void polyfem::assembler::ElementAssemblyValues::finalize3d ( const basis::ElementBases gbasis,
const std::vector< AssemblyValues > &  gbasis_values 
)
private

Definition at line 63 of file ElementAssemblyValues.cpp.

References polyfem::basis::ElementBases::bases, basis_values, det, polyfem::basis::Basis::global(), polyfem::basis::ElementBases::has_parameterization, jac_it, and val.

Referenced by compute().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ finalize_global_element()

void polyfem::assembler::ElementAssemblyValues::finalize_global_element ( const Eigen::MatrixXd &  v)
private

Definition at line 10 of file ElementAssemblyValues.cpp.

References basis_values, det, has_parameterization, jac_it, and val.

Referenced by compute().

Here is the caller graph for this function:

◆ is_geom_mapping_positive() [1/3]

bool polyfem::assembler::ElementAssemblyValues::is_geom_mapping_positive ( const bool  is_volume,
const basis::ElementBases gbasis 
) const

check if the element is flipped

Definition at line 222 of file ElementAssemblyValues.cpp.

References polyfem::basis::ElementBases::bases, polyfem::basis::ElementBases::compute_quadrature(), polyfem::basis::ElementBases::evaluate_grads(), polyfem::basis::Basis::global(), polyfem::basis::ElementBases::has_parameterization, is_geom_mapping_positive(), and polyfem::quadrature::Quadrature::points.

Referenced by is_geom_mapping_positive().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ is_geom_mapping_positive() [2/3]

bool polyfem::assembler::ElementAssemblyValues::is_geom_mapping_positive ( const Eigen::MatrixXd &  dx,
const Eigen::MatrixXd &  dy 
) const
private

Definition at line 45 of file ElementAssemblyValues.cpp.

◆ is_geom_mapping_positive() [3/3]

bool polyfem::assembler::ElementAssemblyValues::is_geom_mapping_positive ( const Eigen::MatrixXd &  dx,
const Eigen::MatrixXd &  dy,
const Eigen::MatrixXd &  dz 
) const
private

Definition at line 26 of file ElementAssemblyValues.cpp.

Member Data Documentation

◆ basis_values

std::vector<AssemblyValues> polyfem::assembler::ElementAssemblyValues::basis_values

Definition at line 20 of file ElementAssemblyValues.hpp.

Referenced by polyfem::solver::OperatorSplittingSolver::advection_FLIP(), polyfem::solver::OperatorSplittingSolver::advection_PIC(), polyfem::assembler::BilaplacianAux::assemble(), polyfem::assembler::Helmholtz::assemble(), polyfem::assembler::HookeLinearElasticity::assemble(), polyfem::assembler::IncompressibleLinearElasticityDispacement::assemble(), polyfem::assembler::IncompressibleLinearElasticityPressure::assemble(), polyfem::assembler::Laplacian::assemble(), polyfem::assembler::LinearElasticity::assemble(), polyfem::assembler::Mass::assemble(), polyfem::assembler::StokesVelocity::assemble(), polyfem::assembler::BilaplacianMixed::assemble(), polyfem::assembler::IncompressibleLinearElasticityMixed::assemble(), polyfem::assembler::StokesMixed::assemble(), polyfem::assembler::GenericElastic< Derived >::assemble_gradient(), polyfem::assembler::HookeLinearElasticity::assemble_gradient(), polyfem::assembler::LinearElasticity::assemble_gradient(), polyfem::assembler::MooneyRivlin3ParamSymbolic::assemble_gradient(), polyfem::assembler::NavierStokesVelocity::assemble_gradient(), polyfem::assembler::NeoHookeanElasticity::assemble_gradient(), polyfem::assembler::SaintVenantElasticity::assemble_gradient(), polyfem::assembler::ViscousDamping::assemble_gradient(), polyfem::assembler::ViscousDampingPrev::assemble_gradient(), polyfem::assembler::GenericElastic< Derived >::assemble_hessian(), polyfem::assembler::HookeLinearElasticity::assemble_hessian(), polyfem::assembler::LinearElasticity::assemble_hessian(), polyfem::assembler::MooneyRivlin3ParamSymbolic::assemble_hessian(), polyfem::assembler::NeoHookeanElasticity::assemble_hessian(), polyfem::assembler::SaintVenantElasticity::assemble_hessian(), polyfem::assembler::ViscousDamping::assemble_hessian(), polyfem::assembler::ViscousDampingPrev::assemble_hessian(), compute(), polyfem::basis::RBFWithQuadratic::compute_constraints_matrix_2d(), polyfem::basis::RBFWithQuadraticLagrange::compute_constraints_matrix_2d(), polyfem::compute_diplacement_grad(), polyfem::compute_diplacement_grad(), polyfem::compute_disp_grad_at_quad(), polyfem::assembler::ViscousDamping::compute_energy(), polyfem::assembler::MooneyRivlin3ParamSymbolic::compute_energy_aux_gradient_fast(), polyfem::assembler::NeoHookeanElasticity::compute_energy_aux_gradient_fast(), polyfem::assembler::MooneyRivlin3ParamSymbolic::compute_energy_hessian_aux_fast(), polyfem::assembler::NeoHookeanElasticity::compute_energy_hessian_aux_fast(), polyfem::compute_integral_constraints(), polyfem::assembler::NavierStokesVelocity::compute_N(), polyfem::assembler::NavierStokesVelocity::compute_W(), finalize2d(), finalize3d(), finalize_global_element(), polyfem::get_local_disp(), and polyfem::basis::RBFWithQuadratic::setup_monomials_vals_2d().

◆ det

◆ element_id

◆ g_basis_values_cache_

std::vector<AssemblyValues> polyfem::assembler::ElementAssemblyValues::g_basis_values_cache_
private

Definition at line 48 of file ElementAssemblyValues.hpp.

Referenced by compute().

◆ has_parameterization

bool polyfem::assembler::ElementAssemblyValues::has_parameterization = true

◆ jac_it

◆ quadrature

◆ val


The documentation for this class was generated from the following files: