PolyFEM
Loading...
Searching...
No Matches
Mass.cpp
Go to the documentation of this file.
1#include "Mass.hpp"
2
3namespace polyfem::assembler
4{
5 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1> Mass::assemble(const LinearAssemblerData &data) const
6 {
7 double tmp = 0;
8
9 // loop over quadrature points
10 for (int q = 0; q < data.da.size(); ++q)
11 {
12 const double rho = density_(data.vals.quadrature.points.row(q), data.vals.val.row(q), data.t, data.vals.element_id);
13 // phi_i * phi_j weighted by quadrature weights
14 tmp += rho * data.vals.basis_values[data.i].val(q) * data.vals.basis_values[data.j].val(q) * data.da(q);
15 }
16
17 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1> res(size() * size(), 1);
18 res.setZero();
19 for (int i = 0; i < size(); ++i)
20 res(i * size() + i) = tmp;
21
22 return res;
23 }
24
25 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1> Mass::compute_rhs(const AutodiffHessianPt &pt) const
26 {
27 assert(false);
28 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1> result;
29
30 return result;
31 }
32
33 void Mass::add_multimaterial(const int index, const json &params, const Units &units)
34 {
35 assert(size_ == 1 || size_ == 2 || size_ == 3);
36
37 density_.add_multimaterial(index, params, units.density());
38 }
39
40 std::map<std::string, Assembler::ParamFunc> Mass::parameters() const
41 {
42 std::map<std::string, ParamFunc> res;
43 res["rho"] = [this](const RowVectorNd &uv, const RowVectorNd &p, double t, int e) { return this->density_(uv, p, t, e); };
44
45 return res;
46 }
47
48} // namespace polyfem::assembler
std::string density() const
Definition Units.hpp:26
virtual void add_multimaterial(const int index, const json &params, const std::string &density_unit)
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
void add_multimaterial(const int index, const json &params, const Units &units) override
inialize material parameter
Definition Mass.cpp:33
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > compute_rhs(const AutodiffHessianPt &pt) const override
uses autodiff to compute the rhs for a fabricated solution in this case it just return pt....
Definition Mass.cpp:25
virtual std::map< std::string, ParamFunc > parameters() const override
Definition Mass.cpp:40
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
computes and returns local stiffness matrix (1x1) for bases i,j (where i,j is passed in through data)...
Definition Mass.cpp:5
Used for test only.
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13