Flow123d  DF_patch_fevalues-8016b85
fe_values.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file fe_values.cc
15  * @brief Class FEValues calculates finite element data on the actual
16  * cells such as shape function values, gradients, Jacobian of
17  * the mapping from the reference cell etc.
18  * @author Jan Stebel
19  */
20 
21 #include "fem/mapping_p1.hh"
22 #include "quadrature/quadrature.hh"
23 #include "fem/element_values.hh"
24 #include "fem/finite_element.hh"
25 #include "fem/fe_values.hh"
26 #include "fem/patch_fe_values.hh"
27 #include "fem/fe_system.hh"
28 #include "fem/fe_values_map.hh"
29 
30 
31 
32 using namespace arma;
33 using namespace std;
34 
35 
36 
37 
38 
39 
40 
41 template<class FV, unsigned int spacedim>
43 : dim_(-1), n_points_(0), n_dofs_(0)
44 {
45 }
46 
47 
48 template<class FV, unsigned int spacedim>
50  : n_points(np),
51  n_dofs(nd)
52 {
53  ref_shape_values.resize(np, vector<arma::vec>(nd));
54  ref_shape_grads.resize(np, vector<arma::mat>(nd));
55 }
56 
57 
58 template<class FV, unsigned int spacedim>
60  const std::vector<unsigned int> &dof_indices,
61  unsigned int first_component_idx,
62  unsigned int ncomps)
63  : FEInternalData(fe_system_data.n_points, dof_indices.size())
64 {
65  for (unsigned int ip=0; ip<n_points; ip++)
66  for (unsigned int id=0; id<dof_indices.size(); id++)
67  {
68  ref_shape_values[ip][id] = fe_system_data.ref_shape_values[ip][dof_indices[id]].subvec(first_component_idx, first_component_idx+ncomps-1);
69  ref_shape_grads[ip][id] = fe_system_data.ref_shape_grads[ip][dof_indices[id]].cols(first_component_idx, first_component_idx+ncomps-1);
70  }
71 }
72 
73 
74 
75 template<class FV, unsigned int spacedim>
76 template<unsigned int DIM>
78 {
79  scalars.clear();
80  vectors.clear();
81  tensors.clear();
82  switch (fe.type_) {
83  case FEType::FEScalar:
84  scalars.push_back(FEValuesViews::Scalar<FV, spacedim>(fv, 0));
85  break;
86  case FEType::FEVector:
89  vectors.push_back(FEValuesViews::Vector<FV, spacedim>(fv, 0));
90  break;
91  case FEType::FETensor:
92  tensors.push_back(FEValuesViews::Tensor<FV, spacedim>(fv, 0));
93  break;
95  const FESystem<DIM> *fe_sys = dynamic_cast<const FESystem<DIM>*>(&fe);
96  ASSERT(fe_sys != nullptr).error("Mixed system must be represented by FESystem.");
97 
98  // Loop through sub-elements and add views according to their types.
99  // Note that the component index is calculated using fe->n_space_components(),
100  // not fe->n_components()!
101  unsigned int comp_offset = 0;
102  for (auto fe : fe_sys->fe())
103  {
104  switch (fe->type_)
105  {
106  case FEType::FEScalar:
107  scalars.push_back(FEValuesViews::Scalar<FV, spacedim>(fv,comp_offset));
108  break;
109  case FEType::FEVector:
112  vectors.push_back(FEValuesViews::Vector<FV, spacedim>(fv,comp_offset));
113  break;
114  case FEType::FETensor:
115  tensors.push_back(FEValuesViews::Tensor<FV, spacedim>(fv,comp_offset));
116  break;
117  default:
118  ASSERT(false).error("Not implemented.");
119  break;
120  }
121 
122  comp_offset += fe->n_space_components(spacedim);
123  }
124  break;
125  }
126 }
127 
128 
129 
130 template<class FV, unsigned int spacedim>
131 template<unsigned int DIM>
133  Quadrature &q,
134  FiniteElement<DIM> &_fe,
135  UpdateFlags _flags)
136 {
137  if (DIM == 0) //return; // avoid unnecessary allocation of dummy 0 dimensional objects
138  ASSERT(q.size() == 1);
139 
140  this->allocate( q, _fe, _flags);
141  this->initialize_in(q, DIM);
142 
143  // In case of mixed system allocate data for sub-elements.
144  if (fe_type_ == FEMixedSystem)
145  {
146  FESystem<DIM> *fe = dynamic_cast<FESystem<DIM>*>(&_fe);
147  ASSERT(fe != nullptr).error("Mixed system must be represented by FESystem.");
148 
149  fe_values_vec.resize(fe->fe().size());
150  init_fe_val_vec();
151  for (unsigned int f=0; f<fe->fe().size(); f++)
152  fe_values_vec[f].initialize(q, *fe->fe()[f], update_flags);
153  }
154 
155  // precompute finite element data
156  if ( q.dim() == DIM )
157  {
158  fe_data_ = init_fe_data(_fe, q);
159  }
160  else if ( q.dim() + 1 == DIM )
161  {
163  for (unsigned int sid = 0; sid < RefElement<DIM>::n_sides; sid++)
164  {
165  side_fe_data_[sid] = init_fe_data(_fe, q.make_from_side<DIM>(sid));
166  }
167  }
168  else
169  ASSERT(false)(q.dim())(DIM).error("Dimension mismatch in FEValues::initialize().");
170 }
171 
172 
173 
174 template<class FV, unsigned int spacedim>
175 template<unsigned int DIM>
177  Quadrature &_q,
178  FiniteElement<DIM> & _fe,
179  UpdateFlags _flags)
180 {
181  // For FEVector and FETensor check number of components.
182  // This cannot be done in FiniteElement since it does not know spacedim.
183  if (_fe.type_ == FEVector) {
184  ASSERT(_fe.n_components() == spacedim).error("FEVector must have spacedim components.");
185  } else if (_fe.type_ == FETensor) {
186  ASSERT(_fe.n_components() == spacedim*spacedim).error("FETensor must have spacedim*spacedim components.");
187  }
188 
189  fe_sys_dofs_.clear();
190  fe_sys_n_components_.clear();
192 
193  dim_ = DIM;
194  n_points_ = _q.size();
195  n_dofs_ = _fe.n_dofs();
196  n_components_ = _fe.n_space_components(spacedim);
197  fe_type_ = _fe.type_;
198  FESystem<DIM> *fe_sys = dynamic_cast<FESystem<DIM>*>(&_fe);
199  if (fe_sys != nullptr)
200  {
201  for (unsigned int f=0; f<fe_sys->fe().size(); f++)
202  {
203  fe_sys_dofs_.push_back(fe_sys->fe_dofs(f));
204  fe_sys_n_components_.push_back(fe_sys->fe()[f]->n_components());
205  fe_sys_n_space_components_.push_back(fe_sys->fe()[f]->n_space_components(spacedim));
206  }
207  }
208 
209  // add flags required by the finite element or mapping
210  update_flags = _flags | _fe.update_each(_flags);
212  this->allocate_in(_q.dim());
213 
214  views_cache_.initialize(*this->fv_, _fe);
215 }
216 
217 
218 
219 template<class FV, unsigned int spacedim>
220 template<unsigned int DIM>
221 std::shared_ptr<typename FEValuesBase<FV, spacedim>::FEInternalData> FEValuesBase<FV, spacedim>::init_fe_data(const FiniteElement<DIM> &fe, const Quadrature &q)
222 {
223  ASSERT( DIM == dim_ );
224  ASSERT( q.dim() == DIM );
225  std::shared_ptr<FEInternalData> data = std::make_shared<FEInternalData>(q.size(), n_dofs_);
226 
227  arma::mat shape_values(n_dofs_, fe.n_components());
228  for (unsigned int i=0; i<q.size(); i++)
229  {
230  for (unsigned int j=0; j<n_dofs_; j++)
231  {
232  for (unsigned int c=0; c<fe.n_components(); c++)
233  shape_values(j,c) = fe.shape_value(j, q.point<DIM>(i), c);
234 
235  data->ref_shape_values[i][j] = trans(shape_values.row(j));
236  }
237  }
238 
239  arma::mat grad(DIM, fe.n_components());
240  for (unsigned int i=0; i<q.size(); i++)
241  {
242  for (unsigned int j=0; j<n_dofs_; j++)
243  {
244  grad.zeros();
245  for (unsigned int c=0; c<fe.n_components(); c++)
246  grad.col(c) += fe.shape_grad(j, q.point<DIM>(i), c);
247 
248  data->ref_shape_grads[i][j] = grad;
249  }
250  }
251 
252  return data;
253 }
254 
255 
256 template<class FV, unsigned int spacedim>
258 {
259  switch (this->fe_type_) {
260  case FEScalar:
261  this->fill_data_specialized<MapScalar<FV, spacedim>>(elm_values, fe_data);
262  break;
263  case FEVector:
264  this->fill_data_specialized<MapVector<FV, spacedim>>(elm_values, fe_data);
265  break;
267  this->fill_data_specialized<MapContravariant<FV, spacedim>>(elm_values, fe_data);
268  break;
269  case FEVectorPiola:
270  this->fill_data_specialized<MapPiola<FV, spacedim>>(elm_values, fe_data);
271  break;
272  case FETensor:
273  this->fill_data_specialized<MapTensor<FV, spacedim>>(elm_values, fe_data);
274  break;
275  case FEMixedSystem:
276  this->fill_data_specialized<MapSystem<FV, spacedim>>(elm_values, fe_data);
277  break;
278  default:
279  ASSERT_PERMANENT(false).error("Not implemented.");
280  }
281 }
282 
283 
284 
285 template<class FV, unsigned int spacedim>
286 template<class MapType>
288  MapType map_type;
289  map_type.fill_values_vec(*this->fv_, elm_values, fe_data);
290  if (this->update_flags & update_values)
291  map_type.update_values(*this->fv_, elm_values, fe_data);
293  map_type.update_gradients(*this->fv_, elm_values, fe_data);
294 }
295 
296 
297 
298 
299 /*template<unsigned int spacedim>
300 double FEValues<spacedim>::shape_value_component(const unsigned int function_no,
301  const unsigned int point_no,
302  const unsigned int comp) const
303 {
304  ASSERT_LT(function_no, n_dofs_);
305  ASSERT_LT(point_no, n_points_);
306  ASSERT_LT(comp, n_components_);
307  return shape_values[point_no][function_no*n_components_+comp];
308 }*/
309 
310 
311 template<unsigned int spacedim>
313 : FEValuesBase<FEValues<spacedim>, spacedim>() {}
314 
315 
316 template<unsigned int spacedim>
318 }
319 
320 
321 template<unsigned int spacedim>
323  Quadrature &q,
324  unsigned int dim)
325 {
326  elm_values_ = std::make_shared<ElementValues<spacedim> >(q, this->update_flags, dim);
327 }
328 
329 
330 
331 template<unsigned int spacedim>
333 {
334  if (this->update_flags & update_values)
335  shape_values_.resize(this->n_points_, vector<double>(this->n_dofs_*this->n_components_));
336 
337  if (this->update_flags & update_gradients)
338  shape_gradients_.resize(this->n_points_, vector<arma::vec::fixed<spacedim> >(this->n_dofs_*this->n_components_));
339 
340  this->fv_ = this;
341 }
342 
343 
344 
345 template<unsigned int spacedim>
346 arma::vec::fixed<spacedim> FEValues<spacedim>::shape_grad_component(const unsigned int function_no,
347  const unsigned int point_no,
348  const unsigned int comp) const
349 {
350  ASSERT_LT(function_no, this->n_dofs_);
351  ASSERT_LT(point_no, this->n_points_);
352  ASSERT_LT(comp, this->n_components_);
353  return shape_gradients_[point_no][function_no*this->n_components_+comp];
354 }
355 
356 
357 /*template<unsigned int spacedim>
358 void FEValues<spacedim>::fill_scalar_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data)
359 {
360  ASSERT(fe_type_ == FEScalar);
361 
362  // shape values
363  if (update_flags & update_values)
364  for (unsigned int i = 0; i < fe_data.n_points; i++)
365  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
366  shape_values[i][j] = fe_data.ref_shape_values[i][j][0];
367 
368  // shape gradients
369  if (update_flags & update_gradients)
370  for (unsigned int i = 0; i < fe_data.n_points; i++)
371  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
372  shape_gradients[i][j] = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
373 }
374 
375 
376 template<unsigned int spacedim>
377 void FEValues<spacedim>::fill_vec_data(const ElementValues<spacedim> &elm_values,
378  const FEInternalData &fe_data)
379 {
380  ASSERT(fe_type_ == FEVector);
381 
382  // shape values
383  if (update_flags & update_values)
384  {
385  for (unsigned int i = 0; i < fe_data.n_points; i++)
386  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
387  {
388  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
389  for (unsigned int c=0; c<spacedim; c++)
390  shape_values[i][j*spacedim+c] = fv_vec[c];
391  }
392  }
393 
394  // shape gradients
395  if (update_flags & update_gradients)
396  {
397  for (unsigned int i = 0; i < fe_data.n_points; i++)
398  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
399  {
400  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
401  for (unsigned int c=0; c<spacedim; c++)
402  shape_gradients[i][j*spacedim+c] = grads.col(c);
403  }
404  }
405 }
406 
407 
408 template<unsigned int spacedim>
409 void FEValues<spacedim>::fill_vec_contravariant_data(const ElementValues<spacedim> &elm_values,
410  const FEInternalData &fe_data)
411 {
412  ASSERT(fe_type_ == FEVectorContravariant);
413 
414  // shape values
415  if (update_flags & update_values)
416  {
417  for (unsigned int i = 0; i < fe_data.n_points; i++)
418  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
419  {
420  arma::vec fv_vec = elm_values.jacobian(i) * fe_data.ref_shape_values[i][j];
421  for (unsigned int c=0; c<spacedim; c++)
422  shape_values[i][j*spacedim+c] = fv_vec[c];
423  }
424  }
425 
426  // shape gradients
427  if (update_flags & update_gradients)
428  {
429  for (unsigned int i = 0; i < fe_data.n_points; i++)
430  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
431  {
432  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i));
433  for (unsigned int c=0; c<spacedim; c++)
434  shape_gradients[i][j*spacedim+c] = grads.col(c);
435  }
436  }
437 }
438 
439 
440 template<unsigned int spacedim>
441 void FEValues<spacedim>::fill_vec_piola_data(const ElementValues<spacedim> &elm_values,
442  const FEInternalData &fe_data)
443 {
444  ASSERT(fe_type_ == FEVectorPiola);
445 
446  // shape values
447  if (update_flags & update_values)
448  {
449  for (unsigned int i = 0; i < fe_data.n_points; i++)
450  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
451  {
452  arma::vec fv_vec = elm_values.jacobian(i)*fe_data.ref_shape_values[i][j]/elm_values.determinant(i);
453  for (unsigned int c=0; c<spacedim; c++)
454  shape_values[i][j*spacedim+c] = fv_vec(c);
455  }
456  }
457 
458  // shape gradients
459  if (update_flags & update_gradients)
460  {
461  for (unsigned int i = 0; i < fe_data.n_points; i++)
462  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
463  {
464  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j] * trans(elm_values.jacobian(i))
465  / elm_values.determinant(i);
466  for (unsigned int c=0; c<spacedim; c++)
467  shape_gradients[i][j*spacedim+c] = grads.col(c);
468  }
469  }
470 }
471 
472 
473 template<unsigned int spacedim>
474 void FEValues<spacedim>::fill_tensor_data(const ElementValues<spacedim> &elm_values,
475  const FEInternalData &fe_data)
476 {
477  ASSERT(fe_type_ == FETensor);
478 
479  // shape values
480  if (update_flags & update_values)
481  {
482  for (unsigned int i = 0; i < fe_data.n_points; i++)
483  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
484  {
485  arma::vec fv_vec = fe_data.ref_shape_values[i][j];
486  for (unsigned int c=0; c<spacedim*spacedim; c++)
487  shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
488  }
489  }
490 
491  // shape gradients
492  if (update_flags & update_gradients)
493  {
494  for (unsigned int i = 0; i < fe_data.n_points; i++)
495  for (unsigned int j = 0; j < fe_data.n_dofs; j++)
496  {
497  arma::mat grads = trans(elm_values.inverse_jacobian(i)) * fe_data.ref_shape_grads[i][j];
498  for (unsigned int c=0; c<spacedim*spacedim; c++)
499  shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
500  }
501  }
502 }
503 
504 
505 template<unsigned int spacedim>
506 void FEValues<spacedim>::fill_system_data(const ElementValues<spacedim> &elm_values, const FEInternalData &fe_data)
507 {
508  ASSERT(fe_type_ == FEMixedSystem);
509 
510  // for mixed system we first fill data in sub-elements
511  unsigned int comp_offset = 0;
512  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
513  {
514  // fill fe_values for base FE
515  FEInternalData vec_fe_data(fe_data, fe_sys_dofs_[f], comp_offset, fe_sys_n_components_[f]);
516  fe_values_vec[f].fill_data(elm_values, vec_fe_data); // fe_values.fe_values_vec
517 
518  comp_offset += fe_sys_n_components_[f];
519  }
520 
521  // shape values
522  if (update_flags & update_values)
523  {
524  arma::vec fv_vec;
525  unsigned int comp_offset = 0;
526  unsigned int shape_offset = 0;
527  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
528  {
529  // gather fe_values in vectors for FESystem
530  for (unsigned int i=0; i<fe_data.n_points; i++)
531  for (unsigned int n=0; n<fe_sys_dofs_[f].size(); n++)
532  for (unsigned int c=0; c<fe_sys_n_space_components_[f]; c++)
533  shape_values[i][shape_offset+n_components_*n+comp_offset+c] = fe_values_vec[f].shape_values[i][n*fe_sys_n_space_components_[f]+c];
534 
535  comp_offset += fe_sys_n_space_components_[f];
536  shape_offset += fe_sys_dofs_[f].size()*n_components_;
537  }
538  }
539 
540  // shape gradients
541  if (update_flags & update_gradients)
542  {
543  arma::mat grads;
544  unsigned int comp_offset = 0;
545  unsigned int shape_offset = 0;
546  for (unsigned int f=0; f<fe_sys_dofs_.size(); f++)
547  {
548  // gather fe_values in vectors for FESystem
549  for (unsigned int i=0; i<fe_data.n_points; i++)
550  for (unsigned int n=0; n<fe_sys_dofs_[f].size(); n++)
551  for (unsigned int c=0; c<fe_sys_n_space_components_[f]; c++)
552  shape_gradients[i][shape_offset+n_components_*n+comp_offset+c] = fe_values_vec[f].shape_gradients[i][n*fe_sys_n_space_components_[f]+c];
553 
554  comp_offset += fe_sys_n_space_components_[f];
555  shape_offset += fe_sys_dofs_[f].size()*n_components_;
556  }
557  }
558 
559 }*/
560 
561 
562 template<unsigned int spacedim>
564 {
565  ASSERT_EQ( this->dim_, cell.dim() );
566 
567  if (!elm_values_->cell().is_valid() ||
568  elm_values_->cell() != cell)
569  {
570  elm_values_->reinit(cell);
571  }
572 
573  this->fill_data(*elm_values_, *this->fe_data_);
574 }
575 
576 
577 template<unsigned int spacedim>
578 void FEValues<spacedim>::reinit(const Side &cell_side)
579 {
580  ASSERT_EQ( this->dim_, cell_side.dim()+1 );
581 
582  if (!elm_values_->side().is_valid() ||
583  elm_values_->side() != cell_side)
584  {
585  elm_values_->reinit(cell_side);
586  }
587 
588  const LongIdx sid = cell_side.side_idx();
589 
590  // calculation of finite element data
591  this->fill_data(*elm_values_, *(this->side_fe_data_[sid]) );
592 }
593 
594 
595 
596 template<unsigned int spacedim>
597 PatchFEValues<spacedim>::PatchFEValues(unsigned int max_size)
598 : FEValuesBase<PatchFEValues<spacedim>, spacedim>(),
599  patch_data_idx_(-1), used_size_(0), max_n_elem_(max_size) {}
600 
601 
602 template<unsigned int spacedim>
604  element_patch_map_.clear();
605  if (object_type_ == ElementFE)
606  used_size_ = patch_elements.size();
607  else
608  used_size_ = patch_elements.size() * (this->dim_+1);
610 
611  unsigned int i=0;
612  for (auto it=patch_elements.begin(); it!=patch_elements.end(); ++it, ++i) {
613  if (object_type_ == ElementFE) {
614  patch_data_idx_ = i;
615  element_patch_map_[it->second] = i;
616  element_data_[i].elm_values_->reinit(it->first);
617  this->fill_data(*element_data_[i].elm_values_, *this->fe_data_);
618  } else {
619  element_patch_map_[it->second] = i * (this->dim_+1);
620  for (unsigned int sid=0; sid<this->dim_+1; ++sid) {
621  patch_data_idx_ = i * (this->dim_+1) + sid;
622  element_data_[patch_data_idx_].elm_values_->reinit( *it->first.side(sid) );
623  this->fill_data(*element_data_[patch_data_idx_].elm_values_, *this->side_fe_data_[sid]);
624 
625  }
626  }
627  }
628 }
629 
630 
631 template<unsigned int spacedim>
632 void PatchFEValues<spacedim>::allocate_in(unsigned int q_dim)
633 {
635 
636  if ( q_dim == this->dim_ ) {
637  element_data_.resize( this->max_n_elem_ );
639  } else if ( q_dim+1 == this->dim_ ) {
640  element_data_.resize( this->max_n_elem_ * (this->dim_+1) );
642  } else
643  ASSERT(false)(q_dim)(this->dim_).error("Invalid dimension of quadrature!");
644 
645  for (uint i=0; i<max_size(); ++i) {
646  if (this->update_flags & update_values)
647  element_data_[i].shape_values_.resize(this->n_points_, vector<double>(this->n_dofs_*this->n_components_));
648 
649  if (this->update_flags & update_gradients)
650  element_data_[i].shape_gradients_.resize(this->n_points_, vector<arma::vec::fixed<spacedim> >(this->n_dofs_*this->n_components_));
651  }
652 
653  this->fv_ = this;
654 }
655 
656 
657 template<unsigned int spacedim>
659  Quadrature &q,
660  unsigned int dim)
661 {
662  for (uint i=0; i<max_size(); ++i)
663  element_data_[i].elm_values_ = std::make_shared<ElementValues<spacedim> >(q, this->update_flags, dim);
664 }
665 
666 
667 template<unsigned int spacedim>
669 {
670  for (unsigned int i=0; i<this->fe_values_vec.size(); ++i)
671  this->fe_values_vec[i].resize( this->max_size() );
672 }
673 
674 
675 template<unsigned int spacedim>
676 arma::vec::fixed<spacedim> PatchFEValues<spacedim>::shape_grad_component(const unsigned int function_no,
677  const unsigned int point_no,
678  const unsigned int comp) const
679 {
680  ASSERT_LT(function_no, this->n_dofs_);
681  ASSERT_LT(point_no, this->n_points_);
682  ASSERT_LT(comp, this->n_components_);
683  return element_data_[patch_data_idx_].shape_gradients_[point_no][function_no*this->n_components_+comp];
684 }
685 
686 
687 
689  QGauss::array &quadrature,
691  UpdateFlags flags)
692 {
694  fv[0].initialize(quadrature[0], *fe[0_d], flags);
695  fv[1].initialize(quadrature[1], *fe[1_d], flags);
696  fv[2].initialize(quadrature[2], *fe[2_d], flags);
697  fv[3].initialize(quadrature[3], *fe[3_d], flags);
698  return fv;
699 }
700 
701 
702 
703 // explicit instantiation
704 template void FEValuesBase<FEValues<3>, 3>::initialize<0>(Quadrature&, FiniteElement<0>&, UpdateFlags);
705 template void FEValuesBase<FEValues<3>, 3>::initialize<1>(Quadrature&, FiniteElement<1>&, UpdateFlags);
706 template void FEValuesBase<FEValues<3>, 3>::initialize<2>(Quadrature&, FiniteElement<2>&, UpdateFlags);
707 template void FEValuesBase<FEValues<3>, 3>::initialize<3>(Quadrature&, FiniteElement<3>&, UpdateFlags);
708 
709 template void FEValuesBase<PatchFEValues<3>, 3>::initialize<0>(Quadrature&, FiniteElement<0>&, UpdateFlags);
710 template void FEValuesBase<PatchFEValues<3>, 3>::initialize<1>(Quadrature&, FiniteElement<1>&, UpdateFlags);
711 template void FEValuesBase<PatchFEValues<3>, 3>::initialize<2>(Quadrature&, FiniteElement<2>&, UpdateFlags);
712 template void FEValuesBase<PatchFEValues<3>, 3>::initialize<3>(Quadrature&, FiniteElement<3>&, UpdateFlags);
713 
714 template void FEValuesBase<FEValues<3>, 3>::fill_data_specialized<MapScalar<FEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<FEValues<3>, 3>::FEInternalData &);
715 template void FEValuesBase<FEValues<3>, 3>::fill_data_specialized<MapPiola<FEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<FEValues<3>, 3>::FEInternalData &);
716 template void FEValuesBase<FEValues<3>, 3>::fill_data_specialized<MapContravariant<FEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<FEValues<3>, 3>::FEInternalData &);
717 template void FEValuesBase<FEValues<3>, 3>::fill_data_specialized<MapVector<FEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<FEValues<3>, 3>::FEInternalData &);
718 template void FEValuesBase<FEValues<3>, 3>::fill_data_specialized<MapTensor<FEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<FEValues<3>, 3>::FEInternalData &);
719 template void FEValuesBase<FEValues<3>, 3>::fill_data_specialized<MapSystem<FEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<FEValues<3>, 3>::FEInternalData &);
720 
721 template void FEValuesBase<PatchFEValues<3>, 3>::fill_data_specialized<MapScalar<PatchFEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<PatchFEValues<3>, 3>::FEInternalData &);
722 template void FEValuesBase<PatchFEValues<3>, 3>::fill_data_specialized<MapPiola<PatchFEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<PatchFEValues<3>, 3>::FEInternalData &);
723 template void FEValuesBase<PatchFEValues<3>, 3>::fill_data_specialized<MapContravariant<PatchFEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<PatchFEValues<3>, 3>::FEInternalData &);
724 template void FEValuesBase<PatchFEValues<3>, 3>::fill_data_specialized<MapVector<PatchFEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<PatchFEValues<3>, 3>::FEInternalData &);
725 template void FEValuesBase<PatchFEValues<3>, 3>::fill_data_specialized<MapTensor<PatchFEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<PatchFEValues<3>, 3>::FEInternalData &);
726 template void FEValuesBase<PatchFEValues<3>, 3>::fill_data_specialized<MapSystem<PatchFEValues<3>, 3>>(const ElementValues<3> &, const typename FEValuesBase<PatchFEValues<3>, 3>::FEInternalData &);
727 
728 
729 
730 
731 
732 
733 
734 
735 
736 
737 
738 
739 
740 
741 
742 template class FEValues<3>;
743 template class PatchFEValues<3>;
#define ASSERT(expr)
Definition: asserts.hh:351
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:348
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Definition: asserts.hh:301
#define ASSERT_PERMANENT_GT(a, b)
Definition of comparative assert macro (Greater Than)
Definition: asserts.hh:313
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Definition: asserts.hh:333
#define ASSERT_LE(a, b)
Definition of comparative assert macro (Less or Equal) only for debug mode.
Definition: asserts.hh:309
unsigned int dim() const
Definition: accessors.hh:188
Class for computation of data on cell and side.
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:102
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
Definition: fe_system.hh:147
std::vector< unsigned int > fe_dofs(unsigned int fe_index)
Return dof indices belonging to given sub-FE.
Definition: fe_system.cc:267
Structure for storing the precomputed finite element data.
Definition: fe_values.hh:154
FEInternalData(unsigned int np, unsigned int nd)
Definition: fe_values.cc:49
unsigned int n_points
Number of quadrature points.
Definition: fe_values.hh:184
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Definition: fe_values.hh:172
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Definition: fe_values.hh:181
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Definition: fe_values.hh:244
void allocate(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
Definition: fe_values.cc:176
unsigned int n_points_
Number of integration points.
Definition: fe_values.hh:226
FV * fv_
Helper object, we need its for ViewsCache initialization.
Definition: fe_values.hh:262
std::vector< shared_ptr< FEInternalData > > side_fe_data_
Precomputed FE data (shape functions on reference element) for all side quadrature points.
Definition: fe_values.hh:259
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
Definition: fe_values.hh:241
FEType fe_type_
Type of finite element (scalar, vector, tensor).
Definition: fe_values.hh:232
unsigned int n_dofs_
Number of finite element dofs.
Definition: fe_values.hh:229
unsigned int n_points() const
Returns the number of quadrature points.
Definition: fe_values.hh:107
virtual void init_fe_val_vec()=0
Initialize fe_values_vec only in PatchFEValues.
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Definition: fe_values.hh:253
virtual void allocate_in(unsigned int)=0
Initialize vectors declared separately in descendants.
void fill_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure.
Definition: fe_values.cc:257
std::vector< FV > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
Definition: fe_values.hh:247
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
Definition: fe_values.hh:238
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
Definition: fe_values.hh:235
std::shared_ptr< typename FEValuesBase< FV, spacedim >::FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
Definition: fe_values.cc:221
unsigned int n_components_
Number of components of the FE.
Definition: fe_values.hh:250
FEValuesBase()
Default constructor with postponed initialization.
Definition: fe_values.cc:42
void fill_data_specialized(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure....
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
Definition: fe_values.cc:132
virtual void initialize_in(Quadrature &, unsigned int)=0
Initialize ElementValues separately in descendants.
std::shared_ptr< FEInternalData > fe_data_
Precomputed finite element data.
Definition: fe_values.hh:256
unsigned int dim_
Dimension of reference space.
Definition: fe_values.hh:223
Calculates finite element data on the actual cell.
Definition: fe_values.hh:288
~FEValues()
Correct deallocation of objects created by 'initialize' methods.
Definition: fe_values.cc:317
FEValues()
Default constructor with postponed initialization.
Definition: fe_values.cc:312
arma::vec::fixed< spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
Definition: fe_values.cc:346
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Definition: fe_values.cc:563
void initialize_in(Quadrature &q, unsigned int dim) override
Implement FEValuesBase::initialize_in.
Definition: fe_values.cc:322
void allocate_in(unsigned int q_dim) override
Implement FEValuesBase::allocate_in.
Definition: fe_values.cc:332
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
unsigned int n_space_components(unsigned int spacedim)
Number of components of FE in a mapped space with dimension spacedim.
virtual UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
FEType type_
Type of FiniteElement.
double shape_value(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the value of the comp-th component of the i-th shape function at the point p on the refere...
arma::vec::fixed< dim > shape_grad(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the comp-th component of the gradient of the i-th shape function at the point p on the ref...
unsigned int n_components() const
Returns numer of components of the basis function.
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:30
unsigned int patch_data_idx_
Patch index of processed element / side.
Definition: fe_values.hh:680
unsigned int max_size() const
Definition: fe_values.hh:510
void reinit(PatchElementsList patch_elements)
Reinit data.
PatchFEValues(unsigned int max_size=0)
Constructor, set maximal number of elements on patch.
unsigned int max_n_elem_
Maximal number of elements on patch.
Definition: fe_values.hh:692
std::vector< ElementFEData > element_data_
Data of elements / sides on patch.
Definition: fe_values.hh:686
unsigned int used_size_
Number of elements / sides on patch. Must be less or equal to size of element_data vector.
Definition: fe_values.hh:689
void init_fe_val_vec() override
Implement FEValuesBase::initialize_in.
MeshObjectType object_type_
Distinguishes using of PatchFEValues for storing data of elements or sides.
Definition: fe_values.hh:695
arma::vec::fixed< spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
void allocate_in(unsigned int q_dim) override
Implement FEValuesBase::allocate_in.
void initialize_in(Quadrature &q, unsigned int dim) override
Implement FEValuesBase::initialize_in.
std::map< unsigned int, unsigned int > element_patch_map_
Map of element patch indexes to element_data_.
Definition: fe_values.hh:683
void resize(unsigned int max_size)
Set size of ElementFEData. Important: Use only during the initialization of FESystem !
Definition: fe_values.hh:666
std::array< QGauss, 4 > array
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: quadrature.hh:48
Quadrature make_from_side(unsigned int sid) const
Definition: quadrature.cc:47
unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:86
unsigned int dim() const
Definition: quadrature.hh:72
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:151
unsigned int side_idx() const
Returns local index of the side on the element.
Definition: accessors.hh:444
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class FESystem for compound finite elements.
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
Definition: fe_values.cc:688
Class FEValues calculates finite element data on the actual cells such as shape function values,...
Abstract class for description of finite elements.
@ FEScalar
@ FEMixedSystem
@ FEVectorPiola
@ FETensor
@ FEVectorContravariant
@ FEVector
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
unsigned int uint
ArmaMat< double, N, M > mat
Definition: armor.hh:936
Class FEValues calculates finite element data on the actual cells such as shape function values,...
#define FMT_UNUSED
Definition: posix.h:75
Basic definitions of numerical quadrature rules.
void initialize(const FV &fv, const FiniteElement< DIM > &fe)
Definition: fe_values.cc:77
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
Definition: update_flags.hh:68
@ update_values
Shape function values.
Definition: update_flags.hh:87
@ update_gradients
Shape function gradients.
Definition: update_flags.hh:95