Reference documentation for deal.II version 8.4.2
fe_tools.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2016 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 
17 #include <deal.II/base/quadrature_lib.h>
18 #include <deal.II/base/qprojector.h>
19 #include <deal.II/base/thread_management.h>
20 #include <deal.II/base/utilities.h>
21 #include <deal.II/lac/full_matrix.h>
22 #include <deal.II/lac/householder.h>
23 #include <deal.II/lac/constraint_matrix.h>
24 #include <deal.II/grid/tria.h>
25 #include <deal.II/grid/tria_iterator.h>
26 #include <deal.II/grid/grid_generator.h>
27 #include <deal.II/fe/fe_tools.h>
28 #include <deal.II/fe/fe.h>
29 #include <deal.II/fe/fe_q.h>
30 #include <deal.II/fe/fe_bernstein.h>
31 #include <deal.II/fe/fe_q_hierarchical.h>
32 #include <deal.II/fe/fe_dgq.h>
33 #include <deal.II/fe/fe_dgp.h>
34 #include <deal.II/fe/fe_dgp_monomial.h>
35 #include <deal.II/fe/fe_dgp_nonparametric.h>
36 #include <deal.II/fe/fe_nedelec.h>
37 #include <deal.II/fe/fe_abf.h>
38 #include <deal.II/fe/fe_bdm.h>
39 #include <deal.II/fe/fe_raviart_thomas.h>
40 #include <deal.II/fe/fe_nothing.h>
41 #include <deal.II/fe/fe_system.h>
42 #include <deal.II/fe/fe_values.h>
43 #include <deal.II/fe/mapping_cartesian.h>
44 #include <deal.II/dofs/dof_handler.h>
45 #include <deal.II/dofs/dof_accessor.h>
46 #include <deal.II/dofs/dof_tools.h>
47 #include <deal.II/hp/dof_handler.h>
48 
49 #include <deal.II/base/std_cxx11/shared_ptr.h>
50 
51 #include <deal.II/base/index_set.h>
52 
53 #include <cctype>
54 #include <iostream>
55 
56 
57 DEAL_II_NAMESPACE_OPEN
58 
59 namespace FETools
60 {
61  // Not implemented in the general case.
62  template <class FE>
65  {
66  Assert(false, ExcNotImplemented());
67  return 0;
68  }
69 
70  // Specializations for FE_Q.
71  template <>
73  FEFactory<FE_Q<1, 1> >::get (const Quadrature<1> &quad) const
74  {
75  return new FE_Q<1>(quad);
76  }
77  template <>
79  FEFactory<FE_Q<2, 2> >::get (const Quadrature<1> &quad) const
80  {
81  return new FE_Q<2>(quad);
82  }
83  template <>
85  FEFactory<FE_Q<3, 3> >::get (const Quadrature<1> &quad) const
86  {
87  return new FE_Q<3>(quad);
88  }
89 
90 
91  // Specializations for FE_DGQArbitraryNodes.
92  template <>
94  FEFactory<FE_DGQ<1> >::get (const Quadrature<1> &quad) const
95  {
96  return new FE_DGQArbitraryNodes<1>(quad);
97  }
98  template <>
100  FEFactory<FE_DGQ<2> >::get (const Quadrature<1> &quad) const
101  {
102  return new FE_DGQArbitraryNodes<2>(quad);
103  }
104  template <>
106  FEFactory<FE_DGQ<3> >::get (const Quadrature<1> &quad) const
107  {
108  return new FE_DGQArbitraryNodes<3>(quad);
109  }
110 }
111 
112 namespace
113 {
114  // The following three functions serve to fill the maps from element
115  // names to elements fe_name_map below. The first one exists because
116  // we have finite elements which are not implemented for nonzero
117  // codimension. These should be transferred to the second function
118  // eventually.
119 
120  template <int dim>
121  void
122  fill_no_codim_fe_names (std::map<std::string,std_cxx11::shared_ptr<const Subscriptor> > &result)
123  {
124  typedef std_cxx11::shared_ptr<const Subscriptor> FEFactoryPointer;
125 
126  result["FE_Q_Hierarchical"]
127  = FEFactoryPointer(new FETools::FEFactory<FE_Q_Hierarchical<dim> >);
128  result["FE_ABF"]
129  = FEFactoryPointer(new FETools::FEFactory<FE_ABF<dim> >);
130  result["FE_BDM"]
131  = FEFactoryPointer(new FETools::FEFactory<FE_BDM<dim> >);
132  result["FE_RaviartThomas"]
133  = FEFactoryPointer(new FETools::FEFactory<FE_RaviartThomas<dim> >);
134  result["FE_RaviartThomasNodal"]
135  = FEFactoryPointer(new FETools::FEFactory<FE_RaviartThomasNodal<dim> >);
136  result["FE_Nedelec"]
137  = FEFactoryPointer(new FETools::FEFactory<FE_Nedelec<dim> >);
138  result["FE_DGPNonparametric"]
139  = FEFactoryPointer(new FETools::FEFactory<FE_DGPNonparametric<dim> >);
140  result["FE_DGP"]
141  = FEFactoryPointer(new FETools::FEFactory<FE_DGP<dim> >);
142  result["FE_DGPMonomial"]
143  = FEFactoryPointer(new FETools::FEFactory<FE_DGPMonomial<dim> >);
144  result["FE_DGQ"]
145  = FEFactoryPointer(new FETools::FEFactory<FE_DGQ<dim> >);
146  result["FE_DGQArbitraryNodes"]
147  = FEFactoryPointer(new FETools::FEFactory<FE_DGQ<dim> >);
148  result["FE_Q"]
149  = FEFactoryPointer(new FETools::FEFactory<FE_Q<dim> >);
150  result["FE_Bernstein"]
151  = FEFactoryPointer(new FETools::FEFactory<FE_Bernstein<dim> >);
152  result["FE_Nothing"]
153  = FEFactoryPointer(new FETools::FEFactory<FE_Nothing<dim> >);
154  }
155 
156  // This function fills a map from names to finite elements for any
157  // dimension and codimension for those elements which support
158  // nonzero codimension.
159  template <int dim, int spacedim>
160  void
161  fill_codim_fe_names (std::map<std::string,std_cxx11::shared_ptr<const Subscriptor> > &result)
162  {
163  typedef std_cxx11::shared_ptr<const Subscriptor> FEFactoryPointer;
164 
165  result["FE_DGP"]
166  = FEFactoryPointer(new FETools::FEFactory<FE_DGP<dim,spacedim> >);
167  result["FE_DGQ"]
168  = FEFactoryPointer(new FETools::FEFactory<FE_DGQ<dim,spacedim> >);
169  result["FE_DGQArbitraryNodes"]
170  = FEFactoryPointer(new FETools::FEFactory<FE_DGQ<dim,spacedim> >);
171  result["FE_Q"]
172  = FEFactoryPointer(new FETools::FEFactory<FE_Q<dim,spacedim> >);
173  result["FE_Bernstein"]
174  = FEFactoryPointer(new FETools::FEFactory<FE_Bernstein<dim,spacedim> >);
175  }
176 
177  // The function filling the vector fe_name_map below. It iterates
178  // through all legal dimension/spacedimension pairs and fills
179  // fe_name_map[dimension][spacedimension] with the maps generated
180  // by the functions above.
181  std::vector<std::vector<
182  std::map<std::string,
183  std_cxx11::shared_ptr<const Subscriptor> > > >
184  fill_default_map()
185  {
186  std::vector<std::vector<
187  std::map<std::string,
188  std_cxx11::shared_ptr<const Subscriptor> > > >
189  result(4);
190 
191  for (unsigned int d=0; d<4; ++d)
192  result[d].resize(4);
193 
194  fill_no_codim_fe_names<1> (result[1][1]);
195  fill_no_codim_fe_names<2> (result[2][2]);
196  fill_no_codim_fe_names<3> (result[3][3]);
197 
198  fill_codim_fe_names<1,2> (result[1][2]);
199  fill_codim_fe_names<1,3> (result[1][3]);
200  fill_codim_fe_names<2,3> (result[2][3]);
201 
202  return result;
203  }
204 
205 
206  // have a lock that guarantees that at most one thread is changing
207  // and accessing the fe_name_map variable. make this lock local to
208  // this file.
209  //
210  // this and the next variable are declared static (even though
211  // they're in an anonymous namespace) in order to make icc happy
212  // (which otherwise reports a multiply defined symbol when linking
213  // libraries for more than one space dimension together
214  static
215  Threads::Mutex fe_name_map_lock;
216 
217  // This is the map used by FETools::get_fe_from_name and
218  // FETools::add_fe_name. It is only accessed by functions in this
219  // file, so it is safe to make it a static variable here. It must be
220  // static so that we can link several dimensions together.
221 
222  // The organization of this storage is such that
223  // fe_name_map[dim][spacedim][name] points to an
224  // FEFactoryBase<dim,spacedim> with the name given. Since
225  // all entries of this vector are of different type, we store
226  // pointers to generic objects and cast them when needed.
227 
228  // We use a shared pointer to factory objects, to ensure that they
229  // get deleted at the end of the program run and don't end up as
230  // apparent memory leaks to programs like valgrind.
231 
232  // This vector is initialized at program start time using the
233  // function above. because at this time there are no threads
234  // running, there are no thread-safety issues here. since this is
235  // compiled for all dimensions at once, need to create objects for
236  // each dimension and then separate between them further down
237  static
238  std::vector<std::vector<
239  std::map<std::string,
240  std_cxx11::shared_ptr<const Subscriptor> > > >
241  fe_name_map = fill_default_map();
242 }
243 
244 
245 
246 
247 
248 
249 namespace
250 {
251 
252  // forwarder function for
253  // FE::get_interpolation_matrix. we
254  // will want to call that function
255  // for arbitrary FullMatrix<T>
256  // types, but it only accepts
257  // double arguments. since it is a
258  // virtual function, this can also
259  // not be changed. so have a
260  // forwarder function that calls
261  // that function directly if
262  // T==double, and otherwise uses a
263  // temporary
264  template <int dim, int spacedim>
265  inline
266  void gim_forwarder (const FiniteElement<dim,spacedim> &fe1,
267  const FiniteElement<dim,spacedim> &fe2,
268  FullMatrix<double> &interpolation_matrix)
269  {
270  fe2.get_interpolation_matrix (fe1, interpolation_matrix);
271  }
272 
273 
274  template <int dim, typename number, int spacedim>
275  inline
276  void gim_forwarder (const FiniteElement<dim,spacedim> &fe1,
277  const FiniteElement<dim,spacedim> &fe2,
278  FullMatrix<number> &interpolation_matrix)
279  {
280  FullMatrix<double> tmp (interpolation_matrix.m(),
281  interpolation_matrix.n());
282  fe2.get_interpolation_matrix (fe1, tmp);
283  interpolation_matrix = tmp;
284  }
285 
286 
287 
288  // return how many characters
289  // starting at the given position
290  // of the string match either the
291  // generic string "<dim>" or the
292  // specialized string with "dim"
293  // replaced with the numeric value
294  // of the template argument
295  template <int dim, int spacedim>
296  inline
297  unsigned int match_dimension (const std::string &name,
298  const unsigned int position)
299  {
300  if (position >= name.size())
301  return 0;
302 
303  if ((position+5 < name.size())
304  &&
305  (name[position] == '<')
306  &&
307  (name[position+1] == 'd')
308  &&
309  (name[position+2] == 'i')
310  &&
311  (name[position+3] == 'm')
312  &&
313  (name[position+4] == '>'))
314  return 5;
315 
316  Assert (dim<10, ExcNotImplemented());
317  const char dim_char = '0'+dim;
318 
319  if ((position+3 < name.size())
320  &&
321  (name[position] == '<')
322  &&
323  (name[position+1] == dim_char)
324  &&
325  (name[position+2] == '>'))
326  return 3;
327 
328  // some other string that doesn't
329  // match
330  return 0;
331  }
332 }
333 
334 
335 namespace FETools
336 {
337  template <int dim, int spacedim>
339  {}
340 
341 
342  template<int dim, int spacedim>
344  const FiniteElement<dim,spacedim> &element,
345  std::vector<unsigned int> &renumbering,
346  std::vector<std::vector<unsigned int> > &comp_start)
347  {
348  Assert(renumbering.size() == element.dofs_per_cell,
349  ExcDimensionMismatch(renumbering.size(),
350  element.dofs_per_cell));
351 
352  comp_start.resize(element.n_base_elements());
353 
354  unsigned int k=0;
355  for (unsigned int i=0; i<comp_start.size(); ++i)
356  {
357  comp_start[i].resize(element.element_multiplicity(i));
358  const unsigned int increment
359  = element.base_element(i).dofs_per_cell;
360 
361  for (unsigned int j=0; j<comp_start[i].size(); ++j)
362  {
363  comp_start[i][j] = k;
364  k += increment;
365  }
366  }
367 
368  // For each index i of the
369  // unstructured cellwise
370  // numbering, renumbering
371  // contains the index of the
372  // cell-block numbering
373  for (unsigned int i=0; i<element.dofs_per_cell; ++i)
374  {
375  std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
376  indices = element.system_to_base_index(i);
377  renumbering[i] = comp_start[indices.first.first][indices.first.second]
378  +indices.second;
379  }
380  }
381 
382 
383 
384  template<int dim, int spacedim>
386  const FiniteElement<dim,spacedim> &element,
387  std::vector<types::global_dof_index> &renumbering,
388  std::vector<types::global_dof_index> &block_data,
389  bool return_start_indices)
390  {
391  Assert(renumbering.size() == element.dofs_per_cell,
392  ExcDimensionMismatch(renumbering.size(),
393  element.dofs_per_cell));
394  Assert(block_data.size() == element.n_blocks(),
395  ExcDimensionMismatch(block_data.size(),
396  element.n_blocks()));
397 
399  unsigned int count=0;
400  for (unsigned int b=0; b<element.n_base_elements(); ++b)
401  for (unsigned int m=0; m<element.element_multiplicity(b); ++m)
402  {
403  block_data[count++] = (return_start_indices)
404  ? k
405  : (element.base_element(b).n_dofs_per_cell());
406  k += element.base_element(b).n_dofs_per_cell();
407  }
408  Assert (count == element.n_blocks(), ExcInternalError());
409 
410  std::vector<types::global_dof_index> start_indices(block_data.size());
411  k = 0;
412  for (unsigned int i=0; i<block_data.size(); ++i)
413  if (return_start_indices)
414  start_indices[i] = block_data[i];
415  else
416  {
417  start_indices[i] = k;
418  k += block_data[i];
419  }
420 
421  for (unsigned int i=0; i<element.dofs_per_cell; ++i)
422  {
423  std::pair<unsigned int, types::global_dof_index>
424  indices = element.system_to_block_index(i);
425  renumbering[i] = start_indices[indices.first]
426  +indices.second;
427  }
428  }
429 
430 
431 
432  template <int dim, typename number, int spacedim>
434  const FiniteElement<dim,spacedim> &fe2,
435  FullMatrix<number> &interpolation_matrix)
436  {
437  Assert (fe1.n_components() == fe2.n_components(),
438  ExcDimensionMismatch(fe1.n_components(), fe2.n_components()));
439  Assert(interpolation_matrix.m()==fe2.dofs_per_cell &&
440  interpolation_matrix.n()==fe1.dofs_per_cell,
441  ExcMatrixDimensionMismatch(interpolation_matrix.m(),
442  interpolation_matrix.n(),
443  fe2.dofs_per_cell,
444  fe1.dofs_per_cell));
445 
446  // first try the easy way: maybe
447  // the FE wants to implement things
448  // itself:
449  bool fe_implements_interpolation = true;
450  try
451  {
452  gim_forwarder (fe1, fe2, interpolation_matrix);
453  }
455  {
456  // too bad....
457  fe_implements_interpolation = false;
458  }
459  if (fe_implements_interpolation == true)
460  return;
461 
462  // uh, so this was not the
463  // case. hm. then do it the hard
464  // way. note that this will only
465  // work if the element is
466  // primitive, so check this first
467  Assert (fe1.is_primitive() == true, ExcFENotPrimitive());
468  Assert (fe2.is_primitive() == true, ExcFENotPrimitive());
469 
470  // Initialize FEValues for fe1 at
471  // the unit support points of the
472  // fe2 element.
473  const std::vector<Point<dim> > &
474  fe2_support_points = fe2.get_unit_support_points ();
475 
476  typedef FiniteElement<dim,spacedim> FEL;
477  Assert(fe2_support_points.size()==fe2.dofs_per_cell,
478  typename FEL::ExcFEHasNoSupportPoints());
479 
480  for (unsigned int i=0; i<fe2.dofs_per_cell; ++i)
481  {
482  const unsigned int i1 = fe2.system_to_component_index(i).first;
483  for (unsigned int j=0; j<fe1.dofs_per_cell; ++j)
484  {
485  const unsigned int j1 = fe1.system_to_component_index(j).first;
486  if (i1==j1)
487  interpolation_matrix(i,j) = fe1.shape_value (j,fe2_support_points[i]);
488  else
489  interpolation_matrix(i,j)=0.;
490  }
491  }
492  }
493 
494 
495 
496  template <int dim, typename number, int spacedim>
498  const FiniteElement<dim,spacedim> &fe2,
499  FullMatrix<number> &interpolation_matrix)
500  {
501  Assert (fe1.n_components() == fe2.n_components(),
502  ExcDimensionMismatch(fe1.n_components(), fe2.n_components()));
503  Assert(interpolation_matrix.m()==fe1.dofs_per_cell &&
504  interpolation_matrix.n()==fe1.dofs_per_cell,
505  ExcMatrixDimensionMismatch(interpolation_matrix.m(),
506  interpolation_matrix.n(),
507  fe1.dofs_per_cell,
508  fe1.dofs_per_cell));
509 
510  FullMatrix<number> first_matrix (fe2.dofs_per_cell, fe1.dofs_per_cell);
511  FullMatrix<number> second_matrix(fe1.dofs_per_cell, fe2.dofs_per_cell);
512 
513  get_interpolation_matrix(fe1, fe2, first_matrix);
514  get_interpolation_matrix(fe2, fe1, second_matrix);
515 
516  // int_matrix=second_matrix*first_matrix
517  second_matrix.mmult(interpolation_matrix, first_matrix);
518  }
519 
520 
521 
522  template <int dim, typename number, int spacedim>
524  const FiniteElement<dim,spacedim> &fe2,
525  FullMatrix<number> &difference_matrix)
526  {
527  Assert (fe1.n_components() == fe2.n_components(),
528  ExcDimensionMismatch(fe1.n_components(), fe2.n_components()));
529  Assert(difference_matrix.m()==fe1.dofs_per_cell &&
530  difference_matrix.n()==fe1.dofs_per_cell,
531  ExcMatrixDimensionMismatch(difference_matrix.m(),
532  difference_matrix.n(),
533  fe1.dofs_per_cell,
534  fe1.dofs_per_cell));
535 
536  FullMatrix<number> interpolation_matrix(fe1.dofs_per_cell);
537  get_back_interpolation_matrix(fe1, fe2, interpolation_matrix);
538 
539  for (unsigned int i=0; i<fe1.dofs_per_cell; ++i)
540  difference_matrix(i,i) = 1.;
541 
542  // compute difference
543  difference_matrix.add (-1, interpolation_matrix);
544  }
545 
546 
547 
548  template <int dim, typename number, int spacedim>
550  const FiniteElement<dim,spacedim> &fe2,
551  FullMatrix<number> &matrix)
552  {
553  Assert (fe1.n_components() == 1, ExcNotImplemented());
554  Assert (fe1.n_components() == fe2.n_components(),
555  ExcDimensionMismatch(fe1.n_components(), fe2.n_components()));
556  Assert(matrix.m()==fe2.dofs_per_cell && matrix.n()==fe1.dofs_per_cell,
557  ExcMatrixDimensionMismatch(matrix.m(), matrix.n(),
558  fe2.dofs_per_cell,
559  fe1.dofs_per_cell));
560  matrix = 0;
561 
562  unsigned int n1 = fe1.dofs_per_cell;
563  unsigned int n2 = fe2.dofs_per_cell;
564 
565  // First, create a local mass matrix for
566  // the unit cell
569 
570  // Choose a quadrature rule
571  // Gauss is exact up to degree 2n-1
572  const unsigned int degree = std::max(fe1.tensor_degree(), fe2.tensor_degree());
574  ExcNotImplemented());
575 
576  QGauss<dim> quadrature(degree+1);
577  // Set up FEValues.
579  FEValues<dim> val1 (fe1, quadrature, update_values);
580  val1.reinit (tr.begin_active());
581  FEValues<dim> val2 (fe2, quadrature, flags);
582  val2.reinit (tr.begin_active());
583 
584  // Integrate and invert mass matrix
585  // This happens in the target space
586  FullMatrix<double> mass (n2, n2);
587 
588  for (unsigned int k=0; k<quadrature.size(); ++k)
589  {
590  const double w = val2.JxW(k);
591  for (unsigned int i=0; i<n2; ++i)
592  {
593  const double v = val2.shape_value(i,k);
594  for (unsigned int j=0; j<n2; ++j)
595  mass(i,j) += w*v * val2.shape_value(j,k);
596  }
597  }
598  // Gauss-Jordan should be
599  // sufficient since we expect the
600  // mass matrix to be
601  // well-conditioned
602  mass.gauss_jordan();
603 
604  // Now, test every function of fe1
605  // with test functions of fe2 and
606  // compute the projection of each
607  // unit vector.
608  Vector<double> b(n2);
609  Vector<double> x(n2);
610 
611  for (unsigned int j=0; j<n1; ++j)
612  {
613  b = 0.;
614  for (unsigned int i=0; i<n2; ++i)
615  for (unsigned int k=0; k<quadrature.size(); ++k)
616  {
617  const double w = val2.JxW(k);
618  const double u = val1.shape_value(j,k);
619  const double v = val2.shape_value(i,k);
620  b(i) += u*v*w;
621  }
622 
623  // Multiply by the inverse
624  mass.vmult(x,b);
625  for (unsigned int i=0; i<n2; ++i)
626  matrix(i,j) = x(i);
627  }
628  }
629 
630 
631  template<int dim, int spacedim>
632  void
635  const FiniteElement<dim,spacedim> &fe)
636  {
637  const unsigned int n_dofs = fe.dofs_per_cell;
638  Assert (fe.has_generalized_support_points(), ExcNotInitialized());
639  Assert (N.n()==n_dofs, ExcDimensionMismatch(N.n(), n_dofs));
640  Assert (N.m()==n_dofs, ExcDimensionMismatch(N.m(), n_dofs));
641 
642  const std::vector<Point<dim> > &points = fe.get_generalized_support_points();
643 
644  // We need the values of the
645  // polynomials in all generalized
646  // support points.
647  std::vector<std::vector<double> >
648  values (dim, std::vector<double>(points.size()));
649 
650  // In this vector, we store the
651  // result of the interpolation
652  std::vector<double> local_dofs(n_dofs);
653 
654  // One row per shape
655  // function. Remember that these
656  // are the 'raw' shape functions
657  // where the inverse node matrix is
658  // empty. Otherwise, this would
659  // yield identity.
660  for (unsigned int i=0; i<n_dofs; ++i)
661  {
662  for (unsigned int k=0; k<values[0].size(); ++k)
663  for (unsigned int d=0; d<dim; ++d)
664  values[d][k] = fe.shape_value_component(i,points[k],d);
665  fe.interpolate(local_dofs, values);
666  // Enter the interpolated dofs
667  // into the matrix
668  for (unsigned int j=0; j<n_dofs; ++j)
669  N(j,i) = local_dofs[j];
670  }
671  }
672 
673 
674  /*
675  template<>
676  void
677  compute_embedding_matrices(const FiniteElement<1,2> &,
678  std::vector<std::vector<FullMatrix<double> > > &,
679  const bool)
680  {
681  Assert(false, ExcNotImplemented());
682  }
683 
684 
685  template<>
686  void
687  compute_embedding_matrices(const FiniteElement<1,3> &,
688  std::vector<std::vector<FullMatrix<double> > > &,
689  const bool)
690  {
691  Assert(false, ExcNotImplemented());
692  }
693 
694 
695 
696  template<>
697  void
698  compute_embedding_matrices(const FiniteElement<2,3>&,
699  std::vector<std::vector<FullMatrix<double> > >&,
700  const bool)
701  {
702  Assert(false, ExcNotImplemented());
703  }
704 
705  */
706 
707  namespace
708  {
709  template<int dim, typename number, int spacedim>
710  void
711  compute_embedding_for_shape_function (
712  const unsigned int i,
714  const FEValues<dim, spacedim> &coarse,
715  const Householder<double> &H,
716  FullMatrix<number> &this_matrix,
717  const double threshold)
718  {
719  const unsigned int n = fe.dofs_per_cell;
720  const unsigned int nd = fe.n_components ();
721  const unsigned int nq = coarse.n_quadrature_points;
722 
723  Vector<number> v_coarse(nq*nd);
724  Vector<number> v_fine(n);
725 
726  // The right hand side of
727  // the least squares
728  // problem consists of the
729  // function values of the
730  // coarse grid function in
731  // each quadrature point.
732  if (fe.is_primitive ())
733  {
734  const unsigned int
735  d = fe.system_to_component_index (i).first;
736  const double *phi_i = &coarse.shape_value (i, 0);
737 
738  for (unsigned int k = 0; k < nq; ++k)
739  v_coarse (k * nd + d) = phi_i[k];
740  }
741 
742  else
743  for (unsigned int d = 0; d < nd; ++d)
744  for (unsigned int k = 0; k < nq; ++k)
745  v_coarse (k * nd + d) = coarse.shape_value_component (i, k, d);
746 
747  // solve the least squares
748  // problem.
749  const double result = H.least_squares (v_fine, v_coarse);
750  Assert (result <= threshold, ExcLeastSquaresError (result));
751  // Avoid warnings in release mode
752  (void)result;
753  (void)threshold;
754 
755  // Copy into the result
756  // matrix. Since the matrix
757  // maps a coarse grid
758  // function to a fine grid
759  // function, the columns
760  // are fine grid.
761  for (unsigned int j = 0; j < n; ++j)
762  this_matrix(j, i) = v_fine(j);
763  }
764 
765 
766  template<int dim, typename number, int spacedim>
767  void
768  compute_embedding_matrices_for_refinement_case (
770  std::vector<FullMatrix<number> > &matrices,
771  const unsigned int ref_case,
772  const double threshold)
773  {
774  const unsigned int n = fe.dofs_per_cell;
775  const unsigned int nc = GeometryInfo<dim>::n_children(RefinementCase<dim>(ref_case));
776  for (unsigned int i = 0; i < nc; ++i)
777  {
778  Assert(matrices[i].n() == n, ExcDimensionMismatch(matrices[i].n (), n));
779  Assert(matrices[i].m() == n, ExcDimensionMismatch(matrices[i].m (), n));
780  }
781 
782  // Set up meshes, one with a single
783  // reference cell and refine it once
785  GridGenerator::hyper_cube (tria, 0, 1);
786  tria.begin_active()->set_refine_flag (RefinementCase<dim>(ref_case));
788 
789  const unsigned int degree = fe.degree;
790  QGauss<dim> q_fine (degree+1);
791  const unsigned int nq = q_fine.size();
792 
793  FEValues<dim,spacedim> fine (fe, q_fine,
796  update_values);
797 
798  // We search for the polynomial on
799  // the small cell, being equal to
800  // the coarse polynomial in all
801  // quadrature points.
802 
803  // First build the matrix for this
804  // least squares problem. This
805  // contains the values of the fine
806  // cell polynomials in the fine
807  // cell grid points.
808 
809  // This matrix is the same for all
810  // children.
811  fine.reinit (tria.begin_active ());
812  const unsigned int nd = fe.n_components ();
813  FullMatrix<number> A (nq*nd, n);
814 
815  for (unsigned int j = 0; j < n; ++j)
816  for (unsigned int d = 0; d < nd; ++d)
817  for (unsigned int k = 0; k < nq; ++k)
818  A (k * nd + d, j) = fine.shape_value_component (j, k, d);
819 
820  Householder<double> H (A);
821  unsigned int cell_number = 0;
822 
823  Threads::TaskGroup<void> task_group;
824 
826  fine_cell = tria.begin_active (); fine_cell != tria.end ();
827  ++fine_cell, ++cell_number)
828  {
829  fine.reinit (fine_cell);
830 
831  // evaluate on the coarse cell (which
832  // is the first -- inactive -- cell on
833  // the lowest level of the
834  // triangulation we have created)
835  const std::vector<Point<spacedim> > &q_points_fine = fine.get_quadrature_points();
836  std::vector<Point<dim> > q_points_coarse(q_points_fine.size());
837  for (unsigned int i=0; i<q_points_fine.size(); ++i)
838  for (unsigned int j=0; j<dim; ++j)
839  q_points_coarse[i](j) = q_points_fine[i](j);
840  const Quadrature<dim> q_coarse (q_points_coarse,
841  fine.get_JxW_values ());
842  FEValues<dim,spacedim> coarse (fe, q_coarse, update_values);
843 
844  coarse.reinit (tria.begin (0));
845 
846  FullMatrix<double> &this_matrix = matrices[cell_number];
847 
848  // Compute this once for each
849  // coarse grid basis function. can
850  // spawn subtasks if n is
851  // sufficiently large so that there
852  // are more than about 5000
853  // operations in the inner loop
854  // (which is basically const * n^2
855  // operations).
856  if (n > 30)
857  {
858  for (unsigned int i = 0; i < n; ++i)
859  {
860  task_group +=
861  Threads::new_task (&compute_embedding_for_shape_function<dim, number, spacedim>,
862  i, fe, coarse, H, this_matrix, threshold);
863  }
864  task_group.join_all();
865  }
866  else
867  {
868  for (unsigned int i = 0; i < n; ++i)
869  {
870  compute_embedding_for_shape_function<dim, number, spacedim>
871  (i, fe, coarse, H, this_matrix, threshold);
872  }
873  }
874 
875  // Remove small entries from
876  // the matrix
877  for (unsigned int i = 0; i < this_matrix.m (); ++i)
878  for (unsigned int j = 0; j < this_matrix.n (); ++j)
879  if (std::fabs (this_matrix (i, j)) < 1e-12)
880  this_matrix (i, j) = 0.;
881  }
882 
883  Assert (cell_number == GeometryInfo<dim>::n_children (RefinementCase<dim> (ref_case)),
884  ExcInternalError ());
885  }
886  }
887 
888 
889  template <int dim, typename number, int spacedim>
890  void
892  std::vector<std::vector<FullMatrix<number> > > &matrices,
893  const bool isotropic_only,
894  const double threshold)
895  {
896  Threads::TaskGroup<void> task_group;
897 
898  // loop over all possible refinement cases
899  unsigned int ref_case = (isotropic_only)
902 
903  for (; ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
904  task_group += Threads::new_task (&compute_embedding_matrices_for_refinement_case<dim, number, spacedim>,
905  fe, matrices[ref_case-1], ref_case, threshold);
906 
907  task_group.join_all ();
908  }
909 
910 
911 
912  template <int dim, typename number, int spacedim>
913  void
916  const unsigned int face_coarse,
917  const unsigned int face_fine,
918  const double threshold)
919  {
920  Assert(face_coarse==0, ExcNotImplemented());
921  Assert(face_fine==0, ExcNotImplemented());
922 
923  const unsigned int nc = GeometryInfo<dim>::max_children_per_face;
924  const unsigned int n = fe.dofs_per_face;
925  const unsigned int nd = fe.n_components();
926  const unsigned int degree = fe.degree;
927 
928  const bool normal = fe.conforms(FiniteElementData<dim>::Hdiv);
929  const bool tangential = fe.conforms(FiniteElementData<dim>::Hcurl);
930 
931  for (unsigned int i=0; i<nc; ++i)
932  {
933  Assert(matrices[i].n() == n, ExcDimensionMismatch(matrices[i].n(),n));
934  Assert(matrices[i].m() == n, ExcDimensionMismatch(matrices[i].m(),n));
935  }
936 
937  // In order to make the loops below
938  // simpler, we introduce vectors
939  // containing for indices 0-n the
940  // number of the corresponding
941  // shape value on the cell.
942  std::vector<unsigned int> face_c_dofs(n);
943  std::vector<unsigned int> face_f_dofs(n);
944  {
945  unsigned int face_dof=0;
946  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_face; ++i)
947  {
948  const unsigned int offset_c = GeometryInfo<dim>::face_to_cell_vertices(face_coarse, i)
949  *fe.dofs_per_vertex;
950  const unsigned int offset_f = GeometryInfo<dim>::face_to_cell_vertices(face_fine, i)
951  *fe.dofs_per_vertex;
952  for (unsigned int j=0; j<fe.dofs_per_vertex; ++j)
953  {
954  face_c_dofs[face_dof] = offset_c + j;
955  face_f_dofs[face_dof] = offset_f + j;
956  ++face_dof;
957  }
958  }
959  for (unsigned int i=1; i<=GeometryInfo<dim>::lines_per_face; ++i)
960  {
961  const unsigned int offset_c = fe.first_line_index
962  + GeometryInfo<dim>::face_to_cell_lines(face_coarse, i-1)
963  *fe.dofs_per_line;
964  const unsigned int offset_f = fe.first_line_index
965  + GeometryInfo<dim>::face_to_cell_lines(face_fine, i-1)
966  *fe.dofs_per_line;
967  for (unsigned int j=0; j<fe.dofs_per_line; ++j)
968  {
969  face_c_dofs[face_dof] = offset_c + j;
970  face_f_dofs[face_dof] = offset_f + j;
971  ++face_dof;
972  }
973  }
974  for (unsigned int i=1; i<=GeometryInfo<dim>::quads_per_face; ++i)
975  {
976  const unsigned int offset_c = fe.first_quad_index
977  + face_coarse
978  *fe.dofs_per_quad;
979  const unsigned int offset_f = fe.first_quad_index
980  + face_fine
981  *fe.dofs_per_quad;
982  for (unsigned int j=0; j<fe.dofs_per_quad; ++j)
983  {
984  face_c_dofs[face_dof] = offset_c + j;
985  face_f_dofs[face_dof] = offset_f + j;
986  ++face_dof;
987  }
988  }
989  Assert (face_dof == fe.dofs_per_face, ExcInternalError());
990  }
991 
992  // Set up meshes, one with a single
993  // reference cell and refine it once
995  GridGenerator::hyper_cube (tria, 0, 1);
996  tria.refine_global(1);
997  MappingCartesian<dim> mapping;
998 
999  // Setup quadrature and FEValues
1000  // for a face. We cannot use
1001  // FEFaceValues and
1002  // FESubfaceValues because of
1003  // some nifty handling of
1004  // refinement cases. Guido stops
1005  // disliking and instead starts
1006  // hating the anisotropic implementation
1007  QGauss<dim-1> q_gauss(degree+1);
1008  const Quadrature<dim> q_fine = QProjector<dim>::project_to_face(q_gauss, face_fine);
1009  const unsigned int nq = q_fine.size();
1010 
1011  FEValues<dim> fine (mapping, fe, q_fine,
1013 
1014  // We search for the polynomial on
1015  // the small cell, being equal to
1016  // the coarse polynomial in all
1017  // quadrature points.
1018 
1019  // First build the matrix for this
1020  // least squares problem. This
1021  // contains the values of the fine
1022  // cell polynomials in the fine
1023  // cell grid points.
1024 
1025  // This matrix is the same for all
1026  // children.
1027  fine.reinit(tria.begin_active());
1028  FullMatrix<number> A(nq*nd, n);
1029  for (unsigned int j=0; j<n; ++j)
1030  for (unsigned int k=0; k<nq; ++k)
1031  if (nd != dim)
1032  for (unsigned int d=0; d<nd; ++d)
1033  A(k*nd+d,j) = fine.shape_value_component(face_f_dofs[j],k,d);
1034  else
1035  {
1036  if (normal)
1037  A(k*nd,j) = fine.shape_value_component(face_f_dofs[j],k,0);
1038  if (tangential)
1039  for (unsigned int d=1; d<dim; ++d)
1040  A(k*nd+d,j) = fine.shape_value_component(face_f_dofs[j],k,d);
1041  }
1042 
1043  Householder<double> H(A);
1044 
1045  Vector<number> v_coarse(nq*nd);
1046  Vector<number> v_fine(n);
1047 
1048 
1049 
1050  for (unsigned int cell_number = 0; cell_number < GeometryInfo<dim>::max_children_per_face;
1051  ++cell_number)
1052  {
1053  const Quadrature<dim> q_coarse
1054  = QProjector<dim>::project_to_subface(q_gauss, face_coarse, cell_number);
1055  FEValues<dim> coarse (mapping, fe, q_coarse, update_values);
1056 
1059  tria.begin(0)->refinement_case(), face_coarse, cell_number));
1060  fine.reinit(fine_cell);
1061  coarse.reinit(tria.begin(0));
1062 
1063  FullMatrix<double> &this_matrix = matrices[cell_number];
1064 
1065  // Compute this once for each
1066  // coarse grid basis function
1067  for (unsigned int i=0; i<n; ++i)
1068  {
1069  // The right hand side of
1070  // the least squares
1071  // problem consists of the
1072  // function values of the
1073  // coarse grid function in
1074  // each quadrature point.
1075  for (unsigned int k=0; k<nq; ++k)
1076  if (nd != dim)
1077  for (unsigned int d=0; d<nd; ++d)
1078  v_coarse(k*nd+d) = coarse.shape_value_component (face_c_dofs[i],k,d);
1079  else
1080  {
1081  if (normal)
1082  v_coarse(k*nd) = coarse.shape_value_component(face_c_dofs[i],k,0);
1083  if (tangential)
1084  for (unsigned int d=1; d<dim; ++d)
1085  v_coarse(k*nd+d) = coarse.shape_value_component(face_c_dofs[i],k,d);
1086  }
1087  // solve the least squares
1088  // problem.
1089  const double result = H.least_squares(v_fine, v_coarse);
1090  Assert (result <= threshold, ExcLeastSquaresError(result));
1091  // Avoid compiler warnings in Release mode
1092  (void)result;
1093  (void)threshold;
1094 
1095  // Copy into the result
1096  // matrix. Since the matrix
1097  // maps a coarse grid
1098  // function to a fine grid
1099  // function, the columns
1100  // are fine grid.
1101  for (unsigned int j=0; j<n; ++j)
1102  this_matrix(j,i) = v_fine(j);
1103  }
1104  // Remove small entries from
1105  // the matrix
1106  for (unsigned int i=0; i<this_matrix.m(); ++i)
1107  for (unsigned int j=0; j<this_matrix.n(); ++j)
1108  if (std::fabs(this_matrix(i,j)) < 1e-12)
1109  this_matrix(i,j) = 0.;
1110  }
1111  }
1112 
1113 
1114 
1115  template <int dim, typename number, int spacedim>
1116  void
1118  std::vector<std::vector<FullMatrix<number> > > &matrices,
1119  const bool isotropic_only)
1120  {
1121  const unsigned int n = fe.dofs_per_cell;
1122  const unsigned int nd = fe.n_components();
1123  const unsigned int degree = fe.degree;
1124 
1125  // prepare FEValues, quadrature etc on
1126  // coarse cell
1127  QGauss<dim> q_fine(degree+1);
1128  const unsigned int nq = q_fine.size();
1129 
1130  // create mass matrix on coarse cell.
1131  FullMatrix<number> mass(n, n);
1132  {
1133  // set up a triangulation for coarse cell
1135  GridGenerator::hyper_cube (tr, 0, 1);
1136 
1137  FEValues<dim,spacedim> coarse (fe, q_fine,
1139 
1140  typename Triangulation<dim,spacedim>::cell_iterator coarse_cell
1141  = tr.begin(0);
1142  coarse.reinit (coarse_cell);
1143 
1144  const std::vector<double> &JxW = coarse.get_JxW_values();
1145  for (unsigned int i=0; i<n; ++i)
1146  for (unsigned int j=0; j<n; ++j)
1147  if (fe.is_primitive())
1148  {
1149  const double *coarse_i = &coarse.shape_value(i,0);
1150  const double *coarse_j = &coarse.shape_value(j,0);
1151  double mass_ij = 0;
1152  for (unsigned int k=0; k<nq; ++k)
1153  mass_ij += JxW[k] * coarse_i[k] * coarse_j[k];
1154  mass(i,j) = mass_ij;
1155  }
1156  else
1157  {
1158  double mass_ij = 0;
1159  for (unsigned int d=0; d<nd; ++d)
1160  for (unsigned int k=0; k<nq; ++k)
1161  mass_ij += JxW[k] * coarse.shape_value_component(i,k,d)
1162  * coarse.shape_value_component(j,k,d);
1163  mass(i,j) = mass_ij;
1164  }
1165 
1166  // invert mass matrix
1167  mass.gauss_jordan();
1168  }
1169 
1170  // loop over all possible
1171  // refinement cases
1172  unsigned int ref_case = (isotropic_only)
1175  for (; ref_case <= RefinementCase<dim>::isotropic_refinement; ++ref_case)
1176  {
1177  const unsigned int
1179 
1180  for (unsigned int i=0; i<nc; ++i)
1181  {
1182  Assert(matrices[ref_case-1][i].n() == n,
1183  ExcDimensionMismatch(matrices[ref_case-1][i].n(),n));
1184  Assert(matrices[ref_case-1][i].m() == n,
1185  ExcDimensionMismatch(matrices[ref_case-1][i].m(),n));
1186  }
1187 
1188  // create a respective refinement on the
1189  // triangulation
1191  GridGenerator::hyper_cube (tr, 0, 1);
1192  tr.begin_active()->set_refine_flag(RefinementCase<dim>(ref_case));
1194 
1197  update_values);
1198 
1199  typename Triangulation<dim,spacedim>::cell_iterator coarse_cell
1200  = tr.begin(0);
1201 
1202  Vector<number> v_coarse(n);
1203  Vector<number> v_fine(n);
1204 
1205  for (unsigned int cell_number=0; cell_number<nc; ++cell_number)
1206  {
1207  FullMatrix<double> &this_matrix = matrices[ref_case-1][cell_number];
1208 
1209  // Compute right hand side,
1210  // which is a fine level basis
1211  // function tested with the
1212  // coarse level functions.
1213  fine.reinit(coarse_cell->child(cell_number));
1214  const std::vector<Point<spacedim> > &q_points_fine = fine.get_quadrature_points();
1215  std::vector<Point<dim> > q_points_coarse(q_points_fine.size());
1216  for (unsigned int q=0; q<q_points_fine.size(); ++q)
1217  for (unsigned int j=0; j<dim; ++j)
1218  q_points_coarse[q](j) = q_points_fine[q](j);
1219  Quadrature<dim> q_coarse (q_points_coarse,
1220  fine.get_JxW_values());
1222  coarse.reinit(coarse_cell);
1223 
1224  // Build RHS
1225 
1226  const std::vector<double> &JxW = fine.get_JxW_values();
1227 
1228  // Outer loop over all fine
1229  // grid shape functions phi_j
1230  for (unsigned int j=0; j<fe.dofs_per_cell; ++j)
1231  {
1232  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
1233  {
1234  if (fe.is_primitive())
1235  {
1236  const double *coarse_i = &coarse.shape_value(i,0);
1237  const double *fine_j = &fine.shape_value(j,0);
1238 
1239  double update = 0;
1240  for (unsigned int k=0; k<nq; ++k)
1241  update += JxW[k] * coarse_i[k] * fine_j[k];
1242  v_fine(i) = update;
1243  }
1244  else
1245  {
1246  double update = 0;
1247  for (unsigned int d=0; d<nd; ++d)
1248  for (unsigned int k=0; k<nq; ++k)
1249  update += JxW[k] * coarse.shape_value_component(i,k,d)
1250  * fine.shape_value_component(j,k,d);
1251  v_fine(i) = update;
1252  }
1253  }
1254 
1255  // RHS ready. Solve system
1256  // and enter row into
1257  // matrix
1258  mass.vmult (v_coarse, v_fine);
1259  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
1260  this_matrix(i,j) = v_coarse(i);
1261  }
1262 
1263  // Remove small entries from
1264  // the matrix
1265  for (unsigned int i=0; i<this_matrix.m(); ++i)
1266  for (unsigned int j=0; j<this_matrix.n(); ++j)
1267  if (std::fabs(this_matrix(i,j)) < 1e-12)
1268  this_matrix(i,j) = 0.;
1269  }
1270  }
1271  }
1272 
1273 
1274 
1275  template <int dim, int spacedim>
1276  void
1277  add_fe_name(const std::string &parameter_name,
1278  const FEFactoryBase<dim,spacedim> *factory)
1279  {
1280  // Erase everything after the
1281  // actual class name
1282  std::string name = parameter_name;
1283  unsigned int name_end =
1284  name.find_first_not_of(std::string("abcdefghijklmnopqrstuvwxyzABCDEFGHIJKLMNOPQRSTUVWXYZ0123456789_"));
1285  if (name_end < name.size())
1286  name.erase(name_end);
1287  // first make sure that no other
1288  // thread intercepts the
1289  // operation of this function;
1290  // for this, acquire the lock
1291  // until we quit this function
1292  Threads::Mutex::ScopedLock lock(fe_name_map_lock);
1293 
1294  Assert(fe_name_map[dim][spacedim].find(name) == fe_name_map[dim][spacedim].end(),
1295  ExcMessage("Cannot change existing element in finite element name list"));
1296 
1297  // Insert the normalized name into
1298  // the map
1299  fe_name_map[dim][spacedim][name] =
1300  std_cxx11::shared_ptr<const Subscriptor> (factory);
1301  }
1302 
1303 
1304  namespace internal
1305  {
1306  namespace
1307  {
1308  // TODO: this encapsulates the call to the
1309  // dimension-dependent fe_name_map so that we
1310  // have a unique interface. could be done
1311  // smarter?
1312  template <int dim, int spacedim>
1314  get_fe_from_name_ext (std::string &name,
1315  const std::map<std::string,
1316  std_cxx11::shared_ptr<const Subscriptor> >
1317  &fe_name_map)
1318  {
1319  // Extract the name of the
1320  // finite element class, which only
1321  // contains characters, numbers and
1322  // underscores.
1323  unsigned int name_end =
1324  name.find_first_not_of(std::string("abcdefghijklmnopqrstuvwxyzABCDEFGHIJKLMNOPQRSTUVWXYZ0123456789_"));
1325  const std::string name_part(name, 0, name_end);
1326  name.erase(0, name_part.size());
1327 
1328  // now things get a little more
1329  // complicated: FESystem. it's
1330  // more complicated, since we
1331  // have to figure out what the
1332  // base elements are. this can
1333  // only be done recursively
1334  if (name_part == "FESystem")
1335  {
1336  // next we have to get at the
1337  // base elements. start with
1338  // the first. wrap the whole
1339  // block into try-catch to
1340  // make sure we destroy the
1341  // pointers we got from
1342  // recursive calls if one of
1343  // these calls should throw
1344  // an exception
1345  std::vector<const FiniteElement<dim,spacedim>*> base_fes;
1346  std::vector<unsigned int> base_multiplicities;
1347  try
1348  {
1349  // Now, just the [...]
1350  // part should be left.
1351  if (name.size() == 0 || name[0] != '[')
1352  throw (std::string("Invalid first character in ") + name);
1353  do
1354  {
1355  // Erase the
1356  // leading '[' or '-'
1357  name.erase(0,1);
1358  // Now, the name of the
1359  // first base element is
1360  // first... Let's get it
1361  base_fes.push_back (get_fe_from_name_ext<dim,spacedim> (name,
1362  fe_name_map));
1363  // next check whether
1364  // FESystem placed a
1365  // multiplicity after
1366  // the element name
1367  if (name[0] == '^')
1368  {
1369  // yes. Delete the '^'
1370  // and read this
1371  // multiplicity
1372  name.erase(0,1);
1373 
1374  const std::pair<int,unsigned int> tmp
1376  name.erase(0, tmp.second);
1377  // add to length,
1378  // including the '^'
1379  base_multiplicities.push_back (tmp.first);
1380  }
1381  else
1382  // no, so
1383  // multiplicity is
1384  // 1
1385  base_multiplicities.push_back (1);
1386 
1387  // so that's it for
1388  // this base
1389  // element. base
1390  // elements are
1391  // separated by '-',
1392  // and the list is
1393  // terminated by ']',
1394  // so loop while the
1395  // next character is
1396  // '-'
1397  }
1398  while (name[0] == '-');
1399 
1400  // so we got to the end
1401  // of the '-' separated
1402  // list. make sure that
1403  // we actually had a ']'
1404  // there
1405  if (name.size() == 0 || name[0] != ']')
1406  throw (std::string("Invalid first character in ") + name);
1407  name.erase(0,1);
1408  // just one more sanity check
1409  Assert ((base_fes.size() == base_multiplicities.size())
1410  &&
1411  (base_fes.size() > 0),
1412  ExcInternalError());
1413 
1414  // ok, apparently
1415  // everything went ok. so
1416  // generate the composed
1417  // element
1418  FiniteElement<dim,spacedim> *system_element = 0;
1419 
1420  // uses new FESystem constructor
1421  // which is independent of
1422  // the number of FEs in the system
1423  system_element = new FESystem<dim,spacedim>(base_fes, base_multiplicities);
1424 
1425  // now we don't need the
1426  // list of base elements
1427  // any more
1428  for (unsigned int i=0; i<base_fes.size(); ++i)
1429  delete base_fes[i];
1430 
1431  // finally return our
1432  // findings
1433  // Add the closing ']' to
1434  // the length
1435  return system_element;
1436 
1437  }
1438  catch (...)
1439  {
1440  // ups, some exception
1441  // was thrown. prevent a
1442  // memory leak, and then
1443  // pass on the exception
1444  // to the caller
1445  for (unsigned int i=0; i<base_fes.size(); ++i)
1446  delete base_fes[i];
1447  throw;
1448  }
1449 
1450  // this is a place where we
1451  // should really never get,
1452  // since above we have either
1453  // returned from the
1454  // try-clause, or have
1455  // re-thrown in the catch
1456  // clause. check that we
1457  // never get here
1458  Assert (false, ExcInternalError());
1459  }
1460  else if (name_part == "FE_Nothing")
1461  {
1462  // remove the () from FE_Nothing()
1463  name.erase(0,2);
1464 
1465  // this is a bit of a hack, as
1466  // FE_Nothing does not take a
1467  // degree, but it does take an
1468  // argument, which defaults to 1,
1469  // so this properly returns
1470  // FE_Nothing()
1471  const Subscriptor *ptr = fe_name_map.find(name_part)->second.get();
1472  const FEFactoryBase<dim,spacedim> *fef=dynamic_cast<const FEFactoryBase<dim,spacedim>*>(ptr);
1473  return fef->get(1);
1474  }
1475  else
1476  {
1477  // Make sure no other thread
1478  // is just adding an element
1479  Threads::Mutex::ScopedLock lock (fe_name_map_lock);
1480  AssertThrow (fe_name_map.find(name_part) != fe_name_map.end(),
1481  ExcInvalidFEName(name));
1482 
1483  // Now, just the (degree)
1484  // or (Quadrature<1>(degree+1))
1485  // part should be left.
1486  if (name.size() == 0 || name[0] != '(')
1487  throw (std::string("Invalid first character in ") + name);
1488  name.erase(0,1);
1489  if (name[0] != 'Q')
1490  {
1491  const std::pair<int,unsigned int> tmp
1493  name.erase(0, tmp.second+1);
1494  const Subscriptor *ptr = fe_name_map.find(name_part)->second.get();
1495  const FEFactoryBase<dim,spacedim> *fef=dynamic_cast<const FEFactoryBase<dim,spacedim>*>(ptr);
1496  return fef->get(tmp.first);
1497  }
1498  else
1499  {
1500  unsigned int position = name.find('(');
1501  const std::string quadrature_name(name, 0, position);
1502  name.erase(0,position+1);
1503  if (quadrature_name.compare("QGaussLobatto") == 0)
1504  {
1505  const std::pair<int,unsigned int> tmp
1507  // delete "))"
1508  name.erase(0, tmp.second+2);
1509  const Subscriptor *ptr = fe_name_map.find(name_part)->second.get();
1510  const FEFactoryBase<dim,spacedim> *fef=dynamic_cast<const FEFactoryBase<dim,spacedim>*>(ptr);
1511  return fef->get(QGaussLobatto<1>(tmp.first));
1512  }
1513  else
1514  {
1515  AssertThrow (false,ExcNotImplemented());
1516  }
1517  }
1518  }
1519 
1520 
1521  // hm, if we have come thus far, we
1522  // didn't know what to do with the
1523  // string we got. so do as the docs
1524  // say: raise an exception
1525  AssertThrow (false, ExcInvalidFEName(name));
1526 
1527  // make some compilers happy that
1528  // do not realize that we can't get
1529  // here after throwing
1530  return 0;
1531  }
1532 
1533 
1534 
1535  template <int dim,int spacedim>
1536  FiniteElement<dim,spacedim> *get_fe_from_name (std::string &name)
1537  {
1538  return get_fe_from_name_ext<dim,spacedim> (name, fe_name_map[dim][spacedim]);
1539  }
1540  }
1541  }
1542 
1543 
1544 
1545 
1546 
1547  template <int dim, int spacedim>
1549  get_fe_by_name (const std::string &parameter_name)
1550  {
1551  std::string name = Utilities::trim(parameter_name);
1552  std::size_t index = 1;
1553  // remove spaces that are not between two word (things that match the
1554  // regular expression [A-Za-z0-9_]) characters.
1555  while (2 < name.size() && index < name.size() - 1)
1556  {
1557  if (name[index] == ' ' &&
1558  (!(std::isalnum(name[index - 1]) || name[index - 1] == '_') ||
1559  !(std::isalnum(name[index + 1]) || name[index + 1] == '_')))
1560  {
1561  name.erase(index, 1);
1562  }
1563  else
1564  {
1565  ++index;
1566  }
1567  }
1568 
1569  // Create a version of the name
1570  // string where all template
1571  // parameters are eliminated.
1572  for (unsigned int pos1 = name.find('<');
1573  pos1 < name.size();
1574  pos1 = name.find('<'))
1575  {
1576 
1577  const unsigned int pos2 = name.find('>');
1578  // If there is only a single
1579  // character between those two,
1580  // it should be 'd' or the number
1581  // representing the dimension.
1582  if (pos2-pos1 == 2)
1583  {
1584  const char dimchar = '0' + dim;
1585  (void)dimchar;
1586  if (name.at(pos1+1) != 'd')
1587  Assert (name.at(pos1+1) == dimchar,
1588  ExcInvalidFEDimension(name.at(pos1+1), dim));
1589  }
1590  else
1591  Assert(pos2-pos1 == 4, ExcInvalidFEName(name));
1592 
1593  // If pos1==pos2, then we are
1594  // probably at the end of the
1595  // string
1596  if (pos2 != pos1)
1597  name.erase(pos1, pos2-pos1+1);
1598  }
1599  // Replace all occurrences of "^dim"
1600  // by "^d" to be handled by the
1601  // next loop
1602  for (unsigned int pos = name.find("^dim");
1603  pos < name.size();
1604  pos = name.find("^dim"))
1605  name.erase(pos+2, 2);
1606 
1607  // Replace all occurrences of "^d"
1608  // by using the actual dimension
1609  for (unsigned int pos = name.find("^d");
1610  pos < name.size();
1611  pos = name.find("^d"))
1612  name.at(pos+1) = '0' + dim;
1613 
1614  try
1615  {
1616  FiniteElement<dim,spacedim> *fe = internal::get_fe_from_name<dim,spacedim> (name);
1617 
1618  // Make sure the auxiliary function
1619  // ate up all characters of the name.
1620  AssertThrow (name.size() == 0,
1621  ExcInvalidFEName(parameter_name
1622  + std::string(" extra characters after "
1623  "end of name")));
1624  return fe;
1625  }
1626  catch (const std::string &errline)
1627  {
1628  AssertThrow(false, ExcInvalidFEName(parameter_name
1629  + std::string(" at ")
1630  + errline));
1631  return 0;
1632  }
1633  }
1634 
1635 
1636  template <int dim>
1638  get_fe_from_name (const std::string &parameter_name)
1639  {
1640  return get_fe_by_name<dim,dim> (parameter_name);
1641  }
1642 
1643 
1644  template <int dim, int spacedim>
1645  void
1646 
1648  const Quadrature<dim> &lhs_quadrature,
1649  const Quadrature<dim> &rhs_quadrature,
1650  FullMatrix<double> &X)
1651  {
1652  Assert (fe.n_components() == 1, ExcNotImplemented());
1653 
1654  // first build the matrices M and Q
1655  // described in the documentation
1657  FullMatrix<double> Q (fe.dofs_per_cell, rhs_quadrature.size());
1658 
1659  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
1660  for (unsigned int j=0; j<fe.dofs_per_cell; ++j)
1661  for (unsigned int q=0; q<lhs_quadrature.size(); ++q)
1662  M(i,j) += fe.shape_value (i, lhs_quadrature.point(q)) *
1663  fe.shape_value (j, lhs_quadrature.point(q)) *
1664  lhs_quadrature.weight(q);
1665 
1666  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
1667  for (unsigned int q=0; q<rhs_quadrature.size(); ++q)
1668  Q(i,q) += fe.shape_value (i, rhs_quadrature.point(q)) *
1669  rhs_quadrature.weight(q);
1670 
1671  // then invert M
1672  FullMatrix<double> M_inverse (fe.dofs_per_cell, fe.dofs_per_cell);
1673  M_inverse.invert (M);
1674 
1675  // finally compute the result
1676  X.reinit (fe.dofs_per_cell, rhs_quadrature.size());
1677  M_inverse.mmult (X, Q);
1678 
1679  Assert (X.m() == fe.dofs_per_cell, ExcInternalError());
1680  Assert (X.n() == rhs_quadrature.size(), ExcInternalError());
1681  }
1682 
1683 
1684 
1685  template <int dim, int spacedim>
1686  void
1688  const Quadrature<dim> &quadrature,
1689  FullMatrix<double> &I_q)
1690  {
1691  Assert (fe.n_components() == 1, ExcNotImplemented());
1692  Assert (I_q.m() == quadrature.size(),
1693  ExcMessage ("Wrong matrix size"));
1694  Assert (I_q.n() == fe.dofs_per_cell, ExcMessage ("Wrong matrix size"));
1695 
1696  for (unsigned int q=0; q<quadrature.size(); ++q)
1697  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
1698  I_q(q,i) = fe.shape_value (i, quadrature.point(q));
1699  }
1700 
1701 
1702 
1703  template <int dim>
1704  void
1706  const FullMatrix<double> &projection_matrix,
1707  const std::vector< Tensor<1, dim > > &vector_of_tensors_at_qp,
1708  std::vector< Tensor<1, dim > > &vector_of_tensors_at_nodes)
1709  {
1710 
1711  // check that the number columns of the projection_matrix
1712  // matches the size of the vector_of_tensors_at_qp
1713  Assert(projection_matrix.n_cols() == vector_of_tensors_at_qp.size(),
1714  ExcDimensionMismatch(projection_matrix.n_cols(),
1715  vector_of_tensors_at_qp.size()));
1716 
1717  // check that the number rows of the projection_matrix
1718  // matches the size of the vector_of_tensors_at_nodes
1719  Assert(projection_matrix.n_rows() == vector_of_tensors_at_nodes.size(),
1720  ExcDimensionMismatch(projection_matrix.n_rows(),
1721  vector_of_tensors_at_nodes.size()));
1722 
1723  // number of support points (nodes) to project to
1724  const unsigned int n_support_points = projection_matrix.n_rows();
1725  // number of quadrature points to project from
1726  const unsigned int n_quad_points = projection_matrix.n_cols();
1727 
1728  // component projected to the nodes
1729  Vector<double> component_at_node(n_support_points);
1730  // component at the quadrature point
1731  Vector<double> component_at_qp(n_quad_points);
1732 
1733  for (unsigned int ii = 0; ii < dim; ++ii)
1734  {
1735 
1736  component_at_qp = 0;
1737 
1738  // populate the vector of components at the qps
1739  // from vector_of_tensors_at_qp
1740  // vector_of_tensors_at_qp data is in form:
1741  // columns: 0, 1, ..., dim
1742  // rows: 0,1,...., n_quad_points
1743  // so extract the ii'th column of vector_of_tensors_at_qp
1744  for (unsigned int q = 0; q < n_quad_points; ++q)
1745  {
1746  component_at_qp(q) = vector_of_tensors_at_qp[q][ii];
1747  }
1748 
1749  // project from the qps -> nodes
1750  // component_at_node = projection_matrix_u * component_at_qp
1751  projection_matrix.vmult(component_at_node, component_at_qp);
1752 
1753  // rewrite the projection of the components
1754  // back into the vector of tensors
1755  for (unsigned int nn =0; nn <n_support_points; ++nn)
1756  {
1757  vector_of_tensors_at_nodes[nn][ii] = component_at_node(nn);
1758  }
1759  }
1760  }
1761 
1762 
1763 
1764  template <int dim>
1765  void
1767  const FullMatrix<double> &projection_matrix,
1768  const std::vector< SymmetricTensor<2, dim > > &vector_of_tensors_at_qp,
1769  std::vector< SymmetricTensor<2, dim > > &vector_of_tensors_at_nodes)
1770  {
1771 
1772  // check that the number columns of the projection_matrix
1773  // matches the size of the vector_of_tensors_at_qp
1774  Assert(projection_matrix.n_cols() == vector_of_tensors_at_qp.size(),
1775  ExcDimensionMismatch(projection_matrix.n_cols(),
1776  vector_of_tensors_at_qp.size()));
1777 
1778  // check that the number rows of the projection_matrix
1779  // matches the size of the vector_of_tensors_at_nodes
1780  Assert(projection_matrix.n_rows() == vector_of_tensors_at_nodes.size(),
1781  ExcDimensionMismatch(projection_matrix.n_rows(),
1782  vector_of_tensors_at_nodes.size()));
1783 
1784  // number of support points (nodes)
1785  const unsigned int n_support_points = projection_matrix.n_rows();
1786  // number of quadrature points to project from
1787  const unsigned int n_quad_points = projection_matrix.n_cols();
1788 
1789  // number of unique entries in a symmetric second-order tensor
1790  const unsigned int n_independent_components =
1792 
1793  // component projected to the nodes
1794  Vector<double> component_at_node(n_support_points);
1795  // component at the quadrature point
1796  Vector<double> component_at_qp(n_quad_points);
1797 
1798  // loop over the number of unique dimensions of the tensor
1799  for (unsigned int ii = 0; ii < n_independent_components; ++ii)
1800  {
1801 
1802  component_at_qp = 0;
1803 
1804  // row-column entry of tensor corresponding the unrolled index
1806  const unsigned int row = row_column_index[0];
1807  const unsigned int column = row_column_index[1];
1808 
1809  // populate the vector of components at the qps
1810  // from vector_of_tensors_at_qp
1811  // vector_of_tensors_at_qp is in form:
1812  // columns: 0, 1, ..., n_independent_components
1813  // rows: 0,1,...., n_quad_points
1814  // so extract the ii'th column of vector_of_tensors_at_qp
1815  for (unsigned int q = 0; q < n_quad_points; ++q)
1816  {
1817  component_at_qp(q) = (vector_of_tensors_at_qp[q])[row][column];
1818  }
1819 
1820  // project from the qps -> nodes
1821  // component_at_node = projection_matrix_u * component_at_qp
1822  projection_matrix.vmult(component_at_node, component_at_qp);
1823 
1824  // rewrite the projection of the components back into the vector of tensors
1825  for (unsigned int nn =0; nn <n_support_points; ++nn)
1826  {
1827  (vector_of_tensors_at_nodes[nn])[row][column] = component_at_node(nn);
1828  }
1829  }
1830  }
1831 
1832 
1833 
1834  template <int dim, int spacedim>
1835  void
1837  const Quadrature<dim-1> &lhs_quadrature,
1838  const Quadrature<dim-1> &rhs_quadrature,
1840  const unsigned int face,
1841  FullMatrix<double> &X)
1842  {
1843  Assert (fe.n_components() == 1, ExcNotImplemented());
1844  Assert (lhs_quadrature.size () > fe.degree, ExcNotGreaterThan (lhs_quadrature.size (), fe.degree));
1845 
1846 
1847 
1848  // build the matrices M and Q
1849  // described in the documentation
1851  FullMatrix<double> Q (fe.dofs_per_cell, rhs_quadrature.size());
1852 
1853  {
1854  // need an FEFaceValues object to evaluate shape function
1855  // values on the specified face.
1856  FEFaceValues <dim> fe_face_values (fe, lhs_quadrature, update_values);
1857  fe_face_values.reinit (cell, face); // setup shape_value on this face.
1858 
1859  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
1860  for (unsigned int j=0; j<fe.dofs_per_cell; ++j)
1861  for (unsigned int q=0; q<lhs_quadrature.size(); ++q)
1862  M(i,j) += fe_face_values.shape_value (i, q) *
1863  fe_face_values.shape_value (j, q) *
1864  lhs_quadrature.weight(q);
1865  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
1866  {
1867  M(i,i) = (M(i,i) == 0 ? 1 : M(i,i));
1868  }
1869  }
1870 
1871  {
1872  FEFaceValues <dim> fe_face_values (fe, rhs_quadrature, update_values);
1873  fe_face_values.reinit (cell, face); // setup shape_value on this face.
1874 
1875  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
1876  for (unsigned int q=0; q<rhs_quadrature.size(); ++q)
1877  Q(i,q) += fe_face_values.shape_value (i, q) *
1878  rhs_quadrature.weight(q);
1879  }
1880  // then invert M
1881  FullMatrix<double> M_inverse (fe.dofs_per_cell, fe.dofs_per_cell);
1882  M_inverse.invert (M);
1883 
1884  // finally compute the result
1885  X.reinit (fe.dofs_per_cell, rhs_quadrature.size());
1886  M_inverse.mmult (X, Q);
1887 
1888  Assert (X.m() == fe.dofs_per_cell, ExcInternalError());
1889  Assert (X.n() == rhs_quadrature.size(), ExcInternalError());
1890  }
1891 
1892 
1893 
1894  template <int dim>
1895  void
1896  hierarchic_to_lexicographic_numbering (unsigned int degree, std::vector<unsigned int> &h2l)
1897  {
1898  // number of support points in each
1899  // direction
1900  const unsigned int n = degree+1;
1901 
1902  unsigned int dofs_per_cell = n;
1903  for (unsigned int i=1; i<dim; ++i)
1904  dofs_per_cell *= n;
1905 
1906  // Assert size maches degree
1907  AssertDimension (h2l.size(), dofs_per_cell);
1908 
1909  // polynomial degree
1910  const unsigned int dofs_per_line = degree - 1;
1911 
1912  // the following lines of code are somewhat odd, due to the way the
1913  // hierarchic numbering is organized. if someone would really want to
1914  // understand these lines, you better draw some pictures where you
1915  // indicate the indices and orders of vertices, lines, etc, along with the
1916  // numbers of the degrees of freedom in hierarchical and lexicographical
1917  // order
1918  switch (dim)
1919  {
1920  case 1:
1921  {
1922  h2l[0] = 0;
1923  h2l[1] = dofs_per_cell-1;
1924  for (unsigned int i=2; i<dofs_per_cell; ++i)
1925  h2l[i] = i-1;
1926 
1927  break;
1928  }
1929 
1930  case 2:
1931  {
1932  unsigned int next_index = 0;
1933  // first the four vertices
1934  h2l[next_index++] = 0;
1935  h2l[next_index++] = n-1;
1936  h2l[next_index++] = n*(n-1);
1937  h2l[next_index++] = n*n-1;
1938 
1939  // left line
1940  for (unsigned int i=0; i<dofs_per_line; ++i)
1941  h2l[next_index++] = (1+i)*n;
1942 
1943  // right line
1944  for (unsigned int i=0; i<dofs_per_line; ++i)
1945  h2l[next_index++] = (2+i)*n-1;
1946 
1947  // bottom line
1948  for (unsigned int i=0; i<dofs_per_line; ++i)
1949  h2l[next_index++] = 1+i;
1950 
1951  // top line
1952  for (unsigned int i=0; i<dofs_per_line; ++i)
1953  h2l[next_index++] = n*(n-1)+i+1;
1954 
1955  // inside quad
1956  for (unsigned int i=0; i<dofs_per_line; ++i)
1957  for (unsigned int j=0; j<dofs_per_line; ++j)
1958  h2l[next_index++] = n*(i+1)+j+1;
1959 
1960  Assert (next_index == dofs_per_cell, ExcInternalError());
1961 
1962  break;
1963  }
1964 
1965  case 3:
1966  {
1967  unsigned int next_index = 0;
1968  // first the eight vertices
1969  h2l[next_index++] = 0; // 0
1970  h2l[next_index++] = ( 1)*degree; // 1
1971  h2l[next_index++] = ( n )*degree; // 2
1972  h2l[next_index++] = ( n+1)*degree; // 3
1973  h2l[next_index++] = (n*n )*degree; // 4
1974  h2l[next_index++] = (n*n +1)*degree; // 5
1975  h2l[next_index++] = (n*n+n )*degree; // 6
1976  h2l[next_index++] = (n*n+n+1)*degree; // 7
1977 
1978  // line 0
1979  for (unsigned int i=0; i<dofs_per_line; ++i)
1980  h2l[next_index++] = (i+1)*n;
1981  // line 1
1982  for (unsigned int i=0; i<dofs_per_line; ++i)
1983  h2l[next_index++] = n-1+(i+1)*n;
1984  // line 2
1985  for (unsigned int i=0; i<dofs_per_line; ++i)
1986  h2l[next_index++] = 1+i;
1987  // line 3
1988  for (unsigned int i=0; i<dofs_per_line; ++i)
1989  h2l[next_index++] = 1+i+n*(n-1);
1990 
1991  // line 4
1992  for (unsigned int i=0; i<dofs_per_line; ++i)
1993  h2l[next_index++] = (n-1)*n*n+(i+1)*n;
1994  // line 5
1995  for (unsigned int i=0; i<dofs_per_line; ++i)
1996  h2l[next_index++] = (n-1)*(n*n+1)+(i+1)*n;
1997  // line 6
1998  for (unsigned int i=0; i<dofs_per_line; ++i)
1999  h2l[next_index++] = n*n*(n-1)+i+1;
2000  // line 7
2001  for (unsigned int i=0; i<dofs_per_line; ++i)
2002  h2l[next_index++] = n*n*(n-1)+i+1+n*(n-1);
2003 
2004  // line 8
2005  for (unsigned int i=0; i<dofs_per_line; ++i)
2006  h2l[next_index++] = (i+1)*n*n;
2007  // line 9
2008  for (unsigned int i=0; i<dofs_per_line; ++i)
2009  h2l[next_index++] = n-1+(i+1)*n*n;
2010  // line 10
2011  for (unsigned int i=0; i<dofs_per_line; ++i)
2012  h2l[next_index++] = (i+1)*n*n+n*(n-1);
2013  // line 11
2014  for (unsigned int i=0; i<dofs_per_line; ++i)
2015  h2l[next_index++] = n-1+(i+1)*n*n+n*(n-1);
2016 
2017 
2018  // inside quads
2019  // face 0
2020  for (unsigned int i=0; i<dofs_per_line; ++i)
2021  for (unsigned int j=0; j<dofs_per_line; ++j)
2022  h2l[next_index++] = (i+1)*n*n+n*(j+1);
2023  // face 1
2024  for (unsigned int i=0; i<dofs_per_line; ++i)
2025  for (unsigned int j=0; j<dofs_per_line; ++j)
2026  h2l[next_index++] = (i+1)*n*n+n-1+n*(j+1);
2027  // face 2, note the orientation!
2028  for (unsigned int i=0; i<dofs_per_line; ++i)
2029  for (unsigned int j=0; j<dofs_per_line; ++j)
2030  h2l[next_index++] = (j+1)*n*n+i+1;
2031  // face 3, note the orientation!
2032  for (unsigned int i=0; i<dofs_per_line; ++i)
2033  for (unsigned int j=0; j<dofs_per_line; ++j)
2034  h2l[next_index++] = (j+1)*n*n+n*(n-1)+i+1;
2035  // face 4
2036  for (unsigned int i=0; i<dofs_per_line; ++i)
2037  for (unsigned int j=0; j<dofs_per_line; ++j)
2038  h2l[next_index++] = n*(i+1)+j+1;
2039  // face 5
2040  for (unsigned int i=0; i<dofs_per_line; ++i)
2041  for (unsigned int j=0; j<dofs_per_line; ++j)
2042  h2l[next_index++] = (n-1)*n*n+n*(i+1)+j+1;
2043 
2044  // inside hex
2045  for (unsigned int i=0; i<dofs_per_line; ++i)
2046  for (unsigned int j=0; j<dofs_per_line; ++j)
2047  for (unsigned int k=0; k<dofs_per_line; ++k)
2048  h2l[next_index++] = n*n*(i+1)+n*(j+1)+k+1;
2049 
2050  Assert (next_index == dofs_per_cell, ExcInternalError());
2051 
2052  break;
2053  }
2054 
2055  default:
2056  Assert (false, ExcNotImplemented());
2057  }
2058  }
2059 
2060 
2061 
2062  template <int dim>
2063  void
2065  std::vector<unsigned int> &h2l)
2066  {
2067  Assert (h2l.size() == fe.dofs_per_cell,
2068  ExcDimensionMismatch (h2l.size(), fe.dofs_per_cell));
2069  hierarchic_to_lexicographic_numbering<dim> (fe.dofs_per_line+1, h2l);
2070  }
2071 
2072 
2073 
2074  template <int dim>
2075  std::vector<unsigned int>
2077  {
2078  Assert (fe.n_components() == 1, ExcInvalidFE());
2079  std::vector<unsigned int> h2l(fe.dofs_per_cell);
2080  hierarchic_to_lexicographic_numbering<dim> (fe.dofs_per_line+1, h2l);
2081  return (h2l);
2082  }
2083 
2084  template <int dim>
2085  void
2087  std::vector<unsigned int> &l2h)
2088  {
2090  }
2091 
2092 
2093 
2094  template <int dim>
2095  std::vector<unsigned int>
2097  {
2099  }
2100 
2101 } // end of namespace FETools
2102 
2103 
2104 
2105 /*-------------- Explicit Instantiations -------------------------------*/
2106 #include "fe_tools.inst"
2107 
2108 
2109 /*---------------------------- fe_tools.cc ---------------------------*/
2110 
2111 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
Shape function values.
size_type m() const
static const unsigned int invalid_unsigned_int
Definition: types.h:164
void compute_interpolation_to_quadrature_points_matrix(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, FullMatrix< double > &I_q)
Definition: fe_tools.cc:1687
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
void compute_component_wise(const FiniteElement< dim, spacedim > &fe, std::vector< unsigned int > &renumbering, std::vector< std::vector< unsigned int > > &start_indices)
Definition: fe_tools.cc:343
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:175
void compute_block_renumbering(const FiniteElement< dim, spacedim > &fe, std::vector< types::global_dof_index > &renumbering, std::vector< types::global_dof_index > &block_data, bool return_start_indices=true)
Definition: fe_tools.cc:385
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe.cc:851
void gauss_jordan()
const std::vector< Point< dim > > & get_generalized_support_points() const
Definition: fe.cc:981
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
virtual void interpolate(std::vector< double > &local_dofs, const std::vector< double > &values) const
Definition: fe.cc:1103
Definition: fe_dgq.h:81
::ExceptionBase & ExcMessage(std::string arg1)
std::string trim(const std::string &input)
Definition: utilities.cc:131
void compute_face_embedding_matrices(const FiniteElement< dim, spacedim > &fe, FullMatrix< number >(&matrices)[GeometryInfo< dim >::max_children_per_face], const unsigned int face_coarse, const unsigned int face_fine, const double threshold=1.e-12)
Definition: fe_tools.cc:914
const unsigned int dofs_per_quad
Definition: fe_base.h:238
static TableIndices< rank > unrolled_to_component_indices(const unsigned int i)
std::pair< int, unsigned int > get_integer_at_position(const std::string &name, const unsigned int position)
Definition: utilities.cc:406
void hierarchic_to_lexicographic_numbering(unsigned int degree, std::vector< unsigned int > &h2l)
Definition: fe_tools.cc:1896
const unsigned int degree
Definition: fe_base.h:299
const Point< dim > & point(const unsigned int i) const
bool has_generalized_support_points() const
Definition: fe.cc:996
void invert(const FullMatrix< number2 > &M)
ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:221
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: tria.cc:10397
void compute_embedding_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false, const double threshold=1.e-12)
Definition: fe_tools.cc:891
Transformed quadrature points.
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
const unsigned int dofs_per_line
Definition: fe_base.h:232
cell_iterator begin(const unsigned int level=0) const
Definition: tria.cc:10377
Definition: fe_bdm.h:56
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
unsigned int n_blocks() const
virtual FiniteElement< dim, spacedim > * get(const unsigned int degree) const =0
void get_projection_matrix(const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &matrix)
Definition: fe_tools.cc:549
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
Definition: fe_abf.h:103
cell_iterator end() const
Definition: tria.cc:10465
const std::vector< Point< spacedim > > & get_quadrature_points() const
virtual void execute_coarsening_and_refinement()
Definition: tria.cc:11602
Definition: fe_dgp.h:311
Definition: fe_q.h:522
static unsigned int n_children(const RefinementCase< dim > &refinement_case)
const unsigned int first_quad_index
Definition: fe_base.h:254
size_type n() const
void get_interpolation_matrix(const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &interpolation_matrix)
Definition: fe_tools.cc:433
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
void compute_projection_from_face_quadrature_points_matrix(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim-1 > &lhs_quadrature, const Quadrature< dim-1 > &rhs_quadrature, const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face, FullMatrix< double > &X)
Definition: fe_tools.cc:1836
unsigned int global_dof_index
Definition: types.h:88
static unsigned int face_to_cell_lines(const unsigned int face, const unsigned int line, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
#define Assert(cond, exc)
Definition: exceptions.h:294
unsigned int element_multiplicity(const unsigned int index) const
Definition: fe.h:2766
UpdateFlags
static void project_to_subface(const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim-1 > &ref_case=RefinementCase< dim-1 >::isotropic_refinement)
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell, const unsigned int face_no)
Definition: fe_values.cc:3929
bool is_primitive(const unsigned int i) const
Definition: fe.h:2937
Definition: fe.h:34
void compute_projection_from_quadrature_points_matrix(const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &lhs_quadrature, const Quadrature< dim > &rhs_quadrature, FullMatrix< double > &X)
Definition: fe_tools.cc:1647
void get_back_interpolation_matrix(const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &interpolation_matrix)
Definition: fe_tools.cc:497
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
void compute_projection_from_quadrature_points(const FullMatrix< double > &projection_matrix, const std::vector< Tensor< 1, dim > > &vector_of_tensors_at_qp, std::vector< Tensor< 1, dim > > &vector_of_tensors_at_nodes)
Definition: fe_tools.cc:1705
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1236
virtual ~FEFactoryBase()
Definition: fe_tools.cc:338
std::vector< unsigned int > invert_permutation(const std::vector< unsigned int > &permutation)
Definition: utilities.cc:483
virtual FiniteElement< FE::dimension, FE::space_dimension > * get(const unsigned int degree) const
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:283
void lexicographic_to_hierarchic_numbering(const FiniteElementData< dim > &fe_data, std::vector< unsigned int > &l2h)
Definition: fe_tools.cc:2086
const unsigned int n_quadrature_points
Definition: fe_values.h:1457
const std::vector< Point< dim > > & get_unit_support_points() const
Definition: fe.cc:956
FiniteElement< dim, dim > * get_fe_from_name(const std::string &name)
Definition: fe_tools.cc:1638
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
Definition: fe.h:2742
void add_fe_name(const std::string &name, const FEFactoryBase< dim, spacedim > *factory)
Definition: fe_tools.cc:1277
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:164
unsigned int n_components() const
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > system_to_base_index(const unsigned int index) const
Definition: fe.h:2827
const std::vector< double > & get_JxW_values() const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
Definition: fe.h:2885
double JxW(const unsigned int quadrature_point) const
static void project_to_face(const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
bool conforms(const Conformity) const
void get_interpolation_difference_matrix(const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &difference_matrix)
Definition: fe_tools.cc:523
void add(const number a, const FullMatrix< number2 > &A)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
const unsigned int dofs_per_face
Definition: fe_base.h:276
const unsigned int first_line_index
Definition: fe_base.h:249
double least_squares(Vector< number2 > &dst, const Vector< number2 > &src) const
void refine_global(const unsigned int times=1)
Definition: tria.cc:9337
void compute_projection_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false)
Definition: fe_tools.cc:1117
const unsigned int dofs_per_vertex
Definition: fe_base.h:226
FiniteElement< dim, spacedim > * get_fe_by_name(const std::string &name)
Definition: fe_tools.cc:1549
unsigned int size(const unsigned int i) const
unsigned int n_base_elements() const
Definition: fe.h:2756
void compute_node_matrix(FullMatrix< double > &M, const FiniteElement< dim, spacedim > &fe)
Definition: fe_tools.cc:633
void reinit(const TriaIterator< DoFCellAccessor< DoFHandlerType< dim, spacedim >, level_dof_access > > &cell)
Definition: fe_values.cc:3739
double weight(const unsigned int i) const
Task< RT > new_task(const std_cxx11::function< RT()> &function)
unsigned int tensor_degree() const