Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
HookeLinearElasticity.cpp
Go to the documentation of this file.
2
4
5namespace polyfem::assembler
6{
7 using namespace basis;
8 namespace
9 {
10 template <class Matrix>
11 Matrix strain_from_disp_grad(const Matrix &disp_grad)
12 {
13 Matrix mat = (disp_grad + disp_grad.transpose());
14
15 for (int i = 0; i < mat.size(); ++i)
16 mat(i) *= 0.5;
17
18 return mat;
19 }
20
21 template <int dim>
22 Eigen::Matrix<double, dim, dim> strain(
23 const Eigen::MatrixXd &grad,
24 const Eigen::MatrixXd &jac_it,
25 int k, int coo)
26 {
27 Eigen::Matrix<double, dim, dim> jac;
28 jac.setZero();
29 jac.row(coo) = grad.row(k);
30 jac = jac * jac_it;
31
32 return strain_from_disp_grad(jac);
33 }
34
35 template <typename T, unsigned long N>
36 T stress(const ElasticityTensor &elasticity_tensor, const std::array<T, N> &strain, const int j)
37 {
38 T res = elasticity_tensor(j, 0) * strain[0];
39
40 for (unsigned long k = 1; k < N; ++k)
41 res += elasticity_tensor(j, k) * strain[k];
42
43 return res;
44 }
45 } // namespace
46
50
51 void HookeLinearElasticity::add_multimaterial(const int index, const json &params, const Units &units)
52 {
53 assert(size() == 2 || size() == 3);
54
55 if (!params.contains("elasticity_tensor") || params["elasticity_tensor"].empty())
56 {
57 if (params.count("young"))
58 {
59 if (params["young"].is_number() && params["nu"].is_number())
60 elasticity_tensor_.set_from_young_poisson(params["young"], params["nu"], units.stress());
61 }
62 else if (params.count("E"))
63 {
64 if (params["E"].is_number() && params["nu"].is_number())
65 elasticity_tensor_.set_from_young_poisson(params["E"], params["nu"], units.stress());
66 }
67 else if (params.count("lambda"))
68 {
69 if (params["lambda"].is_number() && params["mu"].is_number())
70 elasticity_tensor_.set_from_lambda_mu(params["lambda"], params["mu"], units.stress());
71 }
72 }
73 else
74 {
75 std::vector<double> entries = params["elasticity_tensor"];
77 }
78 if (params.contains("fiber_direction"))
79 fiber_direction_.add_multimaterial(index, params["fiber_direction"], units.length());
80 }
81
88
89 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>
91 {
92 const Eigen::MatrixXd &gradi = data.vals.basis_values[data.i].grad;
93 const Eigen::MatrixXd &gradj = data.vals.basis_values[data.j].grad;
94
95 // (C : gradi) : gradj
96 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1> res(size() * size());
97 res.setZero();
98 assert(gradi.cols() == size());
99 assert(gradj.cols() == size());
100 assert(size_t(gradi.rows()) == data.vals.jac_it.size());
101
102 for (long k = 0; k < gradi.rows(); ++k)
103 {
104 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1> res_k(size() * size());
105
108 {
109 const auto fdir_mtx_voigt = fiber_direction_.stiffness_rotation_voigt(data.vals.quadrature.points.row(k), data.vals.val.row(k), data.t, data.vals.element_id);
110
112 elasticity_tensor.rotate_stiffness(fdir_mtx_voigt);
113 }
114 else
115 {
117 }
118
119 if (size() == 2)
120 {
121 const Eigen::Matrix2d eps_x_i = strain<2>(gradi, data.vals.jac_it[k], k, 0);
122 const Eigen::Matrix2d eps_y_i = strain<2>(gradi, data.vals.jac_it[k], k, 1);
123
124 const Eigen::Matrix2d eps_x_j = strain<2>(gradj, data.vals.jac_it[k], k, 0);
125 const Eigen::Matrix2d eps_y_j = strain<2>(gradj, data.vals.jac_it[k], k, 1);
126
127 std::array<double, 3> e_x, e_y;
128 e_x[0] = eps_x_i(0, 0);
129 e_x[1] = eps_x_i(1, 1);
130 e_x[2] = 2 * eps_x_i(0, 1);
131
132 e_y[0] = eps_y_i(0, 0);
133 e_y[1] = eps_y_i(1, 1);
134 e_y[2] = 2 * eps_y_i(0, 1);
135
136 Eigen::Matrix2d sigma_x;
137 sigma_x << elasticity_tensor.compute_stress<3>(e_x, 0), elasticity_tensor.compute_stress<3>(e_x, 2),
139
140 Eigen::Matrix2d sigma_y;
141 sigma_y << elasticity_tensor.compute_stress<3>(e_y, 0), elasticity_tensor.compute_stress<3>(e_y, 2),
143
144 res_k(0) = (sigma_x * eps_x_j).trace();
145 res_k(2) = (sigma_x * eps_y_j).trace();
146
147 res_k(1) = (sigma_y * eps_x_j).trace();
148 res_k(3) = (sigma_y * eps_y_j).trace();
149 }
150 else
151 {
152 const Eigen::Matrix3d eps_x_i = strain<3>(gradi, data.vals.jac_it[k], k, 0);
153 const Eigen::Matrix3d eps_y_i = strain<3>(gradi, data.vals.jac_it[k], k, 1);
154 const Eigen::Matrix3d eps_z_i = strain<3>(gradi, data.vals.jac_it[k], k, 2);
155
156 const Eigen::Matrix3d eps_x_j = strain<3>(gradj, data.vals.jac_it[k], k, 0);
157 const Eigen::Matrix3d eps_y_j = strain<3>(gradj, data.vals.jac_it[k], k, 1);
158 const Eigen::Matrix3d eps_z_j = strain<3>(gradj, data.vals.jac_it[k], k, 2);
159
160 std::array<double, 6> e_x, e_y, e_z;
161 e_x[0] = eps_x_i(0, 0);
162 e_x[1] = eps_x_i(1, 1);
163 e_x[2] = eps_x_i(2, 2);
164 e_x[3] = 2 * eps_x_i(1, 2);
165 e_x[4] = 2 * eps_x_i(0, 2);
166 e_x[5] = 2 * eps_x_i(0, 1);
167
168 e_y[0] = eps_y_i(0, 0);
169 e_y[1] = eps_y_i(1, 1);
170 e_y[2] = eps_y_i(2, 2);
171 e_y[3] = 2 * eps_y_i(1, 2);
172 e_y[4] = 2 * eps_y_i(0, 2);
173 e_y[5] = 2 * eps_y_i(0, 1);
174
175 e_z[0] = eps_z_i(0, 0);
176 e_z[1] = eps_z_i(1, 1);
177 e_z[2] = eps_z_i(2, 2);
178 e_z[3] = 2 * eps_z_i(1, 2);
179 e_z[4] = 2 * eps_z_i(0, 2);
180 e_z[5] = 2 * eps_z_i(0, 1);
181
182 Eigen::Matrix3d sigma_x;
186
187 Eigen::Matrix3d sigma_y;
191
192 Eigen::Matrix3d sigma_z;
196
197 res_k(0) = (sigma_x * eps_x_j).trace();
198 res_k(3) = (sigma_x * eps_y_j).trace();
199 res_k(6) = (sigma_x * eps_z_j).trace();
200
201 res_k(1) = (sigma_y * eps_x_j).trace();
202 res_k(4) = (sigma_y * eps_y_j).trace();
203 res_k(7) = (sigma_y * eps_z_j).trace();
204
205 res_k(2) = (sigma_z * eps_x_j).trace();
206 res_k(5) = (sigma_z * eps_y_j).trace();
207 res_k(8) = (sigma_z * eps_z_j).trace();
208 }
209
210 res += res_k * data.da(k);
211 }
212
213 return res;
214 }
215
217 const OutputData &data,
218 const int all_size,
219 const ElasticityTensorType &type,
220 Eigen::MatrixXd &all,
221 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const
222 {
223 const auto &displacement = data.fun;
224 const auto &local_pts = data.local_pts;
225 const auto &bs = data.bs;
226 const auto &gbs = data.gbs;
227 const auto el_id = data.el_id;
228
229 Eigen::MatrixXd displacement_grad(size(), size());
230
231 assert(displacement.cols() == 1);
232
233 all.resize(local_pts.rows(), all_size);
234
236 vals.compute(el_id, size() == 3, local_pts, bs, gbs);
237
238 for (long p = 0; p < local_pts.rows(); ++p)
239 {
240 compute_diplacement_grad(size(), bs, vals, local_pts, p, displacement, displacement_grad);
241
244 {
245 const auto fdir_mtx_voigt = fiber_direction_.stiffness_rotation_voigt(local_pts.row(p), vals.val.row(p), data.t, el_id);
246
248 elasticity_tensor.rotate_stiffness(fdir_mtx_voigt);
249 }
250 else
251 {
253 }
254
255 if (type == ElasticityTensorType::F)
256 {
257 all.row(p) = fun(displacement_grad + Eigen::MatrixXd::Identity(size(), size()));
258 continue;
259 }
260
261 Eigen::MatrixXd strain = (displacement_grad + displacement_grad.transpose()) / 2;
262 Eigen::MatrixXd sigma(size(), size());
263
264 if (size() == 2)
265 {
266 std::array<double, 3> eps;
267 eps[0] = strain(0, 0);
268 eps[1] = strain(1, 1);
269 eps[2] = 2 * strain(0, 1);
270
273 }
274 else
275 {
276 std::array<double, 6> eps;
277 eps[0] = strain(0, 0);
278 eps[1] = strain(1, 1);
279 eps[2] = strain(2, 2);
280 eps[3] = 2 * strain(1, 2);
281 eps[4] = 2 * strain(0, 2);
282 eps[5] = 2 * strain(0, 1);
283
287 }
288
289 if (type == ElasticityTensorType::PK1)
290 sigma = pk1_from_cauchy(sigma, displacement_grad + Eigen::MatrixXd::Identity(size(), size()));
291 else if (type == ElasticityTensorType::PK2)
292 sigma = pk2_from_cauchy(sigma, displacement_grad + Eigen::MatrixXd::Identity(size(), size()));
293
294 all.row(p) = fun(sigma);
295 }
296 }
297
298 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1>
300 {
301 assert(pt.size() == size());
302 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1> res;
303
304 if (size() == 2)
306 else if (size() == 3)
308 else
309 assert(false);
310
311 return res;
312 }
313
314 std::map<std::string, Assembler::ParamFunc> HookeLinearElasticity::parameters() const
315 {
316 std::map<std::string, ParamFunc> res;
317 const auto &elast_tensor = elasticity_tensor();
318 const int size = this->size() == 2 ? 3 : 6;
319
320 for (int i = 0; i < size; ++i)
321 {
322 for (int j = i; j < size; ++j)
323 {
324 res[fmt::format("C_{}{}", i, j)] = [&elast_tensor, i, j](const RowVectorNd &, const RowVectorNd &, double, int) {
325 return elast_tensor(i, j);
326 };
327 }
328 }
329 return res;
330 }
331
333 {
334 return compute_energy_aux<double>(data);
335 }
336
338 {
339 const int n_bases = data.vals.basis_values.size();
341 size(), n_bases, data,
342 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 6, 1>>>(data); },
343 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 8, 1>>>(data); },
344 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 12, 1>>>(data); },
345 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 18, 1>>>(data); },
346 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 24, 1>>>(data); },
347 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 30, 1>>>(data); },
348 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 60, 1>>>(data); },
349 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, 81, 1>>>(data); },
350 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, SMALL_N, 1>>>(data); },
351 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, BIG_N, 1>>>(data); },
352 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar1<double, Eigen::VectorXd>>(data); });
353 }
354
356 {
357 const int n_bases = data.vals.basis_values.size();
359 size(), n_bases, data,
360 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 6, 1>, Eigen::Matrix<double, 6, 6>>>(data); },
361 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 8, 1>, Eigen::Matrix<double, 8, 8>>>(data); },
362 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 12, 1>, Eigen::Matrix<double, 12, 12>>>(data); },
363 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 18, 1>, Eigen::Matrix<double, 18, 18>>>(data); },
364 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 24, 1>, Eigen::Matrix<double, 24, 24>>>(data); },
365 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 30, 1>, Eigen::Matrix<double, 30, 30>>>(data); },
366 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 60, 1>, Eigen::Matrix<double, 60, 60>>>(data); },
367 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, 81, 1>, Eigen::Matrix<double, 81, 81>>>(data); },
368 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::Matrix<double, Eigen::Dynamic, 1, 0, SMALL_N, 1>, Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, 0, SMALL_N, SMALL_N>>>(data); },
369 [&](const NonLinearAssemblerData &data) { return compute_energy_aux<DScalar2<double, Eigen::VectorXd, Eigen::MatrixXd>>(data); });
370 }
371
372 template <typename T>
374 {
375 typedef Eigen::Matrix<T, Eigen::Dynamic, 1> AutoDiffVect;
376 typedef Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3> AutoDiffGradMat;
377
378 AutoDiffVect local_disp;
379 get_local_disp(data, size(), local_disp);
380
381 AutoDiffGradMat disp_grad(size(), size());
382
383 T energy = T(0.0);
384
385 const int n_pts = data.da.size();
386 for (long p = 0; p < n_pts; ++p)
387 {
388 compute_disp_grad_at_quad(data, local_disp, p, size(), disp_grad);
389
390 AutoDiffGradMat strain = strain_from_disp_grad(disp_grad);
391 AutoDiffGradMat stress_tensor(size(), size());
392
393 if (size() == 2)
394 {
395 std::array<T, 3> eps;
396 eps[0] = strain(0, 0);
397 eps[1] = strain(1, 1);
398 eps[2] = 2 * strain(0, 1);
399
400 stress_tensor << stress(elasticity_tensor_, eps, 0), stress(elasticity_tensor_, eps, 2),
401 stress(elasticity_tensor_, eps, 2), stress(elasticity_tensor_, eps, 1);
402 }
403 else
404 {
405 std::array<T, 6> eps;
406 eps[0] = strain(0, 0);
407 eps[1] = strain(1, 1);
408 eps[2] = strain(2, 2);
409 eps[3] = 2 * strain(1, 2);
410 eps[4] = 2 * strain(0, 2);
411 eps[5] = 2 * strain(0, 1);
412
413 stress_tensor << stress(elasticity_tensor_, eps, 0), stress(elasticity_tensor_, eps, 5), stress(elasticity_tensor_, eps, 4),
414 stress(elasticity_tensor_, eps, 5), stress(elasticity_tensor_, eps, 1), stress(elasticity_tensor_, eps, 3),
415 stress(elasticity_tensor_, eps, 4), stress(elasticity_tensor_, eps, 3), stress(elasticity_tensor_, eps, 2);
416 }
417
418 energy += (stress_tensor * strain).trace() * data.da(p);
419 }
420
421 return energy * 0.5;
422 }
423} // namespace polyfem::assembler
ElementAssemblyValues vals
Definition Assembler.cpp:22
std::vector< Eigen::Triplet< double > > entries
std::string stress() const
Definition Units.hpp:27
const std::string & length() const
Definition Units.hpp:19
virtual void set_size(const int size)
Definition Assembler.hpp:64
void rotate_stiffness(const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 6, 6 > &rotation_mtx_voigt)
void set_from_entries(const std::vector< double > &entries, const std::string &stress_unit)
void set_from_young_poisson(const double young, const double poisson, const std::string &stress_unit)
double compute_stress(const std::array< double, DIM > &strain, const int j) const
void set_from_lambda_mu(const double lambda, const double mu, const std::string &stress_unit)
stores per element basis values at given quadrature points and geometric mapping
void compute(const int el_index, const bool is_volume, const Eigen::MatrixXd &pts, const basis::ElementBases &basis, const basis::ElementBases &gbasis)
computes the per element values at the local (ref el) points (pts) sets basis_values,...
std::vector< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > > jac_it
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 1, 6, 6 > stiffness_rotation_voigt(double px, double py, double pz, double x, double y, double z, double t, int el_id) const
void add_multimaterial(const int index, const json &params, const std::string &unit)
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
local assembly function that defines the bilinear form (LHS) computes and returns a single local stif...
double compute_energy(const NonLinearAssemblerData &data) const override
std::map< std::string, ParamFunc > parameters() const override
const ElasticityTensor & elasticity_tensor() const
void add_multimaterial(const int index, const json &params, const Units &units) override
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
T compute_energy_aux(const NonLinearAssemblerData &data) const
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
void assign_stress_tensor(const OutputData &data, const int all_size, const ElasticityTensorType &type, Eigen::MatrixXd &all, const std::function< Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const override
const ElementAssemblyValues & vals
stores the evaluation for that element
const QuadratureVector & da
contains both the quadrature weight and the change of metric in the integral
const basis::ElementBases & bs
const Eigen::MatrixXd & fun
const Eigen::MatrixXd & local_pts
const basis::ElementBases & gbs
Used for test only.
Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > AutoDiffGradMat
Eigen::Matrix< T, Eigen::Dynamic, 1 > AutoDiffVect
void hooke_3d_function(const AutodiffHessianPt &pt, const assembler::ElasticityTensor &C, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > &res)
void hooke_2d_function(const AutodiffHessianPt &pt, const assembler::ElasticityTensor &C, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > &res)
Eigen::MatrixXd pk2_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
void compute_disp_grad_at_quad(const assembler::NonLinearAssemblerData &data, const AutoDiffVect &local_disp, const int p, const int size, AutoDiffGradMat &def_grad)
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
Eigen::MatrixXd hessian_from_energy(const int size, const int n_bases, const assembler::NonLinearAssemblerData &data, const std::function< DScalar2< double, Eigen::Matrix< double, 6, 1 >, Eigen::Matrix< double, 6, 6 > >(const assembler::NonLinearAssemblerData &)> &fun6, const std::function< DScalar2< double, Eigen::Matrix< double, 8, 1 >, Eigen::Matrix< double, 8, 8 > >(const assembler::NonLinearAssemblerData &)> &fun8, const std::function< DScalar2< double, Eigen::Matrix< double, 12, 1 >, Eigen::Matrix< double, 12, 12 > >(const assembler::NonLinearAssemblerData &)> &fun12, const std::function< DScalar2< double, Eigen::Matrix< double, 18, 1 >, Eigen::Matrix< double, 18, 18 > >(const assembler::NonLinearAssemblerData &)> &fun18, const std::function< DScalar2< double, Eigen::Matrix< double, 24, 1 >, Eigen::Matrix< double, 24, 24 > >(const assembler::NonLinearAssemblerData &)> &fun24, const std::function< DScalar2< double, Eigen::Matrix< double, 30, 1 >, Eigen::Matrix< double, 30, 30 > >(const assembler::NonLinearAssemblerData &)> &fun30, const std::function< DScalar2< double, Eigen::Matrix< double, 60, 1 >, Eigen::Matrix< double, 60, 60 > >(const assembler::NonLinearAssemblerData &)> &fun60, const std::function< DScalar2< double, Eigen::Matrix< double, 81, 1 >, Eigen::Matrix< double, 81, 81 > >(const assembler::NonLinearAssemblerData &)> &fun81, const std::function< DScalar2< double, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, SMALL_N, 1 >, Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, SMALL_N, SMALL_N > >(const assembler::NonLinearAssemblerData &)> &funN, const std::function< DScalar2< double, Eigen::VectorXd, Eigen::MatrixXd >(const assembler::NonLinearAssemblerData &)> &funn)
nlohmann::json json
Definition Common.hpp:9
void compute_diplacement_grad(const int size, const ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const int p, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &displacement_grad)
void get_local_disp(const assembler::NonLinearAssemblerData &data, const int size, AutoDiffVect &local_disp)
Eigen::MatrixXd pk1_from_cauchy(const Eigen::MatrixXd &stress, const Eigen::MatrixXd &F)
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
Eigen::VectorXd gradient_from_energy(const int size, const int n_bases, const assembler::NonLinearAssemblerData &data, const std::function< DScalar1< double, Eigen::Matrix< double, 6, 1 > >(const assembler::NonLinearAssemblerData &)> &fun6, const std::function< DScalar1< double, Eigen::Matrix< double, 8, 1 > >(const assembler::NonLinearAssemblerData &)> &fun8, const std::function< DScalar1< double, Eigen::Matrix< double, 12, 1 > >(const assembler::NonLinearAssemblerData &)> &fun12, const std::function< DScalar1< double, Eigen::Matrix< double, 18, 1 > >(const assembler::NonLinearAssemblerData &)> &fun18, const std::function< DScalar1< double, Eigen::Matrix< double, 24, 1 > >(const assembler::NonLinearAssemblerData &)> &fun24, const std::function< DScalar1< double, Eigen::Matrix< double, 30, 1 > >(const assembler::NonLinearAssemblerData &)> &fun30, const std::function< DScalar1< double, Eigen::Matrix< double, 60, 1 > >(const assembler::NonLinearAssemblerData &)> &fun60, const std::function< DScalar1< double, Eigen::Matrix< double, 81, 1 > >(const assembler::NonLinearAssemblerData &)> &fun81, const std::function< DScalar1< double, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, SMALL_N, 1 > >(const assembler::NonLinearAssemblerData &)> &funN, const std::function< DScalar1< double, Eigen::Matrix< double, Eigen::Dynamic, 1, 0, BIG_N, 1 > >(const assembler::NonLinearAssemblerData &)> &funBigN, const std::function< DScalar1< double, Eigen::VectorXd >(const assembler::NonLinearAssemblerData &)> &funn)