PolyFEM
Loading...
Searching...
No Matches
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 {
23 SPLINE,
24 POLY
25 };
26
27 AssemblerUtils() = delete;
28
32 static void merge_mixed_matrices(
33 const int n_bases, const int n_pressure_bases, const int problem_dim, const bool add_average,
34 const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness,
35 StiffnessMatrix &stiffness);
36
38 static int quadrature_order(const std::string &assembler, const int basis_degree, const BasisType &b_type, const int dim);
39
41 static bool is_elastic_material(const std::string &material);
43 static std::vector<std::string> elastic_materials();
44 };
45
48 {
49
50 public:
52
53 void set_size(const int size);
54 void add_multimaterial(const int index, const json &params, const Units &units);
55 std::shared_ptr<assembler::ElasticityNLAssembler> get_assembler(const std::string &name) const;
56
57 std::map<std::string, Assembler::ParamFunc> parameters() const;
58
59 private:
60 std::unordered_map<std::string, std::shared_ptr<assembler::ElasticityNLAssembler>> elastic_material_map_;
61 };
62} // 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::shared_ptr< assembler::ElasticityNLAssembler > get_assembler(const std::string &name) const
std::map< std::string, Assembler::ParamFunc > parameters() const
std::unordered_map< std::string, std::shared_ptr< assembler::ElasticityNLAssembler > > elastic_material_map_
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:24