PolyFEM
Loading...
Searching...
No Matches
Assembler.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Units.hpp>
4
7
12
13// this casses are instantiated in the cpp, cannot be used with generic assembler
14// without adding template instantiation
15namespace polyfem::assembler
17 // mixed formulation assembler
19 {
20 public:
22 virtual ~MixedAssembler() = default;
23
24 // this assembler takes two bases: psi_bases are the scalar ones, phi_bases are the tensor ones
25 // both have the same geometric mapping
27 const bool is_volume,
28 const int n_psi_basis,
29 const int n_phi_basis,
30 const std::vector<basis::ElementBases> &psi_bases,
31 const std::vector<basis::ElementBases> &phi_bases,
32 const std::vector<basis::ElementBases> &gbases,
33 const AssemblyValsCache &psi_cache,
34 const AssemblyValsCache &phi_cache,
35 const double t,
36 StiffnessMatrix &stiffness) const;
37
38 virtual std::string name() const = 0;
39
40 int size() const { return size_; }
41 virtual void set_size(const int size) { size_ = size; }
42
43 protected:
44 int size_ = -1;
45
46 virtual int rows() const = 0;
47 virtual int cols() const = 0;
48
49 virtual Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1> assemble(const MixedAssemblerData &data) const = 0;
50 };
51
54 {
55 public:
56 typedef std::pair<std::string, Eigen::MatrixXd> NamedMatrix;
57 typedef std::function<double(const RowVectorNd &, const RowVectorNd &, double, int)> ParamFunc;
58
59 virtual ~Assembler() = default;
60
61 virtual std::string name() const = 0;
62
63 int size() const { return size_; }
64 virtual void set_size(const int size) { size_ = size; }
65
66 // assembler stiffness matrix, is the mesh is volumetric, number of bases and bases (FE and geom)
67 // gbases and bases can be the same (ie isoparametric)
68 virtual void assemble(
69 const bool is_volume,
70 const int n_basis,
71 const std::vector<basis::ElementBases> &bases,
72 const std::vector<basis::ElementBases> &gbases,
74 const double t,
75 StiffnessMatrix &stiffness,
76 const bool is_mass = false) const { log_and_throw_error("Assembler not implemented by {}!", name()); }
77
78 // assemble energy
79 virtual double assemble_energy(
80 const bool is_volume,
81 const std::vector<basis::ElementBases> &bases,
82 const std::vector<basis::ElementBases> &gbases,
84 const double t,
85 const double dt,
86 const Eigen::MatrixXd &displacement,
87 const Eigen::MatrixXd &displacement_prev) const { log_and_throw_error("Assemble energy not implemented by {}!", name()); }
88
89 virtual Eigen::VectorXd assemble_energy_per_element(
90 const bool is_volume,
91 const std::vector<basis::ElementBases> &bases,
92 const std::vector<basis::ElementBases> &gbases,
94 const double t,
95 const double dt,
96 const Eigen::MatrixXd &displacement,
97 const Eigen::MatrixXd &displacement_prev) const { log_and_throw_error("Assemble energy not implemented by {}!", name()); }
98
99 // assemble gradient of energy (rhs)
100 virtual void assemble_gradient(
101 const bool is_volume,
102 const int n_basis,
103 const std::vector<basis::ElementBases> &bases,
104 const std::vector<basis::ElementBases> &gbases,
106 const double t,
107 const double dt,
108 const Eigen::MatrixXd &displacement,
109 const Eigen::MatrixXd &displacement_prev,
110 Eigen::MatrixXd &rhs) const { log_and_throw_error("Assemble grad not implemented by {}!", name()); }
111
112 // assemble hessian of energy (grad)
113 virtual void assemble_hessian(
114 const bool is_volume,
115 const int n_basis,
116 const bool project_to_psd,
117 const std::vector<basis::ElementBases> &bases,
118 const std::vector<basis::ElementBases> &gbases,
120 const double t,
121 const double dt,
122 const Eigen::MatrixXd &displacement,
123 const Eigen::MatrixXd &displacement_prev,
124 utils::MatrixCache &mat_cache,
125 StiffnessMatrix &grad) const { log_and_throw_error("Assemble hessian not implemented by {}!", name()); }
126
127 // plotting (eg von mises), assembler is the name of the formulation
129 const OutputData &data,
130 std::vector<NamedMatrix> &result) const {}
131
132 // computes tensor, assembler is the name of the formulation
134 const OutputData &data,
135 std::vector<NamedMatrix> &result) const
136 {
137 }
138
139 // computes tensor, assembler is the name of the formulation
141 const double t,
143 const Eigen::MatrixXd &local_pts,
144 const Eigen::MatrixXd &displacement,
145 Eigen::MatrixXd &tensor) const { log_and_throw_error("Not implemented!"); }
146
148 const OptAssemblerData &data,
149 Eigen::MatrixXd &dstress_dmu,
150 Eigen::MatrixXd &dstress_dlambda) const { log_and_throw_adjoint_error("Not implemented!"); }
151
153 const OptAssemblerData &data,
154 const Eigen::MatrixXd &mat,
155 Eigen::MatrixXd &stress,
156 Eigen::MatrixXd &result) const { log_and_throw_adjoint_error("Not implemented!"); }
157
159 const OptAssemblerData &data,
160 Eigen::MatrixXd &stress,
161 Eigen::MatrixXd &result) const
162 {
163 Eigen::MatrixXd unused;
164 compute_stress_grad_multiply_mat(data, Eigen::MatrixXd::Zero(data.grad_u_i.rows(), data.grad_u_i.cols()), stress, unused);
165 compute_stress_grad_multiply_mat(data, stress, unused, result);
166 }
167
169 const OptAssemblerData &data,
170 const Eigen::MatrixXd &vect,
171 Eigen::MatrixXd &stress,
172 Eigen::MatrixXd &result) const { log_and_throw_error("Not implemented!"); }
173
175 const OptAssemblerData &data,
176 const Eigen::MatrixXd &prev_grad_u_i,
177 Eigen::MatrixXd &stress,
178 Eigen::MatrixXd &result) const { log_and_throw_adjoint_error("Not implemented!"); }
180 const OptAssemblerData &data,
181 const Eigen::MatrixXd &prev_grad_u_i,
182 Eigen::MatrixXd &result) const { log_and_throw_adjoint_error("Not implemented!"); }
183
184 virtual std::map<std::string, ParamFunc> parameters() const = 0;
185 virtual VectorNd compute_rhs(const AutodiffHessianPt &pt) const { log_and_throw_error("Rhs not supported by {}!", name()); }
186
187 virtual Eigen::Matrix<AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1> kernel(const int dim, const AutodiffGradPt &rvect, const AutodiffScalarGrad &r) const { log_and_throw_error("Kernel not supported by {}!", name()); }
188
189 void set_materials(const std::vector<int> &body_ids, const json &body_params, const Units &units);
190 virtual void add_multimaterial(const int index, const json &params, const Units &units) {}
191
192 virtual void update_lame_params(const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus)
193 {
194 log_and_throw_error("Not implemented!");
195 }
196
197 virtual bool is_linear() const = 0;
198 virtual bool is_solution_displacement() const { return false; }
199 virtual bool is_fluid() const { return false; }
200 virtual bool is_tensor() const { return false; }
201
202 protected:
203 int size_ = -1;
204 };
205
208 class LinearAssembler : virtual public Assembler
209 {
210 public:
212 virtual ~LinearAssembler() = default;
213
221 const bool is_volume,
222 const int n_basis,
223 const std::vector<basis::ElementBases> &bases,
224 const std::vector<basis::ElementBases> &gbases,
226 const double t,
227 StiffnessMatrix &stiffness,
228 const bool is_mass = false) const override;
229
230 virtual bool is_linear() const override { return true; }
231
234 virtual Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1> assemble(const LinearAssemblerData &data) const = 0;
235 };
236
237 // non-linear assembler (eg neohookean elasticity)
238 class NLAssembler : virtual public Assembler
239 {
240 public:
241 virtual ~NLAssembler() = default;
242
243 // assemble energy
244 double assemble_energy(
245 const bool is_volume,
246 const std::vector<basis::ElementBases> &bases,
247 const std::vector<basis::ElementBases> &gbases,
249 const double t,
250 const double dt,
251 const Eigen::MatrixXd &displacement,
252 const Eigen::MatrixXd &displacement_prev) const override;
253
254 // assemble the energy per element
255 Eigen::VectorXd assemble_energy_per_element(
256 const bool is_volume,
257 const std::vector<basis::ElementBases> &bases,
258 const std::vector<basis::ElementBases> &gbases,
260 const double t,
261 const double dt,
262 const Eigen::MatrixXd &displacement,
263 const Eigen::MatrixXd &displacement_prev) const override;
264
265 // assemble gradient of energy (rhs)
267 const bool is_volume,
268 const int n_basis,
269 const std::vector<basis::ElementBases> &bases,
270 const std::vector<basis::ElementBases> &gbases,
272 const double t,
273 const double dt,
274 const Eigen::MatrixXd &displacement,
275 const Eigen::MatrixXd &displacement_prev,
276 Eigen::MatrixXd &rhs) const override;
277
278 // assemble hessian of energy (grad)
280 const bool is_volume,
281 const int n_basis,
282 const bool project_to_psd,
283 const std::vector<basis::ElementBases> &bases,
284 const std::vector<basis::ElementBases> &gbases,
286 const double t,
287 const double dt,
288 const Eigen::MatrixXd &displacement,
289 const Eigen::MatrixXd &displacement_prev,
290 utils::MatrixCache &mat_cache,
291 StiffnessMatrix &grad) const override;
292
293 virtual bool is_linear() const override { return false; }
294
295 protected:
296 // energy, gradient, and hessian used in newton method
297 virtual double compute_energy(const NonLinearAssemblerData &data) const = 0;
298 virtual Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const = 0;
299 virtual Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const = 0;
300 };
301
302 class ElasticityAssembler : virtual public Assembler
303 {
304 public:
306 virtual ~ElasticityAssembler() = default;
307
308 // plotting (eg von mises), assembler is the name of the formulation
310 const OutputData &data,
311 std::vector<NamedMatrix> &result) const override
312 {
313 result.clear();
314 Eigen::MatrixXd tmp;
316 result.emplace_back("von_mises", tmp);
317 }
318
319 // computes tensor, assembler is the name of the formulation
321 const OutputData &data,
322 std::vector<NamedMatrix> &result) const override
323 {
324 result.clear();
325 Eigen::MatrixXd cauchy, pk1, pk2, F;
326
331
332 result.emplace_back("cauchy_stess", cauchy);
333 result.emplace_back("pk1_stess", pk1);
334 result.emplace_back("pk2_stess", pk2);
335 result.emplace_back("F", F);
336 }
337
339 const ElasticityTensorType &type,
340 Eigen::MatrixXd &stresses) const
341 {
342 assign_stress_tensor(data, size() * size(), type, stresses, [&](const Eigen::MatrixXd &stress) {
343 Eigen::MatrixXd tmp = stress;
344 auto a = Eigen::Map<Eigen::MatrixXd>(tmp.data(), 1, size() * size());
345 return Eigen::MatrixXd(a);
346 });
347 }
348
350 Eigen::MatrixXd &stresses) const
351 {
352 assign_stress_tensor(data, 1, ElasticityTensorType::CAUCHY, stresses, [&](const Eigen::MatrixXd &stress) {
353 Eigen::Matrix<double, 1, 1> res;
354 res.setConstant(von_mises_stress_for_stress_tensor(stress));
355 return res;
356 });
357 }
358
359 bool is_solution_displacement() const override { return true; }
360 bool is_tensor() const override { return true; }
361
362 protected:
363 virtual void assign_stress_tensor(const OutputData &data,
364 const int all_size,
365 const ElasticityTensorType &type,
366 Eigen::MatrixXd &all,
367 const std::function<Eigen::MatrixXd(const Eigen::MatrixXd &)> &fun) const = 0;
368 };
369} // namespace polyfem::assembler
std::unique_ptr< MatrixCache > cache
Definition Assembler.cpp:21
ElementAssemblyValues vals
Definition Assembler.cpp:22
Eigen::MatrixXd mat
virtual void compute_stress_grad_multiply_stress(const OptAssemblerData &data, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
virtual void assemble(const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, StiffnessMatrix &stiffness, const bool is_mass=false) const
Definition Assembler.hpp:68
void set_materials(const std::vector< int > &body_ids, const json &body_params, const Units &units)
Definition Assembler.cpp:97
virtual void compute_stiffness_value(const double t, const assembler::ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &tensor) const
std::pair< std::string, Eigen::MatrixXd > NamedMatrix
Definition Assembler.hpp:56
virtual ~Assembler()=default
virtual void compute_dstress_dmu_dlambda(const OptAssemblerData &data, Eigen::MatrixXd &dstress_dmu, Eigen::MatrixXd &dstress_dlambda) const
virtual void compute_stress_grad_multiply_vect(const OptAssemblerData &data, const Eigen::MatrixXd &vect, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
virtual bool is_solution_displacement() const
virtual void update_lame_params(const Eigen::MatrixXd &lambdas, const Eigen::MatrixXd &mus)
virtual void compute_stress_grad_multiply_mat(const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
virtual Eigen::Matrix< AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1 > kernel(const int dim, const AutodiffGradPt &rvect, const AutodiffScalarGrad &r) const
virtual bool is_tensor() const
std::function< double(const RowVectorNd &, const RowVectorNd &, double, int)> ParamFunc
Definition Assembler.hpp:57
virtual std::string name() const =0
virtual Eigen::VectorXd assemble_energy_per_element(const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const
Definition Assembler.hpp:89
virtual void assemble_gradient(const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, Eigen::MatrixXd &rhs) const
virtual void compute_stress_grad(const OptAssemblerData &data, const Eigen::MatrixXd &prev_grad_u_i, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const
virtual void set_size(const int size)
Definition Assembler.hpp:64
virtual void compute_scalar_value(const OutputData &data, std::vector< NamedMatrix > &result) const
virtual void assemble_hessian(const bool is_volume, const int n_basis, const bool project_to_psd, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, utils::MatrixCache &mat_cache, StiffnessMatrix &grad) const
virtual std::map< std::string, ParamFunc > parameters() const =0
virtual VectorNd compute_rhs(const AutodiffHessianPt &pt) const
virtual double assemble_energy(const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const
Definition Assembler.hpp:79
virtual bool is_fluid() const
virtual bool is_linear() const =0
virtual void compute_tensor_value(const OutputData &data, std::vector< NamedMatrix > &result) const
virtual void add_multimaterial(const int index, const json &params, const Units &units)
virtual void compute_stress_prev_grad(const OptAssemblerData &data, const Eigen::MatrixXd &prev_grad_u_i, Eigen::MatrixXd &result) const
Caches basis evaluation and geometric mapping at every element.
void compute_von_mises_stresses(const OutputData &data, Eigen::MatrixXd &stresses) const
void compute_scalar_value(const OutputData &data, std::vector< NamedMatrix > &result) const override
bool is_solution_displacement() const override
void compute_tensor_value(const OutputData &data, std::vector< NamedMatrix > &result) const override
void compute_stress_tensor(const OutputData &data, const ElasticityTensorType &type, Eigen::MatrixXd &stresses) const
virtual 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 =0
stores per element basis values at given quadrature points and geometric mapping
assemble matrix based on the local assembler local assembler is eg Laplace, LinearElasticity etc
virtual Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const =0
local assembly function that defines the bilinear form (LHS) computes and returns a single local stif...
void assemble(const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, StiffnessMatrix &stiffness, const bool is_mass=false) const override
assembles the stiffness matrix for the given basis the bilinear form (local assembler) is encoded by ...
virtual bool is_linear() const override
virtual int rows() const =0
virtual void set_size(const int size)
Definition Assembler.hpp:41
virtual Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > assemble(const MixedAssemblerData &data) const =0
void assemble(const bool is_volume, const int n_psi_basis, const int n_phi_basis, const std::vector< basis::ElementBases > &psi_bases, const std::vector< basis::ElementBases > &phi_bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &psi_cache, const AssemblyValsCache &phi_cache, const double t, StiffnessMatrix &stiffness) const
virtual int cols() const =0
virtual std::string name() const =0
virtual bool is_linear() const override
double assemble_energy(const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const override
Eigen::VectorXd assemble_energy_per_element(const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const override
virtual double compute_energy(const NonLinearAssemblerData &data) const =0
virtual ~NLAssembler()=default
void assemble_gradient(const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, Eigen::MatrixXd &rhs) const override
void assemble_hessian(const bool is_volume, const int n_basis, const bool project_to_psd, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, utils::MatrixCache &mat_cache, StiffnessMatrix &grad) const override
virtual Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const =0
virtual Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const =0
abstract class used for caching
Used for test only.
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:9
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
double von_mises_stress_for_stress_tensor(const Eigen::MatrixXd &stress)
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:11
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
Eigen::Matrix< AutodiffScalarGrad, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffGradPt
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20
Automatic differentiation scalar with first-order derivatives.
Definition autodiff.h:112