Loading [MathJax]/jax/input/TeX/config.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
AssemblerUtils.hpp
Go to the documentation of this file.
1#pragma once
2
7
8namespace polyfem::assembler
9{
11 {
12 public:
13 static std::string other_assembler_name(const std::string &formulation);
14
15 static std::shared_ptr<Assembler> make_assembler(const std::string &formulation);
16 static std::shared_ptr<MixedAssembler> make_mixed_assembler(const std::string &formulation);
17
18 enum class BasisType
19 {
22 SPLINE,
23 POLY
24 };
25
26 AssemblerUtils() = delete;
27
31 static void merge_mixed_matrices(
32 const int n_bases, const int n_pressure_bases, const int problem_dim, const bool add_average,
33 const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness,
34 StiffnessMatrix &stiffness);
35
37 static int quadrature_order(const std::string &assembler, const int basis_degree, const BasisType &b_type, const int dim);
38
40 static bool is_elastic_material(const std::string &material);
42 static std::vector<std::string> elastic_materials();
43 };
44
47 {
48
49 public:
51
52 void set_size(const int size);
53 void add_multimaterial(const int index, const json &params, const Units &units);
54 std::shared_ptr<assembler::NLAssembler> get_assembler(const std::string &name) const;
55
56 std::map<std::string, Assembler::ParamFunc> parameters() const;
57
58 private:
59 std::unordered_map<std::string, std::shared_ptr<assembler::NLAssembler>> elastic_material_map_;
60 };
61} // namespace polyfem::assembler
utility to create a map of all elastic materials
void add_multimaterial(const int index, const json &params, const Units &units)
std::unordered_map< std::string, std::shared_ptr< assembler::NLAssembler > > elastic_material_map_
std::map< std::string, Assembler::ParamFunc > parameters() const
std::shared_ptr< assembler::NLAssembler > get_assembler(const std::string &name) const
static std::shared_ptr< MixedAssembler > make_mixed_assembler(const std::string &formulation)
static std::string other_assembler_name(const std::string &formulation)
static int quadrature_order(const std::string &assembler, const int basis_degree, const BasisType &b_type, const int dim)
utility for retrieving the needed quadrature order to precisely integrate the given form on the given...
static void merge_mixed_matrices(const int n_bases, const int n_pressure_bases, const int problem_dim, const bool add_average, const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness, StiffnessMatrix &stiffness)
utility to merge 3 blocks of mixed matrices, A=velocity_stiffness, B=mixed_stiffness,...
static bool is_elastic_material(const std::string &material)
utility to check if material is one of the elastic materials
static std::vector< std::string > elastic_materials()
list of all elastic materials
static std::shared_ptr< Assembler > make_assembler(const std::string &formulation)
Used for test only.
nlohmann::json json
Definition Common.hpp:9
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22