Reference documentation for deal.II version 8.5.1
fe_system.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/memory_consumption.h>
17 #include <deal.II/base/quadrature.h>
18 #include <deal.II/base/qprojector.h>
19 #include <deal.II/grid/tria.h>
20 #include <deal.II/grid/tria_iterator.h>
21 #include <deal.II/dofs/dof_accessor.h>
22 #include <deal.II/fe/mapping.h>
23 #include <deal.II/fe/fe_system.h>
24 #include <deal.II/fe/fe_values.h>
25 #include <deal.II/fe/fe_tools.h>
26 
27 #include <sstream>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 namespace
33 {
34  bool IsNonZero (unsigned int i)
35  {
36  return i>0;
37  }
38 
39  unsigned int count_nonzeros(const std::vector<unsigned int> &vec)
40  {
41  return std::count_if(vec.begin(), vec.end(), IsNonZero);
42  }
43 
44 }
45 /* ----------------------- FESystem::InternalData ------------------- */
46 
47 
48 template <int dim, int spacedim>
49 FESystem<dim,spacedim>::InternalData::InternalData(const unsigned int n_base_elements)
50  :
51  base_fe_datas(n_base_elements),
52  base_fe_output_objects(n_base_elements)
53 {}
54 
55 
56 
57 template <int dim, int spacedim>
59 {
60  // delete pointers and set them to zero to avoid inadvertent use
61  for (unsigned int i=0; i<base_fe_datas.size(); ++i)
62  if (base_fe_datas[i])
63  {
64  delete base_fe_datas[i];
65  base_fe_datas[i] = 0;
66  }
67 }
68 
69 
70 
71 template <int dim, int spacedim>
74 InternalData::get_fe_data (const unsigned int base_no) const
75 {
76  Assert(base_no<base_fe_datas.size(),
77  ExcIndexRange(base_no,0,base_fe_datas.size()));
78  return *base_fe_datas[base_no];
79 }
80 
81 
82 
83 template <int dim, int spacedim>
84 void
86 InternalData::set_fe_data (const unsigned int base_no,
88 {
89  Assert(base_no<base_fe_datas.size(),
90  ExcIndexRange(base_no,0,base_fe_datas.size()));
91  base_fe_datas[base_no]=ptr;
92 }
93 
94 
95 
96 template <int dim, int spacedim>
99 InternalData::get_fe_output_object (const unsigned int base_no) const
100 {
101  Assert(base_no<base_fe_output_objects.size(),
102  ExcIndexRange(base_no,0,base_fe_output_objects.size()));
103  return base_fe_output_objects[base_no];
104 }
105 
106 
107 
108 /* ---------------------------------- FESystem ------------------- */
109 
110 
111 template <int dim, int spacedim>
113 
114 
115 template <int dim, int spacedim>
117  const unsigned int n_elements) :
118  FiniteElement<dim,spacedim> (FETools::Compositing::multiply_dof_numbers(&fe, n_elements),
119  FETools::Compositing::compute_restriction_is_additive_flags (&fe, n_elements),
120  FETools::Compositing::compute_nonzero_components(&fe, n_elements)),
121  base_elements((n_elements>0))
122 {
123  std::vector<const FiniteElement<dim,spacedim>*> fes;
124  fes.push_back(&fe);
125  std::vector<unsigned int> multiplicities;
126  multiplicities.push_back(n_elements);
127  initialize(fes, multiplicities);
128 }
129 
130 
131 
132 template <int dim, int spacedim>
134  const unsigned int n1,
135  const FiniteElement<dim,spacedim> &fe2,
136  const unsigned int n2) :
137  FiniteElement<dim,spacedim> (FETools::Compositing::multiply_dof_numbers(&fe1, n1, &fe2, n2),
138  FETools::Compositing::compute_restriction_is_additive_flags (&fe1, n1,
139  &fe2, n2),
140  FETools::Compositing::compute_nonzero_components(&fe1, n1,
141  &fe2, n2)),
142  base_elements((n1>0)+(n2>0))
143 {
144  std::vector<const FiniteElement<dim,spacedim>*> fes;
145  fes.push_back(&fe1);
146  fes.push_back(&fe2);
147  std::vector<unsigned int> multiplicities;
148  multiplicities.push_back(n1);
149  multiplicities.push_back(n2);
150  initialize(fes, multiplicities);
151 }
152 
153 
154 
155 template <int dim, int spacedim>
157  const unsigned int n1,
158  const FiniteElement<dim,spacedim> &fe2,
159  const unsigned int n2,
160  const FiniteElement<dim,spacedim> &fe3,
161  const unsigned int n3) :
162  FiniteElement<dim,spacedim> (FETools::Compositing::multiply_dof_numbers(&fe1, n1,
163  &fe2, n2,
164  &fe3, n3),
165  FETools::Compositing::compute_restriction_is_additive_flags (&fe1, n1,
166  &fe2, n2,
167  &fe3, n3),
168  FETools::Compositing::compute_nonzero_components(&fe1, n1,
169  &fe2, n2,
170  &fe3, n3)),
171  base_elements((n1>0)+(n2>0)+(n3>0))
172 {
173  std::vector<const FiniteElement<dim,spacedim>*> fes;
174  fes.push_back(&fe1);
175  fes.push_back(&fe2);
176  fes.push_back(&fe3);
177  std::vector<unsigned int> multiplicities;
178  multiplicities.push_back(n1);
179  multiplicities.push_back(n2);
180  multiplicities.push_back(n3);
181  initialize(fes, multiplicities);
182 }
183 
184 
185 
186 template <int dim, int spacedim>
188  const unsigned int n1,
189  const FiniteElement<dim,spacedim> &fe2,
190  const unsigned int n2,
191  const FiniteElement<dim,spacedim> &fe3,
192  const unsigned int n3,
193  const FiniteElement<dim,spacedim> &fe4,
194  const unsigned int n4) :
195  FiniteElement<dim,spacedim> (FETools::Compositing::multiply_dof_numbers(&fe1, n1,
196  &fe2, n2,
197  &fe3, n3,
198  &fe4, n4),
199  FETools::Compositing::compute_restriction_is_additive_flags (&fe1, n1,
200  &fe2, n2,
201  &fe3, n3,
202  &fe4, n4),
203  FETools::Compositing::compute_nonzero_components(&fe1, n1,
204  &fe2, n2,
205  &fe3, n3,
206  &fe4 ,n4)),
207  base_elements((n1>0)+(n2>0)+(n3>0)+(n4>0))
208 {
209  std::vector<const FiniteElement<dim,spacedim>*> fes;
210  fes.push_back(&fe1);
211  fes.push_back(&fe2);
212  fes.push_back(&fe3);
213  fes.push_back(&fe4);
214  std::vector<unsigned int> multiplicities;
215  multiplicities.push_back(n1);
216  multiplicities.push_back(n2);
217  multiplicities.push_back(n3);
218  multiplicities.push_back(n4);
219  initialize(fes, multiplicities);
220 }
221 
222 
223 
224 template <int dim, int spacedim>
226  const unsigned int n1,
227  const FiniteElement<dim,spacedim> &fe2,
228  const unsigned int n2,
229  const FiniteElement<dim,spacedim> &fe3,
230  const unsigned int n3,
231  const FiniteElement<dim,spacedim> &fe4,
232  const unsigned int n4,
233  const FiniteElement<dim,spacedim> &fe5,
234  const unsigned int n5) :
235  FiniteElement<dim,spacedim> (FETools::Compositing::multiply_dof_numbers(&fe1, n1,
236  &fe2, n2,
237  &fe3, n3,
238  &fe4, n4,
239  &fe5, n5),
240  FETools::Compositing::compute_restriction_is_additive_flags (&fe1, n1,
241  &fe2, n2,
242  &fe3, n3,
243  &fe4, n4,
244  &fe5, n5),
245  FETools::Compositing::compute_nonzero_components(&fe1, n1,
246  &fe2, n2,
247  &fe3, n3,
248  &fe4 ,n4,
249  &fe5, n5)),
250  base_elements((n1>0)+(n2>0)+(n3>0)+(n4>0)+(n5>0))
251 {
252  std::vector<const FiniteElement<dim,spacedim>*> fes;
253  fes.push_back(&fe1);
254  fes.push_back(&fe2);
255  fes.push_back(&fe3);
256  fes.push_back(&fe4);
257  fes.push_back(&fe5);
258  std::vector<unsigned int> multiplicities;
259  multiplicities.push_back(n1);
260  multiplicities.push_back(n2);
261  multiplicities.push_back(n3);
262  multiplicities.push_back(n4);
263  multiplicities.push_back(n5);
264  initialize(fes, multiplicities);
265 }
266 
267 
268 
269 template <int dim, int spacedim>
271  const std::vector<const FiniteElement<dim,spacedim>*> &fes,
272  const std::vector<unsigned int> &multiplicities)
273  :
274  FiniteElement<dim,spacedim> (FETools::Compositing::multiply_dof_numbers(fes, multiplicities),
275  FETools::Compositing::compute_restriction_is_additive_flags (fes, multiplicities),
276  FETools::Compositing::compute_nonzero_components(fes, multiplicities)),
277  base_elements(count_nonzeros(multiplicities))
278 {
279  initialize(fes, multiplicities);
280 }
281 
282 
283 
284 template <int dim, int spacedim>
286 {}
287 
288 
289 
290 template <int dim, int spacedim>
291 std::string
293 {
294  // note that the
295  // FETools::get_fe_by_name
296  // function depends on the
297  // particular format of the string
298  // this function returns, so they
299  // have to be kept in synch
300 
301  std::ostringstream namebuf;
302 
303  namebuf << "FESystem<"
304  << Utilities::dim_string(dim,spacedim)
305  << ">[";
306  for (unsigned int i=0; i< this->n_base_elements(); ++i)
307  {
308  namebuf << base_element(i).get_name();
309  if (this->element_multiplicity(i) != 1)
310  namebuf << '^' << this->element_multiplicity(i);
311  if (i != this->n_base_elements()-1)
312  namebuf << '-';
313  }
314  namebuf << ']';
315 
316  return namebuf.str();
317 }
318 
319 
320 
321 template <int dim, int spacedim>
324 {
325  std::vector< const FiniteElement<dim,spacedim>* > fes;
326  std::vector<unsigned int> multiplicities;
327 
328  for (unsigned int i=0; i<this->n_base_elements(); i++)
329  {
330  fes.push_back( & base_element(i) );
331  multiplicities.push_back(this->element_multiplicity(i) );
332  }
333  return new FESystem<dim,spacedim>(fes, multiplicities);
334 }
335 
336 
337 
338 template <int dim, int spacedim>
339 double
341  const Point<dim> &p) const
342 {
343  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
344  Assert (this->is_primitive(i),
346 
347  return (base_element(this->system_to_base_table[i].first.first)
348  .shape_value(this->system_to_base_table[i].second, p));
349 }
350 
351 
352 
353 template <int dim, int spacedim>
354 double
356  const Point<dim> &p,
357  const unsigned int component) const
358 {
359  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
360  Assert (component < this->n_components(),
361  ExcIndexRange (component, 0, this->n_components()));
362 
363  // if this value is supposed to be
364  // zero, then return right away...
365  if (this->nonzero_components[i][component] == false)
366  return 0;
367 
368  // ...otherwise: first find out to
369  // which of the base elements this
370  // desired component belongs, and
371  // which component within this base
372  // element it is
373  const unsigned int base = this->component_to_base_index(component).first;
374  const unsigned int component_in_base = this->component_to_base_index(component).second;
375 
376  // then get value from base
377  // element. note that that will
378  // throw an error should the
379  // respective shape function not be
380  // primitive; thus, there is no
381  // need to check this here
382  return (base_element(base).
383  shape_value_component(this->system_to_base_table[i].second,
384  p,
385  component_in_base));
386 }
387 
388 
389 
390 template <int dim, int spacedim>
393  const Point<dim> &p) const
394 {
395  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
396  Assert (this->is_primitive(i),
398 
399  return (base_element(this->system_to_base_table[i].first.first)
400  .shape_grad(this->system_to_base_table[i].second, p));
401 }
402 
403 
404 
405 template <int dim, int spacedim>
408  const Point<dim> &p,
409  const unsigned int component) const
410 {
411  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
412  Assert (component < this->n_components(),
413  ExcIndexRange (component, 0, this->n_components()));
414 
415  // if this value is supposed to be zero, then return right away...
416  if (this->nonzero_components[i][component] == false)
417  return Tensor<1,dim>();
418 
419  // ...otherwise: first find out to which of the base elements this desired
420  // component belongs, and which component within this base element it is
421  const unsigned int base = this->component_to_base_index(component).first;
422  const unsigned int component_in_base = this->component_to_base_index(component).second;
423 
424  // then get value from base element. note that that will throw an error
425  // should the respective shape function not be primitive; thus, there is no
426  // need to check this here
427  return (base_element(base).
428  shape_grad_component(this->system_to_base_table[i].second,
429  p,
430  component_in_base));
431 }
432 
433 
434 
435 template <int dim, int spacedim>
438  const Point<dim> &p) const
439 {
440  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
441  Assert (this->is_primitive(i),
443 
444  return (base_element(this->system_to_base_table[i].first.first)
445  .shape_grad_grad(this->system_to_base_table[i].second, p));
446 }
447 
448 
449 
450 template <int dim, int spacedim>
453  const Point<dim> &p,
454  const unsigned int component) const
455 {
456  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
457  Assert (component < this->n_components(),
458  ExcIndexRange (component, 0, this->n_components()));
459 
460  // if this value is supposed to be zero, then return right away...
461  if (this->nonzero_components[i][component] == false)
462  return Tensor<2,dim>();
463 
464  // ...otherwise: first find out to which of the base elements this desired
465  // component belongs, and which component within this base element it is
466  const unsigned int base = this->component_to_base_index(component).first;
467  const unsigned int component_in_base = this->component_to_base_index(component).second;
468 
469  // then get value from base element. note that that will throw an error
470  // should the respective shape function not be primitive; thus, there is no
471  // need to check this here
472  return (base_element(base).
473  shape_grad_grad_component(this->system_to_base_table[i].second,
474  p,
475  component_in_base));
476 }
477 
478 
479 
480 template <int dim, int spacedim>
483  const Point<dim> &p) const
484 {
485  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
486  Assert (this->is_primitive(i),
488 
489  return (base_element(this->system_to_base_table[i].first.first)
490  .shape_3rd_derivative(this->system_to_base_table[i].second, p));
491 }
492 
493 
494 
495 template <int dim, int spacedim>
498  const Point<dim> &p,
499  const unsigned int component) const
500 {
501  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
502  Assert (component < this->n_components(),
503  ExcIndexRange (component, 0, this->n_components()));
504 
505  // if this value is supposed to be zero, then return right away...
506  if (this->nonzero_components[i][component] == false)
507  return Tensor<3,dim>();
508 
509  // ...otherwise: first find out to which of the base elements this desired
510  // component belongs, and which component within this base element it is
511  const unsigned int base = this->component_to_base_index(component).first;
512  const unsigned int component_in_base = this->component_to_base_index(component).second;
513 
514  // then get value from base element. note that that will throw an error
515  // should the respective shape function not be primitive; thus, there is no
516  // need to check this here
517  return (base_element(base).
518  shape_3rd_derivative_component(this->system_to_base_table[i].second,
519  p,
520  component_in_base));
521 }
522 
523 
524 
525 template <int dim, int spacedim>
528  const Point<dim> &p) const
529 {
530  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
531  Assert (this->is_primitive(i),
533 
534  return (base_element(this->system_to_base_table[i].first.first)
535  .shape_4th_derivative(this->system_to_base_table[i].second, p));
536 }
537 
538 
539 
540 template <int dim, int spacedim>
543  const Point<dim> &p,
544  const unsigned int component) const
545 {
546  Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
547  Assert (component < this->n_components(),
548  ExcIndexRange (component, 0, this->n_components()));
549 
550  // if this value is supposed to be zero, then return right away...
551  if (this->nonzero_components[i][component] == false)
552  return Tensor<4,dim>();
553 
554  // ...otherwise: first find out to which of the base elements this desired
555  // component belongs, and which component within this base element it is
556  const unsigned int base = this->component_to_base_index(component).first;
557  const unsigned int component_in_base = this->component_to_base_index(component).second;
558 
559  // then get value from base element. note that that will throw an error
560  // should the respective shape function not be primitive; thus, there is no
561  // need to check this here
562  return (base_element(base).
563  shape_4th_derivative_component(this->system_to_base_table[i].second,
564  p,
565  component_in_base));
566 }
567 
568 
569 
570 template <int dim, int spacedim>
571 void
573  const FiniteElement<dim,spacedim> &x_source_fe,
574  FullMatrix<double> &interpolation_matrix) const
575 {
576  // check that the size of the matrices is correct. for historical
577  // reasons, if you call matrix.reinit(8,0), it sets the sizes
578  // to m==n==0 internally. this may happen when we use a FE_Nothing,
579  // so write the test in a more lenient way
580  Assert ((interpolation_matrix.m() == this->dofs_per_cell)
581  ||
582  (x_source_fe.dofs_per_cell == 0),
583  ExcDimensionMismatch (interpolation_matrix.m(),
584  this->dofs_per_cell));
585  Assert ((interpolation_matrix.n() == x_source_fe.dofs_per_cell)
586  ||
587  (this->dofs_per_cell == 0),
588  ExcDimensionMismatch (interpolation_matrix.m(),
589  x_source_fe.dofs_per_cell));
590 
591  // there are certain conditions that the two elements have to satisfy so
592  // that this can work.
593  //
594  // condition 1: the other element must also be a system element
595 
596  typedef FiniteElement<dim,spacedim> FEL;
597  AssertThrow ((x_source_fe.get_name().find ("FESystem<") == 0)
598  ||
599  (dynamic_cast<const FESystem<dim,spacedim>*>(&x_source_fe) != 0),
600  typename FEL::
601  ExcInterpolationNotImplemented());
602 
603  // ok, source is a system element, so we may be able to do the work
604  const FESystem<dim,spacedim> &source_fe
605  = dynamic_cast<const FESystem<dim,spacedim>&>(x_source_fe);
606 
607  // condition 2: same number of basis elements
608  AssertThrow (this->n_base_elements() == source_fe.n_base_elements(),
609  typename FEL::
610  ExcInterpolationNotImplemented());
611 
612  // condition 3: same number of basis elements
613  for (unsigned int i=0; i< this->n_base_elements(); ++i)
614  AssertThrow (this->element_multiplicity(i) ==
615  source_fe.element_multiplicity(i),
616  typename FEL::
617  ExcInterpolationNotImplemented());
618 
619  // ok, so let's try whether it works:
620 
621  // first let's see whether all the basis elements actually generate their
622  // interpolation matrices. if we get past the following loop, then
623  // apparently none of the called base elements threw an exception, so we're
624  // fine continuing and assembling the one big matrix from the small ones of
625  // the base elements
626  std::vector<FullMatrix<double> > base_matrices (this->n_base_elements());
627  for (unsigned int i=0; i<this->n_base_elements(); ++i)
628  {
629  base_matrices[i].reinit (base_element(i).dofs_per_cell,
630  source_fe.base_element(i).dofs_per_cell);
631  base_element(i).get_interpolation_matrix (source_fe.base_element(i),
632  base_matrices[i]);
633  }
634 
635  // first clear big matrix, to make sure that entries that would couple
636  // different bases (or multiplicity indices) are really zero. then assign
637  // entries
638  interpolation_matrix = 0;
639  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
640  for (unsigned int j=0; j<source_fe.dofs_per_cell; ++j)
641  if (this->system_to_base_table[i].first ==
642  source_fe.system_to_base_table[j].first)
643  interpolation_matrix(i,j)
644  = (base_matrices[this->system_to_base_table[i].first.first]
645  (this->system_to_base_table[i].second,
646  source_fe.system_to_base_table[j].second));
647 }
648 
649 
650 
651 template <int dim, int spacedim>
652 const FullMatrix<double> &
654 ::get_restriction_matrix (const unsigned int child,
655  const RefinementCase<dim> &refinement_case) const
656 {
659  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
660  ExcMessage("Restriction matrices are only available for refined cells!"));
661  Assert (child<GeometryInfo<dim>::n_children(refinement_case),
662  ExcIndexRange(child,0,GeometryInfo<dim>::n_children(refinement_case)));
663 
664  // initialization upon first request
665  if (this->restriction[refinement_case-1][child].n() == 0)
666  {
667  Threads::Mutex::ScopedLock lock(this->mutex);
668 
669  // check if updated while waiting for lock
670  if (this->restriction[refinement_case-1][child].n() ==
671  this->dofs_per_cell)
672  return this->restriction[refinement_case-1][child];
673 
674  // Check if some of the matrices of the base elements are void.
675  bool do_restriction = true;
676 
677  // shortcut for accessing local restrictions further down
678  std::vector<const FullMatrix<double> *>
679  base_matrices(this->n_base_elements());
680 
681  for (unsigned int i=0; i<this->n_base_elements(); ++i)
682  {
683  base_matrices[i] =
684  &base_element(i).get_restriction_matrix(child, refinement_case);
685  if (base_matrices[i]->n() != base_element(i).dofs_per_cell)
686  do_restriction = false;
687  }
688  Assert(do_restriction,
690 
691  // if we did not encounter void matrices, initialize the matrix sizes
692  if (do_restriction)
693  {
694  FullMatrix<double> restriction(this->dofs_per_cell,
695  this->dofs_per_cell);
696 
697  // distribute the matrices of the base finite elements to the
698  // matrices of this object. for this, loop over all degrees of
699  // freedom and take the respective entry of the underlying base
700  // element.
701  //
702  // note that we by definition of a base element, they are
703  // independent, i.e. do not couple. only DoFs that belong to the
704  // same instance of a base element may couple
705  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
706  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
707  {
708  // first find out to which base element indices i and j
709  // belong, and which instance thereof in case the base element
710  // has a multiplicity greater than one. if they should not
711  // happen to belong to the same instance of a base element,
712  // then they cannot couple, so go on with the next index
713  if (this->system_to_base_table[i].first !=
714  this->system_to_base_table[j].first)
715  continue;
716 
717  // so get the common base element and the indices therein:
718  const unsigned int
719  base = this->system_to_base_table[i].first.first;
720 
721  const unsigned int
722  base_index_i = this->system_to_base_table[i].second,
723  base_index_j = this->system_to_base_table[j].second;
724 
725  // if we are sure that DoFs i and j may couple, then copy
726  // entries of the matrices:
727  restriction(i,j) = (*base_matrices[base])(base_index_i,base_index_j);
728  }
729 
730  restriction.swap(const_cast<FullMatrix<double> &>
731  (this->restriction[refinement_case-1][child]));
732  }
733  }
734 
735  return this->restriction[refinement_case-1][child];
736 }
737 
738 
739 
740 template <int dim, int spacedim>
741 const FullMatrix<double> &
743 ::get_prolongation_matrix (const unsigned int child,
744  const RefinementCase<dim> &refinement_case) const
745 {
748  Assert (refinement_case!=RefinementCase<dim>::no_refinement,
749  ExcMessage("Restriction matrices are only available for refined cells!"));
750  Assert (child<GeometryInfo<dim>::n_children(refinement_case),
751  ExcIndexRange(child,0,GeometryInfo<dim>::n_children(refinement_case)));
752 
753  // initialization upon first request, construction completely analogous to
754  // restriction matrix
755  if (this->prolongation[refinement_case-1][child].n() == 0)
756  {
757  Threads::Mutex::ScopedLock lock(this->mutex);
758 
759  if (this->prolongation[refinement_case-1][child].n() ==
760  this->dofs_per_cell)
761  return this->prolongation[refinement_case-1][child];
762 
763  bool do_prolongation = true;
764  std::vector<const FullMatrix<double> *>
765  base_matrices(this->n_base_elements());
766  for (unsigned int i=0; i<this->n_base_elements(); ++i)
767  {
768  base_matrices[i] =
769  &base_element(i).get_prolongation_matrix(child, refinement_case);
770  if (base_matrices[i]->n() != base_element(i).dofs_per_cell)
771  do_prolongation = false;
772  }
773  Assert(do_prolongation,
775 
776  if (do_prolongation)
777  {
778  FullMatrix<double> prolongate (this->dofs_per_cell,
779  this->dofs_per_cell);
780 
781  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
782  for (unsigned int j=0; j<this->dofs_per_cell; ++j)
783  {
784  if (this->system_to_base_table[i].first !=
785  this->system_to_base_table[j].first)
786  continue;
787  const unsigned int
788  base = this->system_to_base_table[i].first.first;
789 
790  const unsigned int
791  base_index_i = this->system_to_base_table[i].second,
792  base_index_j = this->system_to_base_table[j].second;
793  prolongate(i,j) = (*base_matrices[base])(base_index_i,base_index_j);
794  }
795  prolongate.swap(const_cast<FullMatrix<double> &>
796  (this->prolongation[refinement_case-1][child]));
797  }
798  }
799 
800  return this->prolongation[refinement_case-1][child];
801 }
802 
803 
804 template <int dim, int spacedim>
805 unsigned int
807 face_to_cell_index (const unsigned int face_dof_index,
808  const unsigned int face,
809  const bool face_orientation,
810  const bool face_flip,
811  const bool face_rotation) const
812 {
813  // we need to ask the base elements how they want to translate
814  // the DoFs within their own numbering. thus, translate to
815  // the base element numbering and then back
816  const std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
817  face_base_index = this->face_system_to_base_index(face_dof_index);
818 
819  const unsigned int
820  base_face_to_cell_index
821  = this->base_element(face_base_index.first.first).face_to_cell_index (face_base_index.second,
822  face,
823  face_orientation,
824  face_flip,
825  face_rotation);
826 
827  // it would be nice if we had a base_to_system_index function, but
828  // all that exists is a component_to_system_index function. we can't do
829  // this here because it won't work for non-primitive elements. consequently,
830  // simply do a loop over all dofs till we find whether it corresponds
831  // to the one we're interested in -- crude, maybe, but works for now
832  const std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
833  target = std::make_pair (face_base_index.first, base_face_to_cell_index);
834  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
835  if (this->system_to_base_index(i) == target)
836  return i;
837 
838  Assert (false, ExcInternalError());
840 }
841 
842 
843 
844 //---------------------------------------------------------------------------
845 // Data field initialization
846 //---------------------------------------------------------------------------
847 
848 
849 
850 template <int dim, int spacedim>
853 {
855  // generate maximal set of flags
856  // that are necessary
857  for (unsigned int base_no=0; base_no<this->n_base_elements(); ++base_no)
858  out |= base_element(base_no).requires_update_flags (flags);
859  return out;
860 }
861 
862 
863 
864 template <int dim, int spacedim>
867 get_data (const UpdateFlags flags,
868  const Mapping<dim,spacedim> &mapping,
869  const Quadrature<dim> &quadrature,
871 {
872  // create an internal data object and set the update flags we will need
873  // to deal with. the current object does not make use of these flags,
874  // but we need to nevertheless set them correctly since we look
875  // into the update_each flag of base elements in fill_fe_values,
876  // and so the current object's update_each flag needs to be
877  // correct in case the current FESystem is a base element for another,
878  // higher-level FESystem itself.
879  InternalData *data = new InternalData(this->n_base_elements());
880  data->update_each = requires_update_flags(flags);
881 
882  // get data objects from each of the base elements and store
883  // them. one might think that doing this in parallel (over the
884  // base elements) would be a good idea, but this turns out to
885  // be wrong because we would then run these jobs on different
886  // threads/processors and this allocates memory in different
887  // NUMA domains; this has large detrimental effects when later
888  // writing into these objects in fill_fe_*_values. all of this
889  // is particularly true when using FEValues objects in
890  // WorkStream contexts where we explicitly make sure that
891  // every function only uses objects previously allocated
892  // in the same NUMA context and on the same thread as the
893  // function is called
894  for (unsigned int base_no=0; base_no<this->n_base_elements(); ++base_no)
895  {
897  = data->get_fe_output_object(base_no);
898  base_fe_output_object.initialize (quadrature.size(), base_element(base_no),
899  flags | base_element(base_no).requires_update_flags(flags));
900 
901  // let base objects produce their scratch objects. they may
902  // also at this time write into the output objects we provide
903  // for them; it would be nice if we could already copy something
904  // out of the base output object into the system output object,
905  // but we can't because we can't know what the elements already
906  // copied and/or will want to update on every cell
907  typename FiniteElement<dim,spacedim>::InternalDataBase *base_fe_data =
908  base_element(base_no).get_data (flags, mapping, quadrature,
909  base_fe_output_object);
910 
911  data->set_fe_data(base_no, base_fe_data);
912  }
913 
914  return data;
915 }
916 
917 // The following function is a clone of get_data, with the exception
918 // that get_face_data of the base elements is called.
919 
920 template <int dim, int spacedim>
924  const Mapping<dim,spacedim> &mapping,
925  const Quadrature<dim-1> &quadrature,
927 {
928  // create an internal data object and set the update flags we will need
929  // to deal with. the current object does not make use of these flags,
930  // but we need to nevertheless set them correctly since we look
931  // into the update_each flag of base elements in fill_fe_values,
932  // and so the current object's update_each flag needs to be
933  // correct in case the current FESystem is a base element for another,
934  // higher-level FESystem itself.
935  InternalData *data = new InternalData(this->n_base_elements());
936  data->update_each = requires_update_flags(flags);
937 
938  // get data objects from each of the base elements and store
939  // them. one might think that doing this in parallel (over the
940  // base elements) would be a good idea, but this turns out to
941  // be wrong because we would then run these jobs on different
942  // threads/processors and this allocates memory in different
943  // NUMA domains; this has large detrimental effects when later
944  // writing into these objects in fill_fe_*_values. all of this
945  // is particularly true when using FEValues objects in
946  // WorkStream contexts where we explicitly make sure that
947  // every function only uses objects previously allocated
948  // in the same NUMA context and on the same thread as the
949  // function is called
950  for (unsigned int base_no=0; base_no<this->n_base_elements(); ++base_no)
951  {
953  = data->get_fe_output_object(base_no);
954  base_fe_output_object.initialize (quadrature.size(), base_element(base_no),
955  flags | base_element(base_no).requires_update_flags(flags));
956 
957  // let base objects produce their scratch objects. they may
958  // also at this time write into the output objects we provide
959  // for them; it would be nice if we could already copy something
960  // out of the base output object into the system output object,
961  // but we can't because we can't know what the elements already
962  // copied and/or will want to update on every cell
963  typename FiniteElement<dim,spacedim>::InternalDataBase *base_fe_data =
964  base_element(base_no).get_face_data (flags, mapping, quadrature,
965  base_fe_output_object);
966 
967  data->set_fe_data(base_no, base_fe_data);
968  }
969 
970  return data;
971 }
972 
973 
974 
975 // The following function is a clone of get_data, with the exception
976 // that get_subface_data of the base elements is called.
977 
978 template <int dim, int spacedim>
982  const Mapping<dim,spacedim> &mapping,
983  const Quadrature<dim-1> &quadrature,
985 {
986  // create an internal data object and set the update flags we will need
987  // to deal with. the current object does not make use of these flags,
988  // but we need to nevertheless set them correctly since we look
989  // into the update_each flag of base elements in fill_fe_values,
990  // and so the current object's update_each flag needs to be
991  // correct in case the current FESystem is a base element for another,
992  // higher-level FESystem itself.
993  InternalData *data = new InternalData(this->n_base_elements());
994  data->update_each = requires_update_flags(flags);
995 
996  // get data objects from each of the base elements and store
997  // them. one might think that doing this in parallel (over the
998  // base elements) would be a good idea, but this turns out to
999  // be wrong because we would then run these jobs on different
1000  // threads/processors and this allocates memory in different
1001  // NUMA domains; this has large detrimental effects when later
1002  // writing into these objects in fill_fe_*_values. all of this
1003  // is particularly true when using FEValues objects in
1004  // WorkStream contexts where we explicitly make sure that
1005  // every function only uses objects previously allocated
1006  // in the same NUMA context and on the same thread as the
1007  // function is called
1008  for (unsigned int base_no=0; base_no<this->n_base_elements(); ++base_no)
1009  {
1011  = data->get_fe_output_object(base_no);
1012  base_fe_output_object.initialize (quadrature.size(), base_element(base_no),
1013  flags | base_element(base_no).requires_update_flags(flags));
1014 
1015  // let base objects produce their scratch objects. they may
1016  // also at this time write into the output objects we provide
1017  // for them; it would be nice if we could already copy something
1018  // out of the base output object into the system output object,
1019  // but we can't because we can't know what the elements already
1020  // copied and/or will want to update on every cell
1021  typename FiniteElement<dim,spacedim>::InternalDataBase *base_fe_data =
1022  base_element(base_no).get_subface_data (flags, mapping, quadrature,
1023  base_fe_output_object);
1024 
1025  data->set_fe_data(base_no, base_fe_data);
1026  }
1027 
1028  return data;
1029 }
1030 
1031 
1032 
1033 template <int dim, int spacedim>
1034 void
1037  const CellSimilarity::Similarity cell_similarity,
1038  const Quadrature<dim> &quadrature,
1039  const Mapping<dim,spacedim> &mapping,
1040  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
1041  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
1042  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
1044 {
1045  compute_fill(mapping, cell, invalid_face_number, invalid_face_number,
1046  quadrature, cell_similarity, mapping_internal, fe_internal,
1047  mapping_data, output_data);
1048 }
1049 
1050 
1051 
1052 template <int dim, int spacedim>
1053 void
1056  const unsigned int face_no,
1057  const Quadrature<dim-1> &quadrature,
1058  const Mapping<dim,spacedim> &mapping,
1059  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
1060  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
1061  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
1063 {
1064  compute_fill (mapping, cell, face_no, invalid_face_number, quadrature,
1065  CellSimilarity::none, mapping_internal, fe_internal,
1066  mapping_data, output_data);
1067 }
1068 
1069 
1070 
1071 
1072 template <int dim, int spacedim>
1073 void
1076  const unsigned int face_no,
1077  const unsigned int sub_no,
1078  const Quadrature<dim-1> &quadrature,
1079  const Mapping<dim,spacedim> &mapping,
1080  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
1081  const ::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
1082  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
1084 {
1085  compute_fill (mapping, cell, face_no, sub_no, quadrature,
1086  CellSimilarity::none, mapping_internal, fe_internal,
1087  mapping_data, output_data);
1088 }
1089 
1090 
1091 
1092 template <int dim, int spacedim>
1093 template <int dim_1>
1094 void
1097  const typename Triangulation<dim,spacedim>::cell_iterator &cell,
1098  const unsigned int face_no,
1099  const unsigned int sub_no,
1100  const Quadrature<dim_1> &quadrature,
1101  const CellSimilarity::Similarity cell_similarity,
1102  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
1103  const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
1106 {
1107  // convert data object to internal
1108  // data for this class. fails with
1109  // an exception if that is not
1110  // possible
1111  Assert (dynamic_cast<const InternalData *> (&fe_internal) != 0, ExcInternalError());
1112  const InternalData &fe_data = static_cast<const InternalData &> (fe_internal);
1113 
1114  // Either dim_1==dim
1115  // (fill_fe_values) or dim_1==dim-1
1116  // (fill_fe_(sub)face_values)
1117  Assert(dim_1==dim || dim_1==dim-1, ExcInternalError());
1118  const UpdateFlags flags = fe_data.update_each;
1119 
1120 
1121  // loop over the base elements, let them compute what they need to compute,
1122  // and then copy what is necessary.
1123  //
1124  // one may think that it would be a good idea to parallelize this over
1125  // base elements, but it turns out to be not worthwhile: doing so lets
1126  // multiple threads access data objects that were created by the current
1127  // thread, leading to many NUMA memory access inefficiencies. we specifically
1128  // want to avoid this if this class is called in a WorkStream context where
1129  // we very carefully allocate objects only on the thread where they
1130  // will actually be used; spawning new tasks here would be counterproductive
1131  if (flags & (update_values | update_gradients
1133  for (unsigned int base_no=0; base_no<this->n_base_elements(); ++base_no)
1134  {
1136  base_fe = base_element(base_no);
1138  base_fe_data = fe_data.get_fe_data(base_no);
1140  base_data = fe_data.get_fe_output_object(base_no);
1141 
1142  // fill_fe_face_values needs argument Quadrature<dim-1> for both cases
1143  // dim_1==dim-1 and dim_1=dim. Hence the following workaround
1144  const Quadrature<dim> *cell_quadrature = 0;
1145  const Quadrature<dim-1> *face_quadrature = 0;
1146  const unsigned int n_q_points = quadrature.size();
1147 
1148  // static cast to the common base class of quadrature being either
1149  // Quadrature<dim> or Quadrature<dim-1>:
1150  const Subscriptor *quadrature_base_pointer = &quadrature;
1151 
1152  if (face_no==invalid_face_number)
1153  {
1154  Assert(dim_1==dim, ExcDimensionMismatch(dim_1,dim));
1155  Assert (dynamic_cast<const Quadrature<dim> *>(quadrature_base_pointer) != 0,
1156  ExcInternalError());
1157 
1158  cell_quadrature
1159  = static_cast<const Quadrature<dim> *>(quadrature_base_pointer);
1160  }
1161  else
1162  {
1163  Assert(dim_1==dim-1, ExcDimensionMismatch(dim_1,dim-1));
1164  Assert (dynamic_cast<const Quadrature<dim-1> *>(quadrature_base_pointer) != 0,
1165  ExcInternalError());
1166 
1167  face_quadrature
1168  = static_cast<const Quadrature<dim-1> *>(quadrature_base_pointer);
1169  }
1170 
1171 
1172  // Make sure that in the case of fill_fe_values the data is only
1173  // copied from base_data to data if base_data is changed. therefore
1174  // use fe_fe_data.current_update_flags()
1175  //
1176  // for the case of fill_fe_(sub)face_values the data needs to be
1177  // copied from base_data to data on each face, therefore use
1178  // base_fe_data.update_flags.
1179  if (face_no==invalid_face_number)
1180  base_fe.fill_fe_values(cell, cell_similarity,
1181  *cell_quadrature,
1182  mapping, mapping_internal, mapping_data,
1183  base_fe_data, base_data);
1184  else if (sub_no==invalid_face_number)
1185  base_fe.fill_fe_face_values(cell, face_no,
1186  *face_quadrature,
1187  mapping, mapping_internal, mapping_data,
1188  base_fe_data, base_data);
1189  else
1190  base_fe.fill_fe_subface_values(cell, face_no, sub_no,
1191  *face_quadrature,
1192  mapping, mapping_internal, mapping_data,
1193  base_fe_data, base_data);
1194 
1195  // now data has been generated, so copy it. we used to work by
1196  // looping over all base elements (i.e. this outer loop), then over
1197  // multiplicity, then over the shape functions from that base
1198  // element, but that requires that we can infer the global number of
1199  // a shape function from its number in the base element. for that we
1200  // had the component_to_system_table.
1201  //
1202  // however, this does of course no longer work since we have
1203  // non-primitive elements. so we go the other way round: loop over
1204  // all shape functions of the composed element, and here only treat
1205  // those shape functions that belong to a given base element
1206  //TODO: Introduce the needed table and loop only over base element shape functions. This here is not efficient at all AND very bad style
1207  const UpdateFlags base_flags = base_fe_data.update_each;
1208 
1209  // some base element might involve values that depend on the shape
1210  // of the geometry, so we always need to copy the shape values around
1211  // also in case we detected a cell similarity (but no heavy work will
1212  // be done inside the individual elements in case we have a
1213  // translation and simple elements).
1214  for (unsigned int system_index=0; system_index<this->dofs_per_cell;
1215  ++system_index)
1216  if (this->system_to_base_table[system_index].first.first == base_no)
1217  {
1218  const unsigned int
1219  base_index = this->system_to_base_table[system_index].second;
1220  Assert (base_index<base_fe.dofs_per_cell, ExcInternalError());
1221 
1222  // now copy. if the shape function is primitive, then there
1223  // is only one value to be copied, but for non-primitive
1224  // elements, there might be more values to be copied
1225  //
1226  // so, find out from which index to take this one value, and
1227  // to which index to put
1228  unsigned int out_index = 0;
1229  for (unsigned int i=0; i<system_index; ++i)
1230  out_index += this->n_nonzero_components(i);
1231  unsigned int in_index = 0;
1232  for (unsigned int i=0; i<base_index; ++i)
1233  in_index += base_fe.n_nonzero_components(i);
1234 
1235  // then loop over the number of components to be copied
1236  Assert (this->n_nonzero_components(system_index) ==
1237  base_fe.n_nonzero_components(base_index),
1238  ExcInternalError());
1239 
1240  if (base_flags & update_values)
1241  for (unsigned int s=0; s<this->n_nonzero_components(system_index); ++s)
1242  for (unsigned int q=0; q<n_q_points; ++q)
1243  output_data.shape_values[out_index+s][q] =
1244  base_data.shape_values(in_index+s,q);
1245 
1246  if (base_flags & update_gradients)
1247  for (unsigned int s=0; s<this->n_nonzero_components(system_index); ++s)
1248  for (unsigned int q=0; q<n_q_points; ++q)
1249  output_data.shape_gradients[out_index+s][q] =
1250  base_data.shape_gradients[in_index+s][q];
1251 
1252  if (base_flags & update_hessians)
1253  for (unsigned int s=0; s<this->n_nonzero_components(system_index); ++s)
1254  for (unsigned int q=0; q<n_q_points; ++q)
1255  output_data.shape_hessians[out_index+s][q] =
1256  base_data.shape_hessians[in_index+s][q];
1257 
1258  if (base_flags & update_3rd_derivatives)
1259  for (unsigned int s=0; s<this->n_nonzero_components(system_index); ++s)
1260  for (unsigned int q=0; q<n_q_points; ++q)
1261  output_data.shape_3rd_derivatives[out_index+s][q] =
1262  base_data.shape_3rd_derivatives[in_index+s][q];
1263 
1264  }
1265  }
1266 }
1267 
1268 
1269 template <int dim, int spacedim>
1271 {
1272  // check whether all base elements implement their interface constraint
1273  // matrices. if this is not the case, then leave the interface costraints of
1274  // this composed element empty as well; however, the rest of the element is
1275  // usable
1276  for (unsigned int base=0; base<this->n_base_elements(); ++base)
1277  if (base_element(base).constraints_are_implemented() == false)
1278  return;
1279 
1280  this->interface_constraints.
1281  TableBase<2,double>::reinit (this->interface_constraints_size());
1282 
1283  // the layout of the constraints matrix is described in the FiniteElement
1284  // class. you may want to look there first before trying to understand the
1285  // following, especially the mapping of the @p{m} index.
1286  //
1287  // in order to map it to the fe-system class, we have to know which base
1288  // element a degree of freedom within a vertex, line, etc belongs to. this
1289  // can be accomplished by the system_to_component_index function in
1290  // conjunction with the numbers first_{line,quad,...}_index
1291  for (unsigned int n=0; n<this->interface_constraints.n(); ++n)
1292  for (unsigned int m=0; m<this->interface_constraints.m(); ++m)
1293  {
1294  // for the pair (n,m) find out which base element they belong to and
1295  // the number therein
1296  //
1297  // first for the n index. this is simple since the n indices are in
1298  // the same order as they are usually on a face. note that for the
1299  // data type, first value in pair is (base element,instance of base
1300  // element), second is index within this instance
1301  const std::pair<std::pair<unsigned int,unsigned int>, unsigned int> n_index
1302  = this->face_system_to_base_table[n];
1303 
1304  // likewise for the m index. this is more complicated due to the
1305  // strange ordering we have for the dofs on the refined faces.
1306  std::pair<std::pair<unsigned int,unsigned int>, unsigned int> m_index;
1307  switch (dim)
1308  {
1309  case 1:
1310  {
1311  // we should never get here! (in 1d, the constraints matrix
1312  // should be of size zero)
1313  Assert (false, ExcInternalError());
1314  break;
1315  };
1316 
1317  case 2:
1318  {
1319  // the indices m=0..d_v-1 are from the center vertex. their order
1320  // is the same as for the first vertex of the whole cell, so we
1321  // can use the system_to_base_table variable (using the
1322  // face_s_t_base_t function would yield the same)
1323  if (m < this->dofs_per_vertex)
1324  m_index = this->system_to_base_table[m];
1325  else
1326  // then come the two sets of line indices
1327  {
1328  const unsigned int index_in_line
1329  = (m-this->dofs_per_vertex) % this->dofs_per_line;
1330  const unsigned int sub_line
1331  = (m-this->dofs_per_vertex) / this->dofs_per_line;
1332  Assert (sub_line < 2, ExcInternalError());
1333 
1334  // from this information, try to get base element and instance
1335  // of base element. we do so by constructing the corresponding
1336  // face index of m in the present element, then use
1337  // face_system_to_base_table
1338  const unsigned int tmp1 = 2*this->dofs_per_vertex+index_in_line;
1339  m_index.first = this->face_system_to_base_table[tmp1].first;
1340 
1341  // what we are still missing is the index of m within the base
1342  // elements interface_constraints table
1343  //
1344  // here, the second value of face_system_to_base_table can
1345  // help: it denotes the face index of that shape function
1346  // within the base element. since we know that it is a line
1347  // dof, we can construct the rest: tmp2 will denote the index
1348  // of this shape function among the line shape functions:
1349  Assert (this->face_system_to_base_table[tmp1].second >=
1350  2*base_element(m_index.first.first).dofs_per_vertex,
1351  ExcInternalError());
1352  const unsigned int tmp2 = this->face_system_to_base_table[tmp1].second -
1353  2*base_element(m_index.first.first).dofs_per_vertex;
1354  Assert (tmp2 < base_element(m_index.first.first).dofs_per_line,
1355  ExcInternalError());
1356  m_index.second = base_element(m_index.first.first).dofs_per_vertex +
1357  base_element(m_index.first.first).dofs_per_line*sub_line +
1358  tmp2;
1359  };
1360  break;
1361  };
1362 
1363  case 3:
1364  {
1365  // same way as above, although a little more complicated...
1366 
1367  // the indices m=0..5*d_v-1 are from the center and the four
1368  // subline vertices. their order is the same as for the first
1369  // vertex of the whole cell, so we can use the simple arithmetic
1370  if (m < 5*this->dofs_per_vertex)
1371  m_index = this->system_to_base_table[m];
1372  else
1373  // then come the 12 sets of line indices
1374  if (m < 5*this->dofs_per_vertex + 12*this->dofs_per_line)
1375  {
1376  // for the meaning of all this, see the 2d part
1377  const unsigned int index_in_line
1378  = (m-5*this->dofs_per_vertex) % this->dofs_per_line;
1379  const unsigned int sub_line
1380  = (m-5*this->dofs_per_vertex) / this->dofs_per_line;
1381  Assert (sub_line < 12, ExcInternalError());
1382 
1383  const unsigned int tmp1 = 4*this->dofs_per_vertex+index_in_line;
1384  m_index.first = this->face_system_to_base_table[tmp1].first;
1385 
1386  Assert (this->face_system_to_base_table[tmp1].second >=
1387  4*base_element(m_index.first.first).dofs_per_vertex,
1388  ExcInternalError());
1389  const unsigned int tmp2 = this->face_system_to_base_table[tmp1].second -
1390  4*base_element(m_index.first.first).dofs_per_vertex;
1391  Assert (tmp2 < base_element(m_index.first.first).dofs_per_line,
1392  ExcInternalError());
1393  m_index.second = 5*base_element(m_index.first.first).dofs_per_vertex +
1394  base_element(m_index.first.first).dofs_per_line*sub_line +
1395  tmp2;
1396  }
1397  else
1398  // on one of the four sub-quads
1399  {
1400  // for the meaning of all this, see the 2d part
1401  const unsigned int index_in_quad
1402  = (m-5*this->dofs_per_vertex-12*this->dofs_per_line) %
1403  this->dofs_per_quad;
1404  Assert (index_in_quad < this->dofs_per_quad,
1405  ExcInternalError());
1406  const unsigned int sub_quad
1407  = ((m-5*this->dofs_per_vertex-12*this->dofs_per_line) /
1408  this->dofs_per_quad);
1409  Assert (sub_quad < 4, ExcInternalError());
1410 
1411  const unsigned int tmp1 = 4*this->dofs_per_vertex +
1412  4*this->dofs_per_line +
1413  index_in_quad;
1414  Assert (tmp1 < this->face_system_to_base_table.size(),
1415  ExcInternalError());
1416  m_index.first = this->face_system_to_base_table[tmp1].first;
1417 
1418  Assert (this->face_system_to_base_table[tmp1].second >=
1419  4*base_element(m_index.first.first).dofs_per_vertex +
1420  4*base_element(m_index.first.first).dofs_per_line,
1421  ExcInternalError());
1422  const unsigned int tmp2 = this->face_system_to_base_table[tmp1].second -
1423  4*base_element(m_index.first.first).dofs_per_vertex -
1424  4*base_element(m_index.first.first).dofs_per_line;
1425  Assert (tmp2 < base_element(m_index.first.first).dofs_per_quad,
1426  ExcInternalError());
1427  m_index.second = 5*base_element(m_index.first.first).dofs_per_vertex +
1428  12*base_element(m_index.first.first).dofs_per_line +
1429  base_element(m_index.first.first).dofs_per_quad*sub_quad +
1430  tmp2;
1431  };
1432 
1433  break;
1434  };
1435 
1436  default:
1437  Assert (false, ExcNotImplemented());
1438  };
1439 
1440  // now that we gathered all information: use it to build the
1441  // matrix. note that if n and m belong to different base elements or
1442  // instances, then there definitely will be no coupling
1443  if (n_index.first == m_index.first)
1444  this->interface_constraints(m,n)
1445  = (base_element(n_index.first.first).constraints()(m_index.second,
1446  n_index.second));
1447  };
1448 }
1449 
1450 
1451 
1452 template <int dim, int spacedim>
1454  const std::vector<unsigned int> &multiplicities)
1455 {
1456  Assert (fes.size() == multiplicities.size(),
1457  ExcDimensionMismatch (fes.size(), multiplicities.size()) );
1458  Assert (fes.size() > 0,
1459  ExcMessage ("Need to pass at least one finite element."));
1460  Assert (count_nonzeros(multiplicities) > 0,
1461  ExcMessage("You only passed FiniteElements with multiplicity 0."));
1462 
1463  // Note that we need to skip every fe with multiplicity 0 in the following block of code
1464 
1465  this->base_to_block_indices.reinit(0, 0);
1466 
1467  for (unsigned int i=0; i<fes.size(); i++)
1468  if (multiplicities[i]>0)
1469  this->base_to_block_indices.push_back( multiplicities[i] );
1470 
1471  std::vector<Threads::Task<FiniteElement<dim,spacedim>*> > clone_base_elements;
1472 
1473  for (unsigned int i=0; i<fes.size(); i++)
1474  if (multiplicities[i]>0)
1475  clone_base_elements.push_back (Threads::new_task (&FiniteElement<dim,spacedim>::clone,
1476  *fes[i]));
1477 
1478  unsigned int ind=0;
1479  for (unsigned int i=0; i<fes.size(); i++)
1480  {
1481  if (multiplicities[i]>0)
1482  {
1483  base_elements[ind] =
1484  std::make_pair (std_cxx11::shared_ptr<const FiniteElement<dim,spacedim> >
1485  (clone_base_elements[ind].return_value()),
1486  multiplicities[i]);
1487  ++ind;
1488  }
1489  }
1490 
1491  Assert(ind>0, ExcInternalError());
1492 
1493  {
1494  // If the system is not primitive, these have not been initialized by
1495  // FiniteElement
1496  this->system_to_component_table.resize(this->dofs_per_cell);
1497  this->face_system_to_component_table.resize(this->dofs_per_face);
1498 
1499  FETools::Compositing::build_cell_tables(this->system_to_base_table,
1500  this->system_to_component_table,
1501  this->component_to_base_table,
1502  *this);
1503 
1504  FETools::Compositing::build_face_tables(this->face_system_to_base_table,
1505  this->face_system_to_component_table,
1506  *this);
1507 
1508  }
1509 
1510  // restriction and prolongation matrices are built on demand
1511 
1512  // now set up the interface constraints. this is kind'o hairy, so don't try
1513  // to do it dimension independent
1514  build_interface_constraints ();
1515 
1516  // finally fill in support points on cell and face
1517  initialize_unit_support_points ();
1518  initialize_unit_face_support_points ();
1519 
1520  initialize_quad_dof_index_permutation ();
1521 }
1522 
1523 
1524 
1525 
1526 
1527 
1528 template <int dim, int spacedim>
1529 void
1532 {
1533  // if one of the base elements has no support points, then it makes no sense
1534  // to define support points for the composed element, so return an empty
1535  // array to demonstrate that fact. Note that we ignore FE_Nothing in this logic.
1536  for (unsigned int base_el=0; base_el<this->n_base_elements(); ++base_el)
1537  if (!base_element(base_el).has_support_points() && base_element(base_el).dofs_per_cell!=0)
1538  {
1539  this->unit_support_points.resize(0);
1540  return;
1541  };
1542 
1543  // generate unit support points from unit support points of sub elements
1544  this->unit_support_points.resize(this->dofs_per_cell);
1545 
1546  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
1547  {
1548  const unsigned int
1549  base = this->system_to_base_table[i].first.first,
1550  base_index = this->system_to_base_table[i].second;
1551  Assert (base<this->n_base_elements(), ExcInternalError());
1552  Assert (base_index<base_element(base).unit_support_points.size(),
1553  ExcInternalError());
1554  this->unit_support_points[i] = base_element(base).unit_support_points[base_index];
1555  };
1556 }
1557 
1558 
1559 
1560 template <int dim, int spacedim>
1561 void
1564 {
1565  // Nothing to do in 1D
1566  if (dim == 1)
1567  return;
1568 
1569  // if one of the base elements has no support points, then it makes no sense
1570  // to define support points for the composed element. In that case, return
1571  // an empty array to demonstrate that fact (note that we ask whether the
1572  // base element has no support points at all, not only none on the face!)
1573  //
1574  // on the other hand, if there is an element that simply has no degrees of
1575  // freedom on the face at all, then we don't care whether it has support
1576  // points or not. this is, for example, the case for the stable Stokes
1577  // element Q(p)^dim \times DGP(p-1).
1578  for (unsigned int base_el=0; base_el<this->n_base_elements(); ++base_el)
1579  if (!base_element(base_el).has_support_points()
1580  &&
1581  (base_element(base_el).dofs_per_face > 0))
1582  {
1583  this->unit_face_support_points.resize(0);
1584  return;
1585  }
1586 
1587 
1588  // generate unit face support points from unit support points of sub
1589  // elements
1590  this->unit_face_support_points.resize(this->dofs_per_face);
1591 
1592  for (unsigned int i=0; i<this->dofs_per_face; ++i)
1593  {
1594  const unsigned int base_i = this->face_system_to_base_table[i].first.first;
1595  const unsigned int index_in_base = this->face_system_to_base_table[i].second;
1596 
1597  Assert (index_in_base < base_element(base_i).unit_face_support_points.size(),
1598  ExcInternalError());
1599 
1600  this->unit_face_support_points[i]
1601  = base_element(base_i).unit_face_support_points[index_in_base];
1602  }
1603 }
1604 
1605 
1606 
1607 template <int dim, int spacedim>
1608 void
1610 {
1611  // nothing to do in other dimensions than 3
1612  if (dim < 3)
1613  return;
1614 
1615  // the array into which we want to write should have the correct size
1616  // already.
1617  Assert (this->adjust_quad_dof_index_for_face_orientation_table.n_elements()==
1618  8*this->dofs_per_quad, ExcInternalError());
1619 
1620  // to obtain the shifts for this composed element, copy the shift
1621  // information of the base elements
1622  unsigned int index = 0;
1623  for (unsigned int b=0; b<this->n_base_elements(); ++b)
1624  {
1625  const Table<2,int> &temp
1626  = this->base_element(b).adjust_quad_dof_index_for_face_orientation_table;
1627  for (unsigned int c=0; c<this->element_multiplicity(b); ++c)
1628  {
1629  for (unsigned int i=0; i<temp.size(0); ++i)
1630  for (unsigned int j=0; j<8; ++j)
1631  this->adjust_quad_dof_index_for_face_orientation_table(index+i,j)=
1632  temp(i,j);
1633  index += temp.size(0);
1634  }
1635  }
1636  Assert (index == this->dofs_per_quad,
1637  ExcInternalError());
1638 
1639  // aditionally compose the permutation information for lines
1640  Assert (this->adjust_line_dof_index_for_line_orientation_table.size()==
1641  this->dofs_per_line, ExcInternalError());
1642  index = 0;
1643  for (unsigned int b=0; b<this->n_base_elements(); ++b)
1644  {
1645  const std::vector<int> &temp2
1646  = this->base_element(b).adjust_line_dof_index_for_line_orientation_table;
1647  for (unsigned int c=0; c<this->element_multiplicity(b); ++c)
1648  {
1649  std::copy(temp2.begin(), temp2.end(),
1650  this->adjust_line_dof_index_for_line_orientation_table.begin()
1651  +index);
1652  index += temp2.size();
1653  }
1654  }
1655  Assert (index == this->dofs_per_line,
1656  ExcInternalError());
1657 }
1658 
1659 
1660 
1661 template <int dim, int spacedim>
1662 bool
1665 {
1666  for (unsigned int b=0; b<this->n_base_elements(); ++b)
1667  if (base_element(b).hp_constraints_are_implemented() == false)
1668  return false;
1669 
1670  return true;
1671 }
1672 
1673 
1674 
1675 template <int dim, int spacedim>
1676 void
1679  FullMatrix<double> &interpolation_matrix) const
1680 {
1681  typedef FiniteElement<dim,spacedim> FEL;
1682  AssertThrow ((x_source_fe.get_name().find ("FE_System<") == 0)
1683  ||
1684  (dynamic_cast<const FESystem<dim,spacedim>*>(&x_source_fe) != 0),
1685  typename FEL::
1686  ExcInterpolationNotImplemented());
1687 
1688  Assert (interpolation_matrix.n() == this->dofs_per_face,
1689  ExcDimensionMismatch (interpolation_matrix.n(),
1690  this->dofs_per_face));
1691  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
1692  ExcDimensionMismatch (interpolation_matrix.m(),
1693  x_source_fe.dofs_per_face));
1694 
1695  // since dofs for each base are independent, we only have to stack things up
1696  // from base element to base element
1697  //
1698  // the problem is that we have to work with two FEs (this and
1699  // fe_other). only deal with the case that both are FESystems and that they
1700  // both have the same number of bases (counting multiplicity) each of which
1701  // match in their number of components. this covers
1702  // FESystem(FE_Q(p),1,FE_Q(q),2) vs FESystem(FE_Q(r),2,FE_Q(s),1), but not
1703  // FESystem(FE_Q(p),1,FE_Q(q),2) vs
1704  // FESystem(FESystem(FE_Q(r),2),1,FE_Q(s),1)
1705  const FESystem<dim,spacedim> *fe_other_system
1706  = dynamic_cast<const FESystem<dim,spacedim>*>(&x_source_fe);
1707 
1708  // clear matrix, since we will not get to set all elements
1709  interpolation_matrix = 0;
1710 
1711  // loop over all the base elements of this and the other element, counting
1712  // their multiplicities
1713  unsigned int base_index = 0,
1714  base_index_other = 0;
1715  unsigned int multiplicity = 0,
1716  multiplicity_other = 0;
1717 
1718  FullMatrix<double> base_to_base_interpolation;
1719 
1720  while (true)
1721  {
1723  &base = base_element(base_index),
1724  &base_other = fe_other_system->base_element(base_index_other);
1725 
1726  Assert (base.n_components() == base_other.n_components(),
1727  ExcNotImplemented());
1728 
1729  // get the interpolation from the bases
1730  base_to_base_interpolation.reinit (base_other.dofs_per_face,
1731  base.dofs_per_face);
1732  base.get_face_interpolation_matrix (base_other,
1733  base_to_base_interpolation);
1734 
1735  // now translate entries. we'd like to have something like
1736  // face_base_to_system_index, but that doesn't exist. rather, all we
1737  // have is the reverse. well, use that then
1738  for (unsigned int i=0; i<this->dofs_per_face; ++i)
1739  if (this->face_system_to_base_index(i).first
1740  ==
1741  std::make_pair (base_index, multiplicity))
1742  for (unsigned int j=0; j<fe_other_system->dofs_per_face; ++j)
1743  if (fe_other_system->face_system_to_base_index(j).first
1744  ==
1745  std::make_pair (base_index_other, multiplicity_other))
1746  interpolation_matrix(j, i)
1747  = base_to_base_interpolation(fe_other_system->face_system_to_base_index(j).second,
1748  this->face_system_to_base_index(i).second);
1749 
1750  // advance to the next base element for this and the other fe_system;
1751  // see if we can simply advance the multiplicity by one, or if have to
1752  // move on to the next base element
1753  ++multiplicity;
1754  if (multiplicity == this->element_multiplicity(base_index))
1755  {
1756  multiplicity = 0;
1757  ++base_index;
1758  }
1759  ++multiplicity_other;
1760  if (multiplicity_other ==
1761  fe_other_system->element_multiplicity(base_index_other))
1762  {
1763  multiplicity_other = 0;
1764  ++base_index_other;
1765  }
1766 
1767  // see if we have reached the end of the present element. if so, we
1768  // should have reached the end of the other one as well
1769  if (base_index == this->n_base_elements())
1770  {
1771  Assert (base_index_other == fe_other_system->n_base_elements(),
1772  ExcInternalError());
1773  break;
1774  }
1775 
1776  // if we haven't reached the end of this element, we shouldn't have
1777  // reached the end of the other one either
1778  Assert (base_index_other != fe_other_system->n_base_elements(),
1779  ExcInternalError());
1780  }
1781 }
1782 
1783 
1784 
1785 template <int dim, int spacedim>
1786 void
1789  const unsigned int subface,
1790  FullMatrix<double> &interpolation_matrix) const
1791 {
1792  typedef FiniteElement<dim,spacedim> FEL;
1793  AssertThrow ((x_source_fe.get_name().find ("FE_System<") == 0)
1794  ||
1795  (dynamic_cast<const FESystem<dim,spacedim>*>(&x_source_fe) != 0),
1796  typename FEL::
1797  ExcInterpolationNotImplemented());
1798 
1799  Assert (interpolation_matrix.n() == this->dofs_per_face,
1800  ExcDimensionMismatch (interpolation_matrix.n(),
1801  this->dofs_per_face));
1802  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
1803  ExcDimensionMismatch (interpolation_matrix.m(),
1804  x_source_fe.dofs_per_face));
1805 
1806  // since dofs for each base are independent, we only have to stack things up
1807  // from base element to base element
1808  //
1809  // the problem is that we have to work with two FEs (this and
1810  // fe_other). only deal with the case that both are FESystems and that they
1811  // both have the same number of bases (counting multiplicity) each of which
1812  // match in their number of components. this covers
1813  // FESystem(FE_Q(p),1,FE_Q(q),2) vs FESystem(FE_Q(r),2,FE_Q(s),1), but not
1814  // FESystem(FE_Q(p),1,FE_Q(q),2) vs
1815  // FESystem(FESystem(FE_Q(r),2),1,FE_Q(s),1)
1816  const FESystem<dim,spacedim> *fe_other_system
1817  = dynamic_cast<const FESystem<dim,spacedim>*>(&x_source_fe);
1818 
1819  // clear matrix, since we will not get to set all elements
1820  interpolation_matrix = 0;
1821 
1822  // loop over all the base elements of this and the other element, counting
1823  // their multiplicities
1824  unsigned int base_index = 0,
1825  base_index_other = 0;
1826  unsigned int multiplicity = 0,
1827  multiplicity_other = 0;
1828 
1829  FullMatrix<double> base_to_base_interpolation;
1830 
1831  while (true)
1832  {
1834  &base = base_element(base_index),
1835  &base_other = fe_other_system->base_element(base_index_other);
1836 
1837  Assert (base.n_components() == base_other.n_components(),
1838  ExcNotImplemented());
1839 
1840  // get the interpolation from the bases
1841  base_to_base_interpolation.reinit (base_other.dofs_per_face,
1842  base.dofs_per_face);
1843  base.get_subface_interpolation_matrix (base_other,
1844  subface,
1845  base_to_base_interpolation);
1846 
1847  // now translate entries. we'd like to have something like
1848  // face_base_to_system_index, but that doesn't exist. rather, all we
1849  // have is the reverse. well, use that then
1850  for (unsigned int i=0; i<this->dofs_per_face; ++i)
1851  if (this->face_system_to_base_index(i).first
1852  ==
1853  std::make_pair (base_index, multiplicity))
1854  for (unsigned int j=0; j<fe_other_system->dofs_per_face; ++j)
1855  if (fe_other_system->face_system_to_base_index(j).first
1856  ==
1857  std::make_pair (base_index_other, multiplicity_other))
1858  interpolation_matrix(j, i)
1859  = base_to_base_interpolation(fe_other_system->face_system_to_base_index(j).second,
1860  this->face_system_to_base_index(i).second);
1861 
1862  // advance to the next base element for this and the other fe_system;
1863  // see if we can simply advance the multiplicity by one, or if have to
1864  // move on to the next base element
1865  ++multiplicity;
1866  if (multiplicity == this->element_multiplicity(base_index))
1867  {
1868  multiplicity = 0;
1869  ++base_index;
1870  }
1871  ++multiplicity_other;
1872  if (multiplicity_other ==
1873  fe_other_system->element_multiplicity(base_index_other))
1874  {
1875  multiplicity_other = 0;
1876  ++base_index_other;
1877  }
1878 
1879  // see if we have reached the end of the present element. if so, we
1880  // should have reached the end of the other one as well
1881  if (base_index == this->n_base_elements())
1882  {
1883  Assert (base_index_other == fe_other_system->n_base_elements(),
1884  ExcInternalError());
1885  break;
1886  }
1887 
1888  // if we haven't reached the end of this element, we shouldn't have
1889  // reached the end of the other one either
1890  Assert (base_index_other != fe_other_system->n_base_elements(),
1891  ExcInternalError());
1892  }
1893 }
1894 
1895 
1896 
1897 template <int dim, int spacedim>
1898 template <int structdim>
1899 std::vector<std::pair<unsigned int, unsigned int> >
1901 {
1902  // since dofs on each subobject (vertex, line, ...) are ordered such that
1903  // first come all from the first base element all multiplicities, then
1904  // second base element all multiplicities, etc., we simply have to stack all
1905  // the identities after each other
1906  //
1907  // the problem is that we have to work with two FEs (this and
1908  // fe_other). only deal with the case that both are FESystems and that they
1909  // both have the same number of bases (counting multiplicity) each of which
1910  // match in their number of components. this covers
1911  // FESystem(FE_Q(p),1,FE_Q(q),2) vs FESystem(FE_Q(r),2,FE_Q(s),1), but not
1912  // FESystem(FE_Q(p),1,FE_Q(q),2) vs
1913  // FESystem(FESystem(FE_Q(r),2),1,FE_Q(s),1)
1914  if (const FESystem<dim,spacedim> *fe_other_system
1915  = dynamic_cast<const FESystem<dim,spacedim>*>(&fe_other))
1916  {
1917  // loop over all the base elements of this and the other element,
1918  // counting their multiplicities
1919  unsigned int base_index = 0,
1920  base_index_other = 0;
1921  unsigned int multiplicity = 0,
1922  multiplicity_other = 0;
1923 
1924  // we also need to keep track of the number of dofs already treated for
1925  // each of the elements
1926  unsigned int dof_offset = 0,
1927  dof_offset_other = 0;
1928 
1929  std::vector<std::pair<unsigned int, unsigned int> > identities;
1930 
1931  while (true)
1932  {
1934  &base = base_element(base_index),
1935  &base_other = fe_other_system->base_element(base_index_other);
1936 
1937  Assert (base.n_components() == base_other.n_components(),
1938  ExcNotImplemented());
1939 
1940  // now translate the identities returned by the base elements to the
1941  // indices of this system element
1942  std::vector<std::pair<unsigned int, unsigned int> > base_identities;
1943  switch (structdim)
1944  {
1945  case 0:
1946  base_identities = base.hp_vertex_dof_identities (base_other);
1947  break;
1948  case 1:
1949  base_identities = base.hp_line_dof_identities (base_other);
1950  break;
1951  case 2:
1952  base_identities = base.hp_quad_dof_identities (base_other);
1953  break;
1954  default:
1955  Assert (false, ExcNotImplemented());
1956  }
1957 
1958  for (unsigned int i=0; i<base_identities.size(); ++i)
1959  identities.push_back
1960  (std::make_pair (base_identities[i].first + dof_offset,
1961  base_identities[i].second + dof_offset_other));
1962 
1963  // record the dofs treated above as already taken care of
1964  dof_offset += base.template n_dofs_per_object<structdim>();
1965  dof_offset_other += base_other.template n_dofs_per_object<structdim>();
1966 
1967  // advance to the next base element for this and the other
1968  // fe_system; see if we can simply advance the multiplicity by one,
1969  // or if have to move on to the next base element
1970  ++multiplicity;
1971  if (multiplicity == this->element_multiplicity(base_index))
1972  {
1973  multiplicity = 0;
1974  ++base_index;
1975  }
1976  ++multiplicity_other;
1977  if (multiplicity_other ==
1978  fe_other_system->element_multiplicity(base_index_other))
1979  {
1980  multiplicity_other = 0;
1981  ++base_index_other;
1982  }
1983 
1984  // see if we have reached the end of the present element. if so, we
1985  // should have reached the end of the other one as well
1986  if (base_index == this->n_base_elements())
1987  {
1988  Assert (base_index_other == fe_other_system->n_base_elements(),
1989  ExcInternalError());
1990  break;
1991  }
1992 
1993  // if we haven't reached the end of this element, we shouldn't have
1994  // reached the end of the other one either
1995  Assert (base_index_other != fe_other_system->n_base_elements(),
1996  ExcInternalError());
1997  }
1998 
1999  return identities;
2000  }
2001  else
2002  {
2003  Assert (false, ExcNotImplemented());
2004  return std::vector<std::pair<unsigned int, unsigned int> >();
2005  }
2006 }
2007 
2008 
2009 
2010 template <int dim, int spacedim>
2011 std::vector<std::pair<unsigned int, unsigned int> >
2013 {
2014  return hp_object_dof_identities<0> (fe_other);
2015 }
2016 
2017 template <int dim, int spacedim>
2018 std::vector<std::pair<unsigned int, unsigned int> >
2020 {
2021  return hp_object_dof_identities<1> (fe_other);
2022 }
2023 
2024 
2025 
2026 template <int dim, int spacedim>
2027 std::vector<std::pair<unsigned int, unsigned int> >
2029 {
2030  return hp_object_dof_identities<2> (fe_other);
2031 }
2032 
2033 
2034 
2035 template <int dim, int spacedim>
2039 {
2040  // at present all we can do is to compare with other FESystems that have the
2041  // same number of components and bases
2042  if (const FESystem<dim,spacedim> *fe_sys_other
2043  = dynamic_cast<const FESystem<dim,spacedim>*>(&fe_other))
2044  {
2045  Assert (this->n_components() == fe_sys_other->n_components(),
2046  ExcNotImplemented());
2047  Assert (this->n_base_elements() == fe_sys_other->n_base_elements(),
2048  ExcNotImplemented());
2049 
2052 
2053  // loop over all base elements and do some sanity checks
2054  for (unsigned int b=0; b<this->n_base_elements(); ++b)
2055  {
2056  Assert (this->base_element(b).n_components() ==
2057  fe_sys_other->base_element(b).n_components(),
2058  ExcNotImplemented());
2059  Assert (this->element_multiplicity(b) ==
2060  fe_sys_other->element_multiplicity(b),
2061  ExcNotImplemented());
2062 
2063  // for this pair of base elements, check who dominates and combine
2064  // with previous result
2066  base_domination = (this->base_element(b)
2067  .compare_for_face_domination (fe_sys_other->base_element(b)));
2068  domination = domination & base_domination;
2069  }
2070 
2071  return domination;
2072  }
2073 
2074  Assert (false, ExcNotImplemented());
2076 }
2077 
2078 
2079 
2080 template <int dim, int spacedim>
2082 FESystem<dim,spacedim>::base_element (const unsigned int index) const
2083 {
2084  Assert (index < base_elements.size(),
2085  ExcIndexRange(index, 0, base_elements.size()));
2086  return *base_elements[index].first;
2087 }
2088 
2089 
2090 
2091 template <int dim, int spacedim>
2092 bool
2093 FESystem<dim,spacedim>::has_support_on_face (const unsigned int shape_index,
2094  const unsigned int face_index) const
2095 {
2096  return (base_element(this->system_to_base_index(shape_index).first.first)
2097  .has_support_on_face(this->system_to_base_index(shape_index).second,
2098  face_index));
2099 }
2100 
2101 
2102 
2103 template <int dim, int spacedim>
2104 Point<dim>
2105 FESystem<dim,spacedim>::unit_support_point (const unsigned int index) const
2106 {
2107  Assert (index < this->dofs_per_cell,
2108  ExcIndexRange (index, 0, this->dofs_per_cell));
2109  typedef FiniteElement<dim,spacedim> FEL;
2110  Assert ((this->unit_support_points.size() == this->dofs_per_cell) ||
2111  (this->unit_support_points.size() == 0),
2112  typename FEL::ExcFEHasNoSupportPoints ());
2113 
2114  // let's see whether we have the information pre-computed
2115  if (this->unit_support_points.size() != 0)
2116  return this->unit_support_points[index];
2117  else
2118  // no. ask the base element whether it would like to provide this
2119  // information
2120  return (base_element(this->system_to_base_index(index).first.first)
2121  .unit_support_point(this->system_to_base_index(index).second));
2122 }
2123 
2124 
2125 
2126 template <int dim, int spacedim>
2127 Point<dim-1>
2128 FESystem<dim,spacedim>::unit_face_support_point (const unsigned int index) const
2129 {
2130  Assert (index < this->dofs_per_face,
2131  ExcIndexRange (index, 0, this->dofs_per_face));
2132  typedef FiniteElement<dim,spacedim> FEL;
2133  Assert ((this->unit_face_support_points.size() == this->dofs_per_face) ||
2134  (this->unit_face_support_points.size() == 0),
2135  typename FEL::ExcFEHasNoSupportPoints ());
2136 
2137  // let's see whether we have the information pre-computed
2138  if (this->unit_face_support_points.size() != 0)
2139  return this->unit_face_support_points[index];
2140  else
2141  // no. ask the base element whether it would like to provide this
2142  // information
2143  return (base_element(this->face_system_to_base_index(index).first.first)
2144  .unit_face_support_point(this->face_system_to_base_index(index).second));
2145 }
2146 
2147 
2148 
2149 template <int dim, int spacedim>
2150 std::pair<Table<2,bool>, std::vector<unsigned int> >
2152 {
2153  // Note that this->n_components() is actually only an estimate of how many
2154  // constant modes we will need. There might be more than one such mode
2155  // (e.g. FE_Q_DG0).
2156  Table<2,bool> constant_modes(this->n_components(), this->dofs_per_cell);
2157  std::vector<unsigned int> components;
2158  for (unsigned int i=0; i<base_elements.size(); ++i)
2159  {
2160  std::pair<Table<2,bool>, std::vector<unsigned int> >
2161  base_table = base_elements[i].first->get_constant_modes();
2162  AssertDimension(base_table.first.n_rows(), base_table.second.size());
2163  const unsigned int element_multiplicity = this->element_multiplicity(i);
2164 
2165  // there might be more than one constant mode for some scalar elements,
2166  // so make sure the table actually fits: Create a new table with more
2167  // rows
2168  const unsigned int comp = components.size();
2169  if (constant_modes.n_rows() < comp+base_table.first.n_rows()*element_multiplicity)
2170  {
2171  Table<2,bool> new_constant_modes(comp+base_table.first.n_rows()*
2172  element_multiplicity,
2173  constant_modes.n_cols());
2174  for (unsigned int r=0; r<comp; ++r)
2175  for (unsigned int c=0; c<this->dofs_per_cell; ++c)
2176  new_constant_modes(r,c) = constant_modes(r,c);
2177  constant_modes.swap(new_constant_modes);
2178  }
2179 
2180  // next, fill the constant modes from the individual components as well
2181  // as the component numbers corresponding to the constant mode rows
2182  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
2183  {
2184  std::pair<std::pair<unsigned int,unsigned int>, unsigned int> ind
2185  = this->system_to_base_index(k);
2186  if (ind.first.first == i)
2187  for (unsigned int c=0; c<base_table.first.n_rows(); ++c)
2188  constant_modes(comp+ind.first.second*base_table.first.n_rows()+c,k)
2189  = base_table.first(c,ind.second);
2190  }
2191  for (unsigned int r=0; r<element_multiplicity; ++r)
2192  for (unsigned int c=0; c<base_table.second.size(); ++c)
2193  components.push_back(comp+r*this->base_elements[i].first->n_components()
2194  +base_table.second[c]);
2195  }
2196  AssertDimension(components.size(), constant_modes.n_rows());
2197  return std::pair<Table<2,bool>, std::vector<unsigned int> >(constant_modes,
2198  components);
2199 }
2200 
2201 
2202 
2203 template <int dim, int spacedim>
2204 std::size_t
2206 {
2207  // neglect size of data stored in @p{base_elements} due to some problems
2208  // with the compiler. should be neglectable after all, considering the size
2209  // of the data of the subelements
2210  std::size_t mem = (FiniteElement<dim,spacedim>::memory_consumption () +
2211  sizeof (base_elements));
2212  for (unsigned int i=0; i<base_elements.size(); ++i)
2213  mem += MemoryConsumption::memory_consumption (*base_elements[i].first);
2214  return mem;
2215 }
2216 
2217 
2218 
2219 
2220 // explicit instantiations
2221 #include "fe_system.inst"
2222 
2223 DEAL_II_NAMESPACE_CLOSE
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
void initialize_unit_support_points()
Definition: fe_system.cc:1531
Shape function values.
size_type m() const
Definition: tria.h:67
static const unsigned int invalid_unsigned_int
Definition: types.h:170
void build_interface_constraints()
Definition: fe_system.cc:1270
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe.cc:867
void initialize(const std::vector< const FiniteElement< dim, spacedim > *> &fes, const std::vector< unsigned int > &multiplicities)
Definition: fe_system.cc:1453
void swap(TableBase< N, T > &v)
std::vector< std::pair< unsigned int, unsigned int > > hp_object_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_system.cc:1900
std::vector< std::pair< std_cxx11::shared_ptr< const FiniteElement< dim, spacedim > >, unsigned int > > base_elements
Definition: fe_system.h:962
unsigned int n_nonzero_components(const unsigned int i) const
Definition: fe.h:3023
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe_system.cc:355
virtual bool hp_constraints_are_implemented() const
Definition: fe_system.cc:1664
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_system.cc:2012
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
void build_face_tables(std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > &face_system_to_base_table, std::vector< std::pair< unsigned int, unsigned int > > &face_system_to_component_table, const FiniteElement< dim, spacedim > &finite_element, const bool do_tensor_product=true)
Definition: fe_tools.cc:776
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:923
virtual FiniteElement< dim, spacedim >::InternalDataBase * get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe_system.cc:923
void initialize_quad_dof_index_permutation()
Definition: fe_system.cc:1609
void set_fe_data(const unsigned int base_no, typename FiniteElement< dim, spacedim >::InternalDataBase *)
Definition: fe_system.cc:86
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_system.cc:2151
virtual Point< dim > unit_support_point(const unsigned int index) const
Definition: fe_system.cc:2105
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
InternalData(const unsigned int n_base_elements)
Definition: fe_system.cc:49
void initialize(const unsigned int n_quadrature_points, const FiniteElement< dim, spacedim > &fe, const UpdateFlags flags)
Definition: fe_values.cc:2181
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
Definition: fe_system.cc:482
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
Definition: fe_system.cc:852
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe_system.cc:497
virtual std::string get_name() const
Definition: fe_system.cc:292
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe_system.cc:542
virtual FiniteElement< dim, spacedim > * clone() const
Definition: fe_system.cc:323
virtual std::size_t memory_consumption() const
Definition: fe.cc:1188
size_type n() const
No update.
Third derivatives of shape functions.
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_system.cc:2028
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:313
unsigned int element_multiplicity(const unsigned int index) const
Definition: fe.h:2863
UpdateFlags
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:46
Definition: fe.h:33
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_system.cc:743
void initialize_unit_face_support_points()
Definition: fe_system.cc:1563
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
Definition: fe_system.cc:437
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe.cc:884
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1251
virtual std::string get_name() const =0
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: fe_system.cc:340
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_system.cc:572
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:901
Second derivatives of shape functions.
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:161
void compute_fill(const Mapping< dim, spacedim > &mapping, const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim_1 > &quadrature, const CellSimilarity::Similarity cell_similarity, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_data, const internal::FEValues::MappingRelatedData< dim, spacedim > &mapping_data, internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe_system.cc:1096
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe_system.cc:407
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_system.cc:1788
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe_system.cc:2082
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_system.cc:1678
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual std::size_t memory_consumption() const
Definition: fe_system.cc:2205
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_system.cc:2019
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:912
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_system_to_base_index(const unsigned int index) const
Definition: fe.h:2937
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_system.cc:2038
unsigned int n_components() const
FiniteElement< dim, spacedim >::InternalDataBase & get_fe_data(const unsigned int base_no) const
Definition: fe_system.cc:74
virtual Point< dim-1 > unit_face_support_point(const unsigned int index) const
Definition: fe_system.cc:2128
Shape function gradients.
void build_cell_tables(std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > &system_to_base_table, std::vector< std::pair< unsigned int, unsigned int > > &system_to_component_table, std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > &component_to_base_table, const FiniteElement< dim, spacedim > &finite_element, const bool do_tensor_product=true)
Definition: fe_tools.cc:581
FESystem(const FiniteElement< dim, spacedim > &fe, const unsigned int n_elements)
Definition: fe_system.cc:116
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
Definition: fe_system.cc:527
const unsigned int dofs_per_face
Definition: fe_base.h:288
virtual FiniteElement< dim, spacedim >::InternalDataBase * get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe_system.cc:981
static ::ExceptionBase & ExcNotImplemented()
virtual FiniteElement< dim, spacedim >::InternalDataBase * get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValues::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe_system.cc:867
unsigned int size(const unsigned int i) const
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
Definition: fe_system.cc:807
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_system.cc:2093
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
internal::FEValues::FiniteElementRelatedData< dim, spacedim > & get_fe_output_object(const unsigned int base_no) const
Definition: fe_system.cc:99
unsigned int n_base_elements() const
Definition: fe.h:2853
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe_system.cc:654
virtual ~FESystem()
Definition: fe_system.cc:285
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe_system.cc:452
UpdateFlags update_each
Definition: fe.h:638
Task< RT > new_task(const std_cxx11::function< RT()> &function)
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
Definition: fe_system.cc:392
static ::ExceptionBase & ExcInternalError()