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