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