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