Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Electrostatics.cpp
Go to the documentation of this file.
1#include "Electrostatics.hpp"
2
3namespace polyfem::assembler
4{
5 namespace
6 {
7 bool delta(int i, int j)
8 {
9 return (i == j) ? true : false;
10 }
11 } // namespace
12
13 void Electrostatics::add_multimaterial(const int index, const json &params, const Units &units)
14 {
15 assert(size() == 1);
16
17 epsilon_.add_multimaterial(index, params, units.permittivity());
18 }
19
20 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1> Electrostatics::assemble(const LinearAssemblerData &data) const
21 {
22 const Eigen::MatrixXd &gradi = data.vals.basis_values[data.i].grad_t_m;
23 const Eigen::MatrixXd &gradj = data.vals.basis_values[data.j].grad_t_m;
24 // return ((gradi.array() * gradj.array()).rowwise().sum().array() * da.array()).colwise().sum();
25
26 double res = 0;
27 assert(gradi.rows() == data.da.size());
28 for (int k = 0; k < gradi.rows(); ++k)
29 {
30 double epsilon = epsilon_(data.vals.val.row(k), data.t, data.vals.element_id);
31 // compute grad(phi_i) dot grad(phi_j) weighted by quadrature weights
32 res += epsilon * gradi.row(k).dot(gradj.row(k)) * data.da(k);
33 }
34 return Eigen::Matrix<double, 1, 1>::Constant(res);
35 }
36
38 const Eigen::MatrixXd &mat,
39 Eigen::MatrixXd &stress,
40 Eigen::MatrixXd &result) const
41 {
42 stress = data.grad_u_i;
43 result = mat;
44 }
45
48 const Eigen::MatrixXd &local_pts,
49 const Eigen::MatrixXd &displacement,
50 Eigen::MatrixXd &tensor) const
51 {
52 const int dim = local_pts.cols();
53 tensor.resize(local_pts.rows(), dim * dim);
54 assert(displacement.cols() == 1);
55
56 for (long p = 0; p < local_pts.rows(); ++p)
57 {
58 for (int i = 0, idx = 0; i < dim; i++)
59 for (int j = 0; j < dim; j++)
60 {
61 tensor(p, idx) = delta(i, j);
62 idx++;
63 }
64 }
65 }
66
67 std::map<std::string, Assembler::ParamFunc> Electrostatics::parameters() const
68 {
69 std::map<std::string, ParamFunc> res;
70
71 res["epsilon"] = [this](const RowVectorNd &, const RowVectorNd &p, double t, int e) {
72 return epsilon_(p, t, e);
73 };
74
75 return res;
76 }
77
79 const bool is_volume,
80 const int n_basis,
81 const std::vector<basis::ElementBases> &bases,
82 const std::vector<basis::ElementBases> &gbases,
84 const double t,
85 const Eigen::MatrixXd &solution)
86 {
88 assemble(is_volume, n_basis, bases, gbases, cache, t, K);
89
90 Eigen::MatrixXd energy = 0.5 * solution.transpose() * (K * solution);
91 assert(energy.size() == 1);
92 return energy(0);
93 }
94
95} // namespace polyfem::assembler
std::unique_ptr< MatrixCache > cache
Definition Assembler.cpp:21
ElementAssemblyValues vals
Definition Assembler.cpp:22
std::string permittivity() const
Definition Units.hpp:25
Caches basis evaluation and geometric mapping at every element.
double compute_stored_energy(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 Eigen::MatrixXd &solution)
void compute_stress_grad_multiply_mat(const OptAssemblerData &data, const Eigen::MatrixXd &mat, Eigen::MatrixXd &stress, Eigen::MatrixXd &result) const override
void compute_stiffness_value(const double t, const assembler::ElementAssemblyValues &vals, const Eigen::MatrixXd &local_pts, const Eigen::MatrixXd &displacement, Eigen::MatrixXd &tensor) const override
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
computes local stiffness matrix (1x1) for bases i,j where i,j is passed in through data ie integral o...
void add_multimaterial(const int index, const json &params, const Units &units) override
std::map< std::string, ParamFunc > parameters() const override
stores per element basis values at given quadrature points and geometric mapping
void add_multimaterial(const int index, const json &params, const std::string &unit_type)
Definition MatParams.cpp:28
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
Used for test only.
nlohmann::json json
Definition Common.hpp:9
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22