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>();
60 else if (formulation ==
"Bilaplacian")
61 return std::make_shared<BilaplacianMain>();
62 else if (formulation ==
"BilaplacianAux")
63 return std::make_shared<BilaplacianAux>();
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>();
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>();
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>();
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>();
125 const int n_bases,
const int n_pressure_bases,
const int problem_dim,
const bool add_average,
129 assert(velocity_stiffness.rows() == velocity_stiffness.cols());
130 assert(velocity_stiffness.rows() == n_bases * problem_dim);
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);
135 assert(pressure_stiffness.size() == 0 || pressure_stiffness.rows() == n_pressure_bases);
136 assert(pressure_stiffness.size() == 0 || pressure_stiffness.cols() == n_pressure_bases);
138 const int avg_offset = add_average ? 1 : 0;
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());
143 for (
int k = 0; k < velocity_stiffness.outerSize(); ++k)
145 for (StiffnessMatrix::InnerIterator it(velocity_stiffness, k); it; ++it)
147 blocks.emplace_back(it.row(), it.col(), it.value());
151 for (
int k = 0; k < mixed_stiffness.outerSize(); ++k)
153 for (StiffnessMatrix::InnerIterator it(mixed_stiffness, k); it; ++it)
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());
160 for (
int k = 0; k < pressure_stiffness.outerSize(); ++k)
162 for (StiffnessMatrix::InnerIterator it(pressure_stiffness, k); it; ++it)
164 blocks.emplace_back(n_bases * problem_dim + it.row(), n_bases * problem_dim + it.col(), it.value());
170 const double val = 1.0 / n_pressure_bases;
171 for (
int i = 0; i < n_pressure_bases; ++i)
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);
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();
192 if (assembler ==
"Mass")
196 return std::max(basis_degree * 2, 1);
198 return basis_degree * 2 + 1;
200 else if (assembler ==
"NavierStokes")
203 return std::max((basis_degree - 1) + basis_degree, 1);
205 return std::max(basis_degree * 2, 1);
207 return basis_degree * 2 + 1;
215 return std::max((basis_degree - 1) * 2, 1);
226 return std::max(basis_degree * 2, 1);
230 return (basis_degree - 1) * 2 + 1;
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,...