Reference documentation for deal.II version 8.5.1
fe_face.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2009 - 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/fe/fe_face.h>
18 #include <deal.II/fe/fe_poly_face.templates.h>
19 #include <deal.II/fe/fe_nothing.h>
20 #include <deal.II/fe/fe_tools.h>
21 #include <deal.II/base/quadrature_lib.h>
22 #include <deal.II/lac/householder.h>
23 #include <sstream>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 
28 namespace
29 {
30  std::vector<Point<1> >
31  get_QGaussLobatto_points (const unsigned int degree)
32  {
33  if (degree > 0)
34  return QGaussLobatto<1>(degree+1).get_points();
35  else
36  return std::vector<Point<1> >(1, Point<1>(0.5));
37  }
38 }
39 
40 template <int dim, int spacedim>
41 FE_FaceQ<dim,spacedim>::FE_FaceQ (const unsigned int degree)
42  :
43  FE_PolyFace<TensorProductPolynomials<dim-1>, dim, spacedim> (
44  TensorProductPolynomials<dim-1>(Polynomials::generate_complete_Lagrange_basis(get_QGaussLobatto_points(degree))),
45  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree, FiniteElementData<dim>::L2),
46  std::vector<bool>(1,true))
47 {
48  // initialize unit face support points
49  const unsigned int codim = dim-1;
50  this->unit_face_support_points.resize(Utilities::fixed_power<codim>(this->degree+1));
51 
52  if (this->degree == 0)
53  for (unsigned int d=0; d<codim; ++d)
54  this->unit_face_support_points[0][d] = 0.5;
55  else
56  {
57  std::vector<Point<1> > points = get_QGaussLobatto_points(degree);
58 
59  unsigned int k=0;
60  for (unsigned int iz=0; iz <= ((codim>2) ? this->degree : 0) ; ++iz)
61  for (unsigned int iy=0; iy <= ((codim>1) ? this->degree : 0) ; ++iy)
62  for (unsigned int ix=0; ix<=this->degree; ++ix)
63  {
64  Point<codim> p;
65 
66  p(0) = points[ix][0];
67  if (codim>1)
68  p(1) = points[iy][0];
69  if (codim>2)
70  p(2) = points[iz][0];
71 
72  this->unit_face_support_points[k++] = p;
73  }
74  AssertDimension (k, this->unit_face_support_points.size());
75  }
76 
77  // initialize unit support points (this makes it possible to assign initial
78  // values to FE_FaceQ)
80  this->unit_face_support_points.size());
81  const unsigned int n_face_dofs = this->unit_face_support_points.size();
82  for (unsigned int i=0; i<n_face_dofs; ++i)
83  for (unsigned int d=0; d<dim; ++d)
84  {
85  for (unsigned int e=0, c=0; e<dim; ++e)
86  if (d!=e)
87  {
88  // faces in y-direction are oriented differently
89  unsigned int renumber = i;
90  if (dim == 3 && d == 1)
91  renumber = i/(degree+1)+(degree+1)*(i%(degree+1));
92  this->unit_support_points[n_face_dofs*2*d+i][e] =
93  this->unit_face_support_points[renumber][c];
94  this->unit_support_points[n_face_dofs*(2*d+1)+i][e] =
95  this->unit_face_support_points[renumber][c];
96  this->unit_support_points[n_face_dofs*(2*d+1)+i][d] = 1;
97  ++c;
98  }
99  }
100 }
101 
102 
103 
104 template <int dim, int spacedim>
107 {
108  return new FE_FaceQ<dim,spacedim>(this->degree);
109 }
110 
111 
112 
113 template <int dim, int spacedim>
114 std::string
116 {
117  // note that the FETools::get_fe_by_name function depends on the
118  // particular format of the string this function returns, so they have to be
119  // kept in synch
120  std::ostringstream namebuf;
121  namebuf << "FE_FaceQ<"
122  << Utilities::dim_string(dim,spacedim)
123  << ">(" << this->degree << ")";
124 
125  return namebuf.str();
126 }
127 
128 
129 
130 template <int dim, int spacedim>
131 void
134  FullMatrix<double> &interpolation_matrix) const
135 {
136  get_subface_interpolation_matrix (source_fe, numbers::invalid_unsigned_int,
137  interpolation_matrix);
138 }
139 
140 
141 
142 template <int dim, int spacedim>
143 void
146  const unsigned int subface,
147  FullMatrix<double> &interpolation_matrix) const
148 {
149  // this function is similar to the respective method in FE_Q
150 
151  Assert (interpolation_matrix.n() == this->dofs_per_face,
152  ExcDimensionMismatch (interpolation_matrix.n(),
153  this->dofs_per_face));
154  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
155  ExcDimensionMismatch (interpolation_matrix.m(),
156  x_source_fe.dofs_per_face));
157 
158  // see if source is a FaceQ element
159  if (const FE_FaceQ<dim,spacedim> *source_fe
160  = dynamic_cast<const FE_FaceQ<dim,spacedim> *>(&x_source_fe))
161  {
162 
163  // Make sure that the element for which the DoFs should be constrained
164  // is the one with the higher polynomial degree. Actually the procedure
165  // will work also if this assertion is not satisfied. But the matrices
166  // produced in that case might lead to problems in the hp procedures,
167  // which use this method.
168  Assert (this->dofs_per_face <= source_fe->dofs_per_face,
169  (typename FiniteElement<dim,spacedim>::
170  ExcInterpolationNotImplemented ()));
171 
172  // generate a quadrature with the unit face support points.
173  const Quadrature<dim-1> face_quadrature (source_fe->get_unit_face_support_points ());
174 
175  // Rule of thumb for FP accuracy, that can be expected for a given
176  // polynomial degree. This value is used to cut off values close to
177  // zero.
178  const double eps = 2e-13*(this->degree+1)*(dim-1);
179 
180  // compute the interpolation matrix by simply taking the value at the
181  // support points.
182  for (unsigned int i=0; i<source_fe->dofs_per_face; ++i)
183  {
184  const Point<dim-1> p =
186  ?
187  face_quadrature.point(i)
188  :
189  GeometryInfo<dim-1>::child_to_cell_coordinates (face_quadrature.point(i),
190  subface);
191 
192  for (unsigned int j=0; j<this->dofs_per_face; ++j)
193  {
194  double matrix_entry = this->poly_space.compute_value (j, p);
195 
196  // Correct the interpolated value. I.e. if it is close to 1 or 0,
197  // make it exactly 1 or 0. Unfortunately, this is required to avoid
198  // problems with higher order elements.
199  if (std::fabs (matrix_entry - 1.0) < eps)
200  matrix_entry = 1.0;
201  if (std::fabs (matrix_entry) < eps)
202  matrix_entry = 0.0;
203 
204  interpolation_matrix(i,j) = matrix_entry;
205  }
206  }
207 
208  // make sure that the row sum of each of the matrices is 1 at this
209  // point. this must be so since the shape functions sum up to 1
210  for (unsigned int j=0; j<source_fe->dofs_per_face; ++j)
211  {
212  double sum = 0.;
213 
214  for (unsigned int i=0; i<this->dofs_per_face; ++i)
215  sum += interpolation_matrix(j,i);
216 
217  Assert (std::fabs(sum-1) < eps, ExcInternalError());
218  }
219  }
220  else if (dynamic_cast<const FE_Nothing<dim> *>(&x_source_fe) != 0)
221  {
222  // nothing to do here, the FE_Nothing has no degrees of freedom anyway
223  }
224  else
225  AssertThrow (false,(typename FiniteElement<dim,spacedim>::
226  ExcInterpolationNotImplemented()));
227 }
228 
229 
230 
231 template <int dim, int spacedim>
232 bool
234  const unsigned int shape_index,
235  const unsigned int face_index) const
236 {
237  return (face_index == (shape_index/this->dofs_per_face));
238 }
239 
240 
241 
242 template <int dim, int spacedim>
243 std::vector<unsigned int>
245 {
246  std::vector<unsigned int> dpo(dim+1, 0U);
247  dpo[dim-1] = deg+1;
248  for (unsigned int i=1; i<dim-1; ++i)
249  dpo[dim-1] *= deg+1;
250  return dpo;
251 }
252 
253 
254 
255 template <int dim, int spacedim>
256 bool
258 {
259  return true;
260 }
261 
262 
263 
264 template <int dim, int spacedim>
268 {
269  if (const FE_FaceQ<dim,spacedim> *fe_q_other
270  = dynamic_cast<const FE_FaceQ<dim,spacedim>*>(&fe_other))
271  {
272  if (this->degree < fe_q_other->degree)
274  else if (this->degree == fe_q_other->degree)
276  else
278  }
279  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
280  {
281  if (fe_nothing->is_dominating())
282  {
284  }
285  else
286  {
287  // the FE_Nothing has no degrees of freedom and it is typically used in
288  // a context where we don't require any continuity along the interface
290  }
291  }
292 
293  Assert (false, ExcNotImplemented());
295 }
296 
297 
298 
299 template <int dim, int spacedim>
300 std::pair<Table<2,bool>, std::vector<unsigned int> >
302 {
303  Table<2,bool> constant_modes(1, this->dofs_per_cell);
304  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
305  constant_modes(0,i) = true;
306  return std::pair<Table<2,bool>, std::vector<unsigned int> >
307  (constant_modes, std::vector<unsigned int>(1, 0));
308 }
309 
310 
311 
312 // ----------------------------- FE_FaceQ<1,spacedim> ------------------------
313 
314 template <int spacedim>
315 FE_FaceQ<1,spacedim>::FE_FaceQ (const unsigned int degree)
316  :
317  FiniteElement<1,spacedim> (FiniteElementData<1>(get_dpo_vector(degree), 1, degree, FiniteElementData<1>::L2),
318  std::vector<bool>(1,true),
319  std::vector<ComponentMask> (1, ComponentMask(1,true)))
320 {
321  this->unit_face_support_points.resize(1);
322 
323  // initialize unit support points (this makes it possible to assign initial
324  // values to FE_FaceQ)
326  this->unit_support_points[1] = Point<1>(1.);
327 }
328 
329 
330 
331 template <int spacedim>
334 {
335  return new FE_FaceQ<1,spacedim>(this->degree);
336 }
337 
338 
339 
340 template <int spacedim>
341 std::string
343 {
344  // note that the FETools::get_fe_by_name function depends on the
345  // particular format of the string this function returns, so they have to be
346  // kept in synch
347  std::ostringstream namebuf;
348  namebuf << "FE_FaceQ<"
349  << Utilities::dim_string(1,spacedim)
350  << ">(" << this->degree << ")";
351 
352  return namebuf.str();
353 }
354 
355 
356 
357 template <int spacedim>
358 void
361  FullMatrix<double> &interpolation_matrix) const
362 {
364  interpolation_matrix);
365 }
366 
367 
368 
369 template <int spacedim>
370 void
373  const unsigned int /*subface*/,
374  FullMatrix<double> &interpolation_matrix) const
375 {
376  (void)x_source_fe;
377  Assert (interpolation_matrix.n() == this->dofs_per_face,
378  ExcDimensionMismatch (interpolation_matrix.n(),
379  this->dofs_per_face));
380  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
381  ExcDimensionMismatch (interpolation_matrix.m(),
382  x_source_fe.dofs_per_face));
383  interpolation_matrix(0,0) = 1.;
384 }
385 
386 
387 
388 template <int spacedim>
389 bool
391  const unsigned int shape_index,
392  const unsigned int face_index) const
393 {
394  AssertIndexRange(shape_index, 2);
395  return (face_index == shape_index);
396 }
397 
398 
399 
400 template <int spacedim>
401 std::vector<unsigned int>
403 {
404  std::vector<unsigned int> dpo(2, 0U);
405  dpo[0] = 1;
406  return dpo;
407 }
408 
409 
410 
411 template <int spacedim>
412 bool
414 {
415  return true;
416 }
417 
418 
419 
420 template <int spacedim>
424 {
426 }
427 
428 
429 
430 template <int spacedim>
431 std::pair<Table<2,bool>, std::vector<unsigned int> >
433 {
434  Table<2,bool> constant_modes(1, this->dofs_per_cell);
435  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
436  constant_modes(0,i) = true;
437  return std::pair<Table<2,bool>, std::vector<unsigned int> >
438  (constant_modes, std::vector<unsigned int>(1,0));
439 }
440 
441 
442 
443 template <int spacedim>
446 {
447  UpdateFlags out = flags & update_values;
448  if (flags & update_gradients)
450  if (flags & update_hessians)
452  if (flags & update_normal_vectors)
454 
455  return out;
456 }
457 
458 
459 template <int spacedim>
460 void
464  const Quadrature<1> &,
465  const Mapping<1,spacedim> &,
466  const typename Mapping<1,spacedim>::InternalDataBase &,
467  const ::internal::FEValues::MappingRelatedData<1, spacedim> &,
470 {
471  // Do nothing, since we do not have values in the interior
472 }
473 
474 
475 
476 template <int spacedim>
477 void
480  const unsigned int face,
481  const Quadrature<0> &,
482  const Mapping<1,spacedim> &,
483  const typename Mapping<1,spacedim>::InternalDataBase &,
484  const ::internal::FEValues::MappingRelatedData<1, spacedim> &,
485  const typename FiniteElement<1,spacedim>::InternalDataBase &fe_internal,
487 {
488  const unsigned int foffset = face;
489  if (fe_internal.update_each & update_values)
490  {
491  for (unsigned int k=0; k<this->dofs_per_cell; ++k)
492  output_data.shape_values(k,0) = 0.;
493  output_data.shape_values(foffset,0) = 1;
494  }
495 }
496 
497 
498 template <int spacedim>
499 void
502  const unsigned int ,
503  const unsigned int ,
504  const Quadrature<0> &,
505  const Mapping<1,spacedim> &,
506  const typename Mapping<1,spacedim>::InternalDataBase &,
507  const ::internal::FEValues::MappingRelatedData<1, spacedim> &,
510 {
511  Assert(false, ExcMessage("There are no sub-face values to fill in 1D!"));
512 }
513 
514 
515 
516 // --------------------------------------- FE_FaceP --------------------------
517 
518 template <int dim, int spacedim>
519 FE_FaceP<dim,spacedim>::FE_FaceP (const unsigned int degree)
520  :
521  FE_PolyFace<PolynomialSpace<dim-1>, dim, spacedim>
522  (PolynomialSpace<dim-1>(Polynomials::Legendre::generate_complete_basis(degree)),
523  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree, FiniteElementData<dim>::L2),
524  std::vector<bool>(1,true))
525 {}
526 
527 
528 
529 template <int dim, int spacedim>
532 {
533  return new FE_FaceP<dim,spacedim>(this->degree);
534 }
535 
536 
537 
538 template <int dim, int spacedim>
539 std::string
541 {
542  // note that the FETools::get_fe_by_name function depends on the
543  // particular format of the string this function returns, so they have to be
544  // kept in synch
545  std::ostringstream namebuf;
546  namebuf << "FE_FaceP<"
547  << Utilities::dim_string(dim,spacedim)
548  << ">(" << this->degree << ")";
549 
550  return namebuf.str();
551 }
552 
553 
554 
555 template <int dim, int spacedim>
556 bool
558  const unsigned int shape_index,
559  const unsigned int face_index) const
560 {
561  return (face_index == (shape_index/this->dofs_per_face));
562 }
563 
564 
565 
566 template <int dim, int spacedim>
567 std::vector<unsigned int>
569 {
570  std::vector<unsigned int> dpo(dim+1, 0U);
571  dpo[dim-1] = deg+1;
572  for (unsigned int i=1; i<dim-1; ++i)
573  {
574  dpo[dim-1] *= deg+1+i;
575  dpo[dim-1] /= i+1;
576  }
577  return dpo;
578 }
579 
580 
581 
582 
583 template <int dim, int spacedim>
584 bool
586 {
587  return true;
588 }
589 
590 
591 
592 template <int dim, int spacedim>
596 {
597  if (const FE_FaceP<dim,spacedim> *fe_q_other
598  = dynamic_cast<const FE_FaceP<dim,spacedim>*>(&fe_other))
599  {
600  if (this->degree < fe_q_other->degree)
602  else if (this->degree == fe_q_other->degree)
604  else
606  }
607  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
608  {
609  if (fe_nothing->is_dominating())
610  {
612  }
613  else
614  {
615  // the FE_Nothing has no degrees of freedom and it is typically used in
616  // a context where we don't require any continuity along the interface
618  }
619  }
620 
621  Assert (false, ExcNotImplemented());
623 }
624 
625 
626 
627 
628 template <int dim, int spacedim>
629 void
632  FullMatrix<double> &interpolation_matrix) const
633 {
634  get_subface_interpolation_matrix (source_fe, numbers::invalid_unsigned_int,
635  interpolation_matrix);
636 }
637 
638 
639 
640 template <int dim, int spacedim>
641 void
644  const unsigned int subface,
645  FullMatrix<double> &interpolation_matrix) const
646 {
647  // this function is similar to the respective method in FE_Q
648 
649  Assert (interpolation_matrix.n() == this->dofs_per_face,
650  ExcDimensionMismatch (interpolation_matrix.n(),
651  this->dofs_per_face));
652  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
653  ExcDimensionMismatch (interpolation_matrix.m(),
654  x_source_fe.dofs_per_face));
655 
656  // see if source is a FaceP element
657  if (const FE_FaceP<dim,spacedim> *source_fe
658  = dynamic_cast<const FE_FaceP<dim,spacedim> *>(&x_source_fe))
659  {
660  // Make sure that the element for which the DoFs should be constrained
661  // is the one with the higher polynomial degree. Actually the procedure
662  // will work also if this assertion is not satisfied. But the matrices
663  // produced in that case might lead to problems in the hp procedures,
664  // which use this method.
665  Assert (this->dofs_per_face <= source_fe->dofs_per_face,
666  (typename FiniteElement<dim,spacedim>::
667  ExcInterpolationNotImplemented ()));
668 
669  // do this as in FETools by solving a least squares problem where we
670  // force the source FE polynomial to be equal the given FE on all
671  // quadrature points
672  const QGauss<dim-1> face_quadrature (source_fe->degree+1);
673 
674  // Rule of thumb for FP accuracy, that can be expected for a given
675  // polynomial degree. This value is used to cut off values close to
676  // zero.
677  const double eps = 2e-13*(this->degree+1)*(dim-1);
678 
679  FullMatrix<double> mass (face_quadrature.size(), source_fe->dofs_per_face);
680 
681  for (unsigned int k = 0; k < face_quadrature.size(); ++k)
682  {
683  const Point<dim-1> p =
684  subface == numbers::invalid_unsigned_int ?
685  face_quadrature.point(k) :
686  GeometryInfo<dim-1>::child_to_cell_coordinates (face_quadrature.point(k),
687  subface);
688 
689  for (unsigned int j = 0; j < source_fe->dofs_per_face; ++j)
690  mass (k , j) = source_fe->poly_space.compute_value(j, p);
691  }
692 
693  Householder<double> H(mass);
694  Vector<double> v_in(face_quadrature.size());
695  Vector<double> v_out(source_fe->dofs_per_face);
696 
697 
698  // compute the interpolation matrix by evaluating on the fine side and
699  // then solving the least squares problem
700  for (unsigned int i=0; i<this->dofs_per_face; ++i)
701  {
702  for (unsigned int k = 0; k < face_quadrature.size(); ++k)
703  {
704  const Point<dim-1> p = numbers::invalid_unsigned_int ?
705  face_quadrature.point(k) :
706  GeometryInfo<dim-1>::child_to_cell_coordinates (face_quadrature.point(k),
707  subface);
708  v_in(k) = this->poly_space.compute_value(i, p);
709  }
710  const double result = H.least_squares(v_out, v_in);
711  (void)result;
712  Assert(result < 1e-12, FETools::ExcLeastSquaresError (result));
713 
714  for (unsigned int j = 0; j < source_fe->dofs_per_face; ++j)
715  {
716  double matrix_entry = v_out(j);
717 
718  // Correct the interpolated value. I.e. if it is close to 1 or 0,
719  // make it exactly 1 or 0. Unfortunately, this is required to avoid
720  // problems with higher order elements.
721  if (std::fabs (matrix_entry - 1.0) < eps)
722  matrix_entry = 1.0;
723  if (std::fabs (matrix_entry) < eps)
724  matrix_entry = 0.0;
725 
726  interpolation_matrix(j,i) = matrix_entry;
727  }
728  }
729  }
730  else if (dynamic_cast<const FE_Nothing<dim> *>(&x_source_fe) != 0)
731  {
732  // nothing to do here, the FE_Nothing has no degrees of freedom anyway
733  }
734  else
735  AssertThrow (false,(typename FiniteElement<dim,spacedim>::
736  ExcInterpolationNotImplemented()));
737 }
738 
739 
740 
741 template <int dim, int spacedim>
742 std::pair<Table<2,bool>, std::vector<unsigned int> >
744 {
745  Table<2,bool> constant_modes(1, this->dofs_per_cell);
746  for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
747  constant_modes(0, face*this->dofs_per_face) = true;
748  return std::pair<Table<2,bool>, std::vector<unsigned int> >
749  (constant_modes, std::vector<unsigned int>(1, 0));
750 }
751 
752 
753 
754 template <int spacedim>
755 FE_FaceP<1,spacedim>::FE_FaceP (const unsigned int degree)
756  :
757  FE_FaceQ<1,spacedim> (degree)
758 {}
759 
760 
761 
762 template <int spacedim>
763 std::string
765 {
766  // note that the FETools::get_fe_by_name function depends on the
767  // particular format of the string this function returns, so they have to be
768  // kept in synch
769  std::ostringstream namebuf;
770  namebuf << "FE_FaceP<"
771  << Utilities::dim_string(1,spacedim)
772  << ">(" << this->degree << ")";
773 
774  return namebuf.str();
775 }
776 
777 
778 
779 // explicit instantiations
780 #include "fe_face.inst"
781 
782 
783 DEAL_II_NAMESPACE_CLOSE
static Point< dim > child_to_cell_coordinates(const Point< dim > &p, const unsigned int child_index, const RefinementCase< dim > refine_case=RefinementCase< dim >::isotropic_refinement)
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_face.cc:557
Transformed quadrature weights.
Shape function values.
size_type m() const
static const unsigned int invalid_unsigned_int
Definition: types.h:170
static ::ExceptionBase & ExcLeastSquaresError(double arg1)
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
virtual FiniteElement< dim, spacedim > * clone() const
Definition: fe_face.cc:531
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
const unsigned int degree
Definition: fe_base.h:311
virtual FiniteElement< dim, spacedim > * clone() const
Definition: fe_face.cc:106
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_face.cc:595
STL namespace.
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_face.cc:133
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_face.cc:145
Definition: point.h:89
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2223
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_face.cc:301
virtual std::string get_name() const
Definition: fe_face.cc:540
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:46
std::vector< Point< dim-1 > > unit_face_support_points
Definition: fe.h:2230
virtual bool hp_constraints_are_implemented() const
Definition: fe_face.cc:585
FE_FaceQ(const unsigned int p)
Definition: fe_face.cc:41
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
Definition: fe_face.cc:568
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_face.cc:233
Second derivatives of shape functions.
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:161
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_face.cc:643
SymmetricTensor< rank, dim, Number > sum(const SymmetricTensor< rank, dim, Number > &local, const MPI_Comm &mpi_communicator)
const unsigned int dofs_per_cell
Definition: fe_base.h:295
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_face.cc:631
Shape function gradients.
Normal vectors.
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
Definition: fe_face.cc:244
const unsigned int dofs_per_face
Definition: fe_base.h:288
double least_squares(Vector< number2 > &dst, const Vector< number2 > &src) const
static ::ExceptionBase & ExcNotImplemented()
FE_FaceP(unsigned int p)
Definition: fe_face.cc:519
virtual std::string get_name() const
Definition: fe_face.cc:115
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_face.cc:743
Covariant transformation.
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_face.cc:267
static ::ExceptionBase & ExcInternalError()
virtual bool hp_constraints_are_implemented() const
Definition: fe_face.cc:257