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>();
59 else if (formulation ==
"Bilaplacian")
60 return std::make_shared<BilaplacianMain>();
61 else if (formulation ==
"BilaplacianAux")
62 return std::make_shared<BilaplacianAux>();
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>();
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>();
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>();
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>();
122 const int n_bases,
const int n_pressure_bases,
const int problem_dim,
const bool add_average,
126 assert(velocity_stiffness.rows() == velocity_stiffness.cols());
127 assert(velocity_stiffness.rows() == n_bases * problem_dim);
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);
132 assert(pressure_stiffness.size() == 0 || pressure_stiffness.rows() == n_pressure_bases);
133 assert(pressure_stiffness.size() == 0 || pressure_stiffness.cols() == n_pressure_bases);
135 const int avg_offset = add_average ? 1 : 0;
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());
140 for (
int k = 0; k < velocity_stiffness.outerSize(); ++k)
142 for (StiffnessMatrix::InnerIterator it(velocity_stiffness, k); it; ++it)
144 blocks.emplace_back(it.row(), it.col(), it.value());
148 for (
int k = 0; k < mixed_stiffness.outerSize(); ++k)
150 for (StiffnessMatrix::InnerIterator it(mixed_stiffness, k); it; ++it)
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());
157 for (
int k = 0; k < pressure_stiffness.outerSize(); ++k)
159 for (StiffnessMatrix::InnerIterator it(pressure_stiffness, k); it; ++it)
161 blocks.emplace_back(n_bases * problem_dim + it.row(), n_bases * problem_dim + it.col(), it.value());
167 const double val = 1.0 / n_pressure_bases;
168 for (
int i = 0; i < n_pressure_bases; ++i)
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);
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();
189 if (assembler ==
"Mass")
193 return std::max(basis_degree * 2, 1);
195 return basis_degree * 2 + 1;
197 else if (assembler ==
"NavierStokes")
200 return std::max((basis_degree - 1) + basis_degree, 1);
202 return std::max(basis_degree * 2, 1);
204 return basis_degree * 2 + 1;
212 return std::max((basis_degree - 1) * 2, 1);
223 return std::max(basis_degree * 2, 1);
227 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,...