PolyFEM
Loading...
Searching...
No Matches
AssemblerUtils.cpp
Go to the documentation of this file.
1
2#include "AssemblerUtils.hpp"
3
22
25
26// #include <unsupported/Eigen/SparseExtra>
27
28namespace polyfem
29{
30 using namespace basis;
31 using namespace utils;
32
33 namespace assembler
34 {
35
36 std::string AssemblerUtils::other_assembler_name(const std::string &formulation)
37 {
38 if (formulation == "Bilaplacian")
39 return "BilaplacianAux";
40 else if (formulation == "Stokes" || formulation == "NavierStokes" || formulation == "OperatorSplitting")
41 return "StokesPressure";
42 else if (formulation == "IncompressibleLinearElasticity")
43 return "IncompressibleLinearElasticityPressure";
44
45 return "";
46 }
47
48 std::shared_ptr<Assembler> AssemblerUtils::make_assembler(const std::string &formulation)
49 {
50 if (formulation == "Helmholtz")
51 return std::make_shared<Helmholtz>();
52 else if (formulation == "Laplacian")
53 return std::make_shared<Laplacian>();
54 else if (formulation == "Bilaplacian")
55 return std::make_shared<BilaplacianMain>();
56 else if (formulation == "BilaplacianAux")
57 return std::make_shared<BilaplacianAux>();
58
59 else if (formulation == "LinearElasticity")
60 return std::make_shared<LinearElasticity>();
61 else if (formulation == "HookeLinearElasticity")
62 return std::make_shared<HookeLinearElasticity>();
63 else if (formulation == "IncompressibleLinearElasticity")
64 return std::make_shared<IncompressibleLinearElasticityDispacement>();
65 else if (formulation == "IncompressibleLinearElasticityPressure")
66 return std::make_shared<IncompressibleLinearElasticityPressure>();
67
68 else if (formulation == "SaintVenant")
69 return std::make_shared<SaintVenantElasticity>();
70 else if (formulation == "NeoHookean")
71 return std::make_shared<NeoHookeanElasticity>();
72 else if (formulation == "MooneyRivlin")
73 return std::make_shared<MooneyRivlinElasticity>();
74 else if (formulation == "MooneyRivlin3Param")
75 return std::make_shared<MooneyRivlin3ParamElasticity>();
76 else if (formulation == "MooneyRivlin3ParamSymbolic")
77 return std::make_shared<MooneyRivlin3ParamSymbolic>();
78 else if (formulation == "MultiModels")
79 return std::make_shared<MultiModel>();
80 else if (formulation == "UnconstrainedOgden")
81 return std::make_shared<UnconstrainedOgdenElasticity>();
82 else if (formulation == "IncompressibleOgden")
83 return std::make_shared<IncompressibleOgdenElasticity>();
84
85 else if (formulation == "Stokes")
86 return std::make_shared<StokesVelocity>();
87 else if (formulation == "StokesPressure")
88 return std::make_shared<StokesPressure>();
89 else if (formulation == "NavierStokes")
90 return std::make_shared<NavierStokesVelocity>();
91 else if (formulation == "OperatorSplitting")
92 return std::make_shared<OperatorSplitting>();
93
94 else if (formulation == "AMIPS")
95 return std::make_shared<AMIPSEnergy>();
96
97 log_and_throw_error("Inavalid assembler name {}", formulation);
98 }
99
100 std::shared_ptr<MixedAssembler> AssemblerUtils::make_mixed_assembler(const std::string &formulation)
101 {
102 if (formulation == "Bilaplacian")
103 return std::make_shared<BilaplacianMixed>();
104 else if (formulation == "IncompressibleLinearElasticity")
105 return std::make_shared<IncompressibleLinearElasticityMixed>();
106 else if (formulation == "Stokes" || formulation == "NavierStokes" || formulation == "OperatorSplitting")
107 return std::make_shared<StokesMixed>();
108
109 log_and_throw_error("Inavalid mixed assembler name {}", formulation);
110 }
111
113 const int n_bases, const int n_pressure_bases, const int problem_dim, const bool add_average,
114 const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness,
115 StiffnessMatrix &stiffness)
116 {
117 assert(velocity_stiffness.rows() == velocity_stiffness.cols());
118 assert(velocity_stiffness.rows() == n_bases * problem_dim);
119
120 assert(mixed_stiffness.size() == 0 || mixed_stiffness.rows() == n_bases * problem_dim);
121 assert(mixed_stiffness.size() == 0 || mixed_stiffness.cols() == n_pressure_bases);
122
123 assert(pressure_stiffness.size() == 0 || pressure_stiffness.rows() == n_pressure_bases);
124 assert(pressure_stiffness.size() == 0 || pressure_stiffness.cols() == n_pressure_bases);
125
126 const int avg_offset = add_average ? 1 : 0;
127
128 std::vector<Eigen::Triplet<double>> blocks;
129 blocks.reserve(velocity_stiffness.nonZeros() + 2 * mixed_stiffness.nonZeros() + pressure_stiffness.nonZeros() + 2 * avg_offset * velocity_stiffness.rows());
130
131 for (int k = 0; k < velocity_stiffness.outerSize(); ++k)
132 {
133 for (StiffnessMatrix::InnerIterator it(velocity_stiffness, k); it; ++it)
134 {
135 blocks.emplace_back(it.row(), it.col(), it.value());
136 }
137 }
138
139 for (int k = 0; k < mixed_stiffness.outerSize(); ++k)
140 {
141 for (StiffnessMatrix::InnerIterator it(mixed_stiffness, k); it; ++it)
142 {
143 blocks.emplace_back(it.row(), n_bases * problem_dim + it.col(), it.value());
144 blocks.emplace_back(it.col() + n_bases * problem_dim, it.row(), it.value());
145 }
146 }
147
148 for (int k = 0; k < pressure_stiffness.outerSize(); ++k)
149 {
150 for (StiffnessMatrix::InnerIterator it(pressure_stiffness, k); it; ++it)
151 {
152 blocks.emplace_back(n_bases * problem_dim + it.row(), n_bases * problem_dim + it.col(), it.value());
153 }
154 }
155
156 if (add_average)
157 {
158 const double val = 1.0 / n_pressure_bases;
159 for (int i = 0; i < n_pressure_bases; ++i)
160 {
161 blocks.emplace_back(n_bases * problem_dim + i, n_bases * problem_dim + n_pressure_bases, val);
162 blocks.emplace_back(n_bases * problem_dim + n_pressure_bases, n_bases * problem_dim + i, val);
163 }
164 }
165
166 stiffness.resize(n_bases * problem_dim + n_pressure_bases + avg_offset, n_bases * problem_dim + n_pressure_bases + avg_offset);
167 stiffness.setFromTriplets(blocks.begin(), blocks.end());
168 stiffness.makeCompressed();
169
170 // static int c = 0;
171 // Eigen::saveMarket(stiffness, "stiffness.txt");
172 // Eigen::saveMarket(velocity_stiffness, "velocity_stiffness.txt");
173 // Eigen::saveMarket(mixed_stiffness, "mixed_stiffness.txt");
174 // Eigen::saveMarket(pressure_stiffness, "pressure_stiffness.txt");
175 }
176
177 int AssemblerUtils::quadrature_order(const std::string &assembler, const int basis_degree, const BasisType &b_type, const int dim)
178 {
179 // note: minimum quadrature order is always 1
180 if (assembler == "Mass")
181 {
182 // multiply by two since we are multiplying phi_i by phi_j
183 if (b_type == BasisType::SIMPLEX_LAGRANGE || b_type == BasisType::CUBE_LAGRANGE)
184 return std::max(basis_degree * 2, 1);
185 else
186 return basis_degree * 2 + 1;
187 }
188 else if (assembler == "NavierStokes")
189 {
190 if (b_type == BasisType::SIMPLEX_LAGRANGE)
191 return std::max((basis_degree - 1) + basis_degree, 1);
192 else if (b_type == BasisType::CUBE_LAGRANGE)
193 return std::max(basis_degree * 2, 1);
194 else
195 return basis_degree * 2 + 1;
196 }
197 else
198 {
199 // subtract one since we take a derivative (lowers polynomial order by 1)
200 // multiply by two since we are multiplying grad phi_i by grad phi_j
201 if (b_type == BasisType::SIMPLEX_LAGRANGE) {
202 return std::max((basis_degree - 1) * 2, 1);
203 }
204 else if (b_type == BasisType::CUBE_LAGRANGE) {
205 // in this case we have a tensor product basis
206 // this computes the quadrature order along a single axis
207 // the Quadrature itself takes a tensor product of the given quadrature points
208 // to form the full quadrature for the basis
209 // taking a gradient leaves at least one variable whose power remains unchanged
210 // thus, we don't subtract 1
211 // note that this is overkill for the variable that was differentiated
212 return std::max(basis_degree * 2, 1);
213 }
214 else {
215 return (basis_degree - 1) * 2 + 1;
216 }
217 }
218 }
219
220 } // namespace assembler
221} // namespace polyfem
double val
Definition Assembler.cpp:86
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 std::shared_ptr< Assembler > make_assembler(const std::string &formulation)
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20