Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
AssemblerUtils.cpp
Go to the documentation of this file.
1
2#include "AssemblerUtils.hpp"
3
24
27
28// #include <unsupported/Eigen/SparseExtra>
29
30namespace polyfem
31{
32 using namespace basis;
33 using namespace utils;
34
35 namespace assembler
36 {
37
38 std::string AssemblerUtils::other_assembler_name(const std::string &formulation)
39 {
40 if (formulation == "Bilaplacian")
41 return "BilaplacianAux";
42 else if (formulation == "Stokes" || formulation == "NavierStokes" || formulation == "OperatorSplitting")
43 return "StokesPressure";
44 else if (formulation == "IncompressibleLinearElasticity")
45 return "IncompressibleLinearElasticityPressure";
46
47 return "";
48 }
49
50 std::shared_ptr<Assembler> AssemblerUtils::make_assembler(const std::string &formulation)
51 {
52 if (formulation == "Helmholtz")
53 return std::make_shared<Helmholtz>();
54 else if (formulation == "Laplacian")
55 return std::make_shared<Laplacian>();
56 else if (formulation == "Electrostatics")
57 return std::make_shared<Electrostatics>();
58
59 else if (formulation == "Bilaplacian")
60 return std::make_shared<BilaplacianMain>();
61 else if (formulation == "BilaplacianAux")
62 return std::make_shared<BilaplacianAux>();
63
64 else if (formulation == "LinearElasticity")
65 return std::make_shared<LinearElasticity>();
66 else if (formulation == "HookeLinearElasticity")
67 return std::make_shared<HookeLinearElasticity>();
68 else if (formulation == "IncompressibleLinearElasticity")
69 return std::make_shared<IncompressibleLinearElasticityDispacement>();
70 else if (formulation == "IncompressibleLinearElasticityPressure")
71 return std::make_shared<IncompressibleLinearElasticityPressure>();
72
73 else if (formulation == "SaintVenant")
74 return std::make_shared<SaintVenantElasticity>();
75 else if (formulation == "NeoHookean")
76 return std::make_shared<NeoHookeanElasticity>();
77 else if (formulation == "MooneyRivlin")
78 return std::make_shared<MooneyRivlinElasticity>();
79 else if (formulation == "MooneyRivlin3Param")
80 return std::make_shared<MooneyRivlin3ParamElasticity>();
81 else if (formulation == "MooneyRivlin3ParamSymbolic")
82 return std::make_shared<MooneyRivlin3ParamSymbolic>();
83 else if (formulation == "MultiModels")
84 return std::make_shared<MultiModel>();
85 else if (formulation == "UnconstrainedOgden")
86 return std::make_shared<UnconstrainedOgdenElasticity>();
87 else if (formulation == "IncompressibleOgden")
88 return std::make_shared<IncompressibleOgdenElasticity>();
89
90 else if (formulation == "Stokes")
91 return std::make_shared<StokesVelocity>();
92 else if (formulation == "StokesPressure")
93 return std::make_shared<StokesPressure>();
94 else if (formulation == "NavierStokes")
95 return std::make_shared<NavierStokesVelocity>();
96 else if (formulation == "OperatorSplitting")
97 return std::make_shared<OperatorSplitting>();
98
99 else if (formulation == "AMIPS")
100 return std::make_shared<AMIPSEnergy>();
101 else if (formulation == "AMIPSAutodiff")
102 return std::make_shared<AMIPSEnergyAutodiff>();
103 else if (formulation == "FixedCorotational")
104 return std::make_shared<FixedCorotational>();
105
106 log_and_throw_error("Inavalid assembler name {}", formulation);
107 }
108
109 std::shared_ptr<MixedAssembler> AssemblerUtils::make_mixed_assembler(const std::string &formulation)
110 {
111 if (formulation == "Bilaplacian")
112 return std::make_shared<BilaplacianMixed>();
113 else if (formulation == "IncompressibleLinearElasticity")
114 return std::make_shared<IncompressibleLinearElasticityMixed>();
115 else if (formulation == "Stokes" || formulation == "NavierStokes" || formulation == "OperatorSplitting")
116 return std::make_shared<StokesMixed>();
117
118 log_and_throw_error("Inavalid mixed assembler name {}", formulation);
119 }
120
122 const int n_bases, const int n_pressure_bases, const int problem_dim, const bool add_average,
123 const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness,
124 StiffnessMatrix &stiffness)
125 {
126 assert(velocity_stiffness.rows() == velocity_stiffness.cols());
127 assert(velocity_stiffness.rows() == n_bases * problem_dim);
128
129 assert(mixed_stiffness.size() == 0 || mixed_stiffness.rows() == n_bases * problem_dim);
130 assert(mixed_stiffness.size() == 0 || mixed_stiffness.cols() == n_pressure_bases);
131
132 assert(pressure_stiffness.size() == 0 || pressure_stiffness.rows() == n_pressure_bases);
133 assert(pressure_stiffness.size() == 0 || pressure_stiffness.cols() == n_pressure_bases);
134
135 const int avg_offset = add_average ? 1 : 0;
136
137 std::vector<Eigen::Triplet<double>> blocks;
138 blocks.reserve(velocity_stiffness.nonZeros() + 2 * mixed_stiffness.nonZeros() + pressure_stiffness.nonZeros() + 2 * avg_offset * velocity_stiffness.rows());
139
140 for (int k = 0; k < velocity_stiffness.outerSize(); ++k)
141 {
142 for (StiffnessMatrix::InnerIterator it(velocity_stiffness, k); it; ++it)
143 {
144 blocks.emplace_back(it.row(), it.col(), it.value());
145 }
146 }
147
148 for (int k = 0; k < mixed_stiffness.outerSize(); ++k)
149 {
150 for (StiffnessMatrix::InnerIterator it(mixed_stiffness, k); it; ++it)
151 {
152 blocks.emplace_back(it.row(), n_bases * problem_dim + it.col(), it.value());
153 blocks.emplace_back(it.col() + n_bases * problem_dim, it.row(), it.value());
154 }
155 }
156
157 for (int k = 0; k < pressure_stiffness.outerSize(); ++k)
158 {
159 for (StiffnessMatrix::InnerIterator it(pressure_stiffness, k); it; ++it)
160 {
161 blocks.emplace_back(n_bases * problem_dim + it.row(), n_bases * problem_dim + it.col(), it.value());
162 }
163 }
164
165 if (add_average)
166 {
167 const double val = 1.0 / n_pressure_bases;
168 for (int i = 0; i < n_pressure_bases; ++i)
169 {
170 blocks.emplace_back(n_bases * problem_dim + i, n_bases * problem_dim + n_pressure_bases, val);
171 blocks.emplace_back(n_bases * problem_dim + n_pressure_bases, n_bases * problem_dim + i, val);
172 }
173 }
174
175 stiffness.resize(n_bases * problem_dim + n_pressure_bases + avg_offset, n_bases * problem_dim + n_pressure_bases + avg_offset);
176 stiffness.setFromTriplets(blocks.begin(), blocks.end());
177 stiffness.makeCompressed();
178
179 // static int c = 0;
180 // Eigen::saveMarket(stiffness, "stiffness.txt");
181 // Eigen::saveMarket(velocity_stiffness, "velocity_stiffness.txt");
182 // Eigen::saveMarket(mixed_stiffness, "mixed_stiffness.txt");
183 // Eigen::saveMarket(pressure_stiffness, "pressure_stiffness.txt");
184 }
185
186 int AssemblerUtils::quadrature_order(const std::string &assembler, const int basis_degree, const BasisType &b_type, const int dim)
187 {
188 // note: minimum quadrature order is always 1
189 if (assembler == "Mass")
190 {
191 // multiply by two since we are multiplying phi_i by phi_j
192 if (b_type == BasisType::SIMPLEX_LAGRANGE || b_type == BasisType::CUBE_LAGRANGE)
193 return std::max(basis_degree * 2, 1);
194 else
195 return basis_degree * 2 + 1;
196 }
197 else if (assembler == "NavierStokes")
198 {
199 if (b_type == BasisType::SIMPLEX_LAGRANGE)
200 return std::max((basis_degree - 1) + basis_degree, 1);
201 else if (b_type == BasisType::CUBE_LAGRANGE)
202 return std::max(basis_degree * 2, 1);
203 else
204 return basis_degree * 2 + 1;
205 }
206 else
207 {
208 // subtract one since we take a derivative (lowers polynomial order by 1)
209 // multiply by two since we are multiplying grad phi_i by grad phi_j
210 if (b_type == BasisType::SIMPLEX_LAGRANGE)
211 {
212 return std::max((basis_degree - 1) * 2, 1);
213 }
214 else if (b_type == BasisType::CUBE_LAGRANGE)
215 {
216 // in this case we have a tensor product basis
217 // this computes the quadrature order along a single axis
218 // the Quadrature itself takes a tensor product of the given quadrature points
219 // to form the full quadrature for the basis
220 // taking a gradient leaves at least one variable whose power remains unchanged
221 // thus, we don't subtract 1
222 // note that this is overkill for the variable that was differentiated
223 return std::max(basis_degree * 2, 1);
224 }
225 else
226 {
227 return (basis_degree - 1) * 2 + 1;
228 }
229 }
230 }
231
232 std::vector<std::string> AssemblerUtils::elastic_materials()
233 {
234 const static std::vector<std::string> elastic_materials = {
235 "LinearElasticity",
236 "HookeLinearElasticity",
237 "SaintVenant",
238 "NeoHookean",
239 "MooneyRivlin",
240 "MooneyRivlin3Param",
241 "UnconstrainedOgden",
242 "IncompressibleOgden",
243 "FixedCorotational",
244 "MultiModels"};
245
246 return elastic_materials;
247 }
248
249 bool AssemblerUtils::is_elastic_material(const std::string &material)
250 {
251 for (const auto &m : elastic_materials())
252 {
253 if (material == m)
254 return true;
255 }
256 return false;
257 }
258
260 {
261 for (const auto &m : AssemblerUtils::elastic_materials())
262 {
263 // skip multimodels
264 // this is a special case where we have multiple models
265 // and we need to create a new assembler for each model
266 // this is handled in the MultiModel class
267 // and not here
268 if (m == "MultiModels")
269 continue;
270 const auto assembler = AssemblerUtils::make_assembler(m);
271 // cast assembler to elasticity assembler
272 elastic_material_map_[m] = std::dynamic_pointer_cast<NLAssembler>(assembler);
273 assert(elastic_material_map_[m] != nullptr);
274 }
275 }
276
277 void AllElasticMaterials::set_size(const int size)
278 {
279 for (auto &it : elastic_material_map_)
280 {
281 it.second->set_size(size);
282 }
283 }
284
285 void AllElasticMaterials::add_multimaterial(const int index, const json &params, const Units &units)
286 {
287 for (auto &it : elastic_material_map_)
288 {
289 it.second->add_multimaterial(index, params, units);
290 }
291 }
292
293 std::shared_ptr<assembler::NLAssembler> AllElasticMaterials::get_assembler(const std::string &name) const
294 {
295 return elastic_material_map_.at(name);
296 }
297
298 std::map<std::string, Assembler::ParamFunc> AllElasticMaterials::parameters() const
299 {
300 std::map<std::string, Assembler::ParamFunc> params;
301 for (const auto &m : elastic_material_map_)
302 {
303 const auto assembler = m.second;
304 auto p = assembler->parameters();
305 for (auto &it : p)
306 {
307 params[m.first + "/" + it.first] = it.second;
308 }
309 }
310 return params;
311 }
312 } // namespace assembler
313} // namespace polyfem
double val
Definition Assembler.cpp:86
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)
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22