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