PolyFEM
Loading...
Searching...
No Matches
ElementAssemblyValues.cpp
Go to the documentation of this file.
4
5namespace polyfem
6{
7 using namespace basis;
8 using namespace quadrature;
9
10 namespace assembler
11 {
13 {
14 val = v;
15
17 det.resize(v.rows(), 1);
18
19 jac_it.resize(v.rows());
20 for (long i = 0; i < v.rows(); ++i)
21 jac_it[i] = Eigen::MatrixXd::Identity(v.cols(), v.cols());
22
23 det.setConstant(1); // volume (det of the geometric mapping)
24 for (std::size_t j = 0; j < basis_values.size(); ++j)
25 basis_values[j].grad_t_m = basis_values[j].grad; // / scaling
26 }
27
28 bool ElementAssemblyValues::is_geom_mapping_positive(const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy, const Eigen::MatrixXd &dz) const
29 {
30 Eigen::Matrix3d tmp;
31 for (long i = 0; i < dx.rows(); ++i)
32 {
33 tmp.row(0) = dx.row(i);
34 tmp.row(1) = dy.row(i);
35 tmp.row(2) = dz.row(i);
36
37 if (tmp.determinant() <= 0)
38 {
39 // std::cout<<tmp.determinant()<<std::endl;
40 return false;
41 }
42 }
43
44 return true;
45 }
46
47 bool ElementAssemblyValues::is_geom_mapping_positive(const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy) const
48 {
49 Eigen::Matrix2d tmp;
50 for (long i = 0; i < dx.rows(); ++i)
51 {
52 tmp.row(0) = dx.row(i);
53 tmp.row(1) = dy.row(i);
54
55 if (tmp.determinant() <= 0)
56 {
57 // std::cout<<tmp.determinant()<<std::endl;
58 return false;
59 }
60 }
61
62 return true;
63 }
64
65 void ElementAssemblyValues::finalize3d(const ElementBases &gbasis, const std::vector<AssemblyValues> &gbasis_values)
66 {
67 // if(det.size() != val.rows())
68 // logger().trace("Reallocating memory");
69 det.resize(val.rows(), 1);
70
71 for (std::size_t j = 0; j < basis_values.size(); ++j)
72 basis_values[j].finalize();
73
74 Eigen::Matrix3d tmp;
75 jac_it.resize(val.rows());
76
77 // loop over points
78 for (long k = 0; k < val.rows(); ++k)
79 {
80 tmp.setZero();
81 for (int j = 0; j < gbasis_values.size(); ++j)
82 {
83 const Basis &b = gbasis.bases[j];
84 assert(gbasis.has_parameterization);
85 assert(gbasis_values[j].grad.rows() == val.rows());
86 assert(gbasis_values[j].grad.cols() == 3);
87
88 for (std::size_t ii = 0; ii < b.global().size(); ++ii)
89 {
90 tmp.row(0) += gbasis_values[j].grad(k, 0) * b.global()[ii].node * b.global()[ii].val;
91 tmp.row(1) += gbasis_values[j].grad(k, 1) * b.global()[ii].node * b.global()[ii].val;
92 tmp.row(2) += gbasis_values[j].grad(k, 2) * b.global()[ii].node * b.global()[ii].val;
93 }
94 }
95
96 det(k) = tmp.determinant();
97 // assert(det(k)>0);
98 // std::cout<<det(k)<<std::endl;
99
100 jac_it[k] = tmp.inverse().transpose();
101 for (std::size_t j = 0; j < basis_values.size(); ++j)
102 basis_values[j].grad_t_m.row(k) = basis_values[j].grad.row(k) * jac_it[k];
103 }
104 }
105
106 void ElementAssemblyValues::finalize2d(const ElementBases &gbasis, const std::vector<AssemblyValues> &gbasis_values)
107 {
108 det.resize(val.rows(), 1);
109
110 for (std::size_t j = 0; j < basis_values.size(); ++j)
111 basis_values[j].finalize();
112
113 Eigen::Matrix2d tmp;
114 jac_it.resize(val.rows());
115
116 // loop over points
117 for (long k = 0; k < val.rows(); ++k)
118 {
119 tmp.setZero();
120 for (int j = 0; j < gbasis_values.size(); ++j)
121 {
122 const Basis &b = gbasis.bases[j];
123 assert(gbasis.has_parameterization);
124 assert(gbasis_values[j].grad.rows() == val.rows());
125 assert(gbasis_values[j].grad.cols() == 2);
126
127 for (std::size_t ii = 0; ii < b.global().size(); ++ii)
128 {
129 // add given geometric basis function + node's contribution to the Jacobian
130 tmp.row(0) += gbasis_values[j].grad(k, 0) * b.global()[ii].node * b.global()[ii].val;
131 tmp.row(1) += gbasis_values[j].grad(k, 1) * b.global()[ii].node * b.global()[ii].val;
132 }
133 }
134
135 // save Jacobian's determinant
136 det(k) = tmp.determinant();
137 // assert(det(k)>0);
138 // std::cout<<det(k)<<std::endl;
139
140 // save Jacobian's inverse transpose
141 jac_it[k] = tmp.inverse().transpose();
142 for (std::size_t j = 0; j < basis_values.size(); ++j)
143 basis_values[j].grad_t_m.row(k) = basis_values[j].grad.row(k) * jac_it[k];
144 }
145 }
146
147 Eigen::VectorXd ElementAssemblyValues::eval_deformed_jacobian_determinant(const Eigen::VectorXd &disp) const
148 {
149 assert(basis_);
150 assert(gbasis_);
151 int order = 0;
152 for (const auto &b : basis_->bases)
153 order = std::max(order, b.order());
154 const Eigen::MatrixXd cp = utils::extract_nodes(is_volume_ ? 3 : 2, *basis_, *gbasis_, disp, order);
156 }
157
158 void ElementAssemblyValues::compute(const int el_index, const bool is_volume, const ElementBases &basis, const ElementBases &gbasis)
159 {
161 compute(el_index, is_volume, quadrature.points, basis, gbasis);
162 }
163
164 void ElementAssemblyValues::compute(const int el_index, const bool is_volume, const Eigen::MatrixXd &pts, const ElementBases &basis, const ElementBases &gbasis)
165 {
166 basis_ = &basis;
167 gbasis_ = &gbasis;
168
169 element_id = el_index;
170 is_volume_ = is_volume;
171 // const bool poly = !gbasis.has_parameterization;
172
173 basis_values.resize(basis.bases.size());
174
175 if (&basis != &gbasis)
176 g_basis_values_cache_.resize(gbasis.bases.size());
177
178 const int n_local_bases = int(basis.bases.size());
179 const int n_local_g_bases = int(gbasis.bases.size());
180
181 // evaluate on reference element
182 basis.evaluate_bases(pts, basis_values);
183 basis.evaluate_grads(pts, basis_values);
184
185 if (&basis != &gbasis)
186 {
189 }
190
191 for (int j = 0; j < n_local_bases; ++j)
192 {
193 AssemblyValues &ass_val = basis_values[j];
194 ass_val.global = basis.bases[j].global();
195 assert(ass_val.val.cols() == 1);
196 assert(ass_val.grad.cols() == pts.cols());
197 }
198
199 if (!gbasis.has_parameterization)
200 {
201 // v = G(pts)
203 return;
204 }
205
206 // compute geometric mapping as linear combination of geometric basis functions
207 const auto &gbasis_values = (&basis == &gbasis) ? basis_values : g_basis_values_cache_;
208 assert(gbasis_values.size() == n_local_g_bases);
209 val.resize(pts.rows(), pts.cols());
210 val.setZero();
211
212 // loop over geometric basis functions
213 for (int j = 0; j < n_local_g_bases; ++j)
214 {
215 const Basis &b = gbasis.bases[j];
216 const auto &tmp = gbasis_values[j].val;
217
218 assert(gbasis.has_parameterization);
219 assert(tmp.size() == val.rows());
220
221 // loop over relevant global nodes
222 for (std::size_t ii = 0; ii < b.global().size(); ++ii)
223 {
224 for (long k = 0; k < val.rows(); ++k)
225 {
226 // add contribution from geometric basis function + node
227 val.row(k) += tmp(k) * b.global()[ii].node * b.global()[ii].val;
228 }
229 }
230 }
231
232 // compute Jacobian
233 if (is_volume)
234 finalize3d(gbasis, gbasis_values);
235 else
236 finalize2d(gbasis, gbasis_values);
237 }
238
239 bool ElementAssemblyValues::is_geom_mapping_positive(const bool is_volume, const ElementBases &gbasis) const
240 {
241 if (!gbasis.has_parameterization)
242 return true;
243
244 const int n_local_bases = int(gbasis.bases.size());
245
246 if (n_local_bases <= 0)
247 return true;
248
249 Quadrature quad;
250 gbasis.compute_quadrature(quad);
251
252 std::vector<AssemblyValues> tmp;
253
254 Eigen::MatrixXd dxmv = Eigen::MatrixXd::Zero(quad.points.rows(), quad.points.cols());
255 Eigen::MatrixXd dymv = Eigen::MatrixXd::Zero(quad.points.rows(), quad.points.cols());
256 Eigen::MatrixXd dzmv;
257 if (is_volume)
258 dzmv = Eigen::MatrixXd::Zero(quad.points.rows(), quad.points.cols());
259
260 gbasis.evaluate_grads(quad.points, tmp);
261
262 for (int j = 0; j < n_local_bases; ++j)
263 {
264 const Basis &b = gbasis.bases[j];
265
266 // b.grad(quad.points, grad);
267 // assert(grad.cols() == quad.points.cols());
268
269 for (std::size_t ii = 0; ii < b.global().size(); ++ii)
270 {
271 for (long k = 0; k < quad.points.rows(); ++k)
272 {
273 dxmv.row(k) += tmp[j].grad(k, 0) * b.global()[ii].node * b.global()[ii].val;
274 dymv.row(k) += tmp[j].grad(k, 1) * b.global()[ii].node * b.global()[ii].val;
275 if (is_volume)
276 dzmv.row(k) += tmp[j].grad(k, 2) * b.global()[ii].node * b.global()[ii].val;
277 }
278 }
279 }
280
281 return is_volume ? is_geom_mapping_positive(dxmv, dymv, dzmv) : is_geom_mapping_positive(dxmv, dymv);
282 }
283 } // namespace assembler
284} // namespace polyfem
Quadrature quadrature
stores per local bases evaluations
std::vector< basis::Local2Global > global
bool is_geom_mapping_positive(const bool is_volume, const basis::ElementBases &gbasis) const
check if the element is flipped
std::vector< AssemblyValues > g_basis_values_cache_
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,...
void finalize2d(const basis::ElementBases &gbasis, const std::vector< AssemblyValues > &gbasis_values)
void finalize(const Eigen::MatrixXd &v, const Eigen::MatrixXd &dx, const Eigen::MatrixXd &dy); void f...
std::vector< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, 0, 3, 3 > > jac_it
Eigen::VectorXd eval_deformed_jacobian_determinant(const Eigen::VectorXd &disp) const
void finalize_global_element(const Eigen::MatrixXd &v)
void finalize3d(const basis::ElementBases &gbasis, const std::vector< AssemblyValues > &gbasis_values)
Represents one basis function and its gradient.
Definition Basis.hpp:44
const std::vector< Local2Global > & global() const
Definition Basis.hpp:104
Stores the basis functions for a given element in a mesh (facet in 2d, cell in 3d).
void compute_quadrature(quadrature::Quadrature &quadrature) const
quadrature points to evaluate the basis functions inside the element
void evaluate_bases(const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values) const
evaluate stored bases at given points on the reference element saves results to basis_values
void evaluate_grads(const Eigen::MatrixXd &uv, std::vector< assembler::AssemblyValues > &basis_values) const
evaluate gradient of stored bases at given points on the reference element saves results to basis_val...
std::vector< Basis > bases
one basis function per node in the element
Eigen::VectorXd robust_evaluate_jacobian(const int order, const Eigen::MatrixXd &cp, const Eigen::MatrixXd &uv)
Definition Jacobian.cpp:124
Eigen::MatrixXd extract_nodes(const int dim, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::VectorXd &u, int order, int n_elem)
Definition Jacobian.cpp:77