Reference documentation for deal.II version 8.4.2
fe_trace.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2015 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/config.h>
17 #include <deal.II/base/tensor_product_polynomials.h>
18 #include <deal.II/base/quadrature_lib.h>
19 #include <deal.II/base/qprojector.h>
20 #include <deal.II/fe/fe_poly_face.h>
21 #include <deal.II/fe/fe_q.h>
22 #include <deal.II/fe/fe_nothing.h>
23 #include <deal.II/fe/fe_tools.h>
24 #include <deal.II/fe/fe_trace.h>
25 
26 #include <sstream>
27 #include <fstream>
28 #include <iostream>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
33 
34 template <int dim, int spacedim>
35 FE_TraceQ<dim,spacedim>::FE_TraceQ (const unsigned int degree)
36  :
37  FE_PolyFace<TensorProductPolynomials<dim-1>, dim, spacedim> (
38  TensorProductPolynomials<dim-1>(Polynomials::LagrangeEquidistant::generate_complete_basis(degree)),
39  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree, FiniteElementData<dim>::L2),
40  std::vector<bool>(1,true)),
41  fe_q (degree)
42 {
43  Assert (degree > 0,
44  ExcMessage ("FE_Trace can only be used for polynomial degrees "
45  "greater than zero"));
46  std::vector<unsigned int> renumber (this->dofs_per_face);
48  this->poly_space.set_numbering(renumber);
49 
50  // Initialize face support points
51  this->unit_face_support_points = fe_q.get_unit_face_support_points();
52  // Initialize constraint matrices
53  this->interface_constraints = fe_q.constraints();
54 }
55 
56 
57 
58 template <int dim, int spacedim>
61 {
62  return new FE_TraceQ<dim,spacedim>(this->degree);
63 }
64 
65 
66 
67 template <int dim, int spacedim>
68 std::string
70 {
71  // note that the FETools::get_fe_from_name function depends on the
72  // particular format of the string this function returns, so they have to be
73  // kept in synch
74 
75  std::ostringstream namebuf;
76  namebuf << "FE_TraceQ<"
77  << Utilities::dim_string(dim,spacedim)
78  << ">(" << this->degree << ")";
79 
80  return namebuf.str();
81 }
82 
83 
84 
85 template <int dim, int spacedim>
86 bool
87 FE_TraceQ<dim,spacedim>::has_support_on_face (const unsigned int shape_index,
88  const unsigned int face_index) const
89 {
90  Assert (shape_index < this->dofs_per_cell,
91  ExcIndexRange (shape_index, 0, this->dofs_per_cell));
93  ExcIndexRange (face_index, 0, GeometryInfo<dim>::faces_per_cell));
94 
95  // FE_TraceQ shares the numbering of elemental degrees of freedom with FE_Q
96  // except for the missing interior ones (quad dofs in 2D and hex dofs in
97  // 3D). Therefore, it is safe to ask fe_q for the corresponding
98  // information. The assertion 'shape_index < this->dofs_per_cell' will make
99  // sure that we only access the trace dofs.
100  return fe_q.has_support_on_face (shape_index, face_index);
101 }
102 
103 
104 
105 template <int dim, int spacedim>
106 std::pair<Table<2,bool>, std::vector<unsigned int> >
108 {
109  Table<2,bool> constant_modes(1, this->dofs_per_cell);
110  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
111  constant_modes(0,i) = true;
112  return std::pair<Table<2,bool>, std::vector<unsigned int> >
113  (constant_modes, std::vector<unsigned int>(1, 0));
114 }
115 
116 
117 
118 template <int dim, int spacedim>
119 std::vector<unsigned int>
121 {
122  // This constructs FE_TraceQ in exactly the same way as FE_Q except for the
123  // interior degrees of freedom that are not present here (line in 1D, quad
124  // in 2D, hex in 3D).
125  AssertThrow(deg>0,ExcMessage("FE_TraceQ needs to be of degree > 0."));
126  std::vector<unsigned int> dpo(dim+1, 1U);
127  dpo[dim]=0;
128  dpo[0]=1;
129  for (unsigned int i=1; i<dim; ++i)
130  dpo[i] = dpo[i-1]*(deg-1);
131  return dpo;
132 }
133 
134 
135 
136 template <int dim, int spacedim>
137 bool
139 {
140  return fe_q.hp_constraints_are_implemented ();
141 }
142 
143 
144 template <int dim, int spacedim>
148 {
149  if (const FE_TraceQ<dim,spacedim> *fe_q_other
150  = dynamic_cast<const FE_TraceQ<dim,spacedim>*>(&fe_other))
151  {
152  if (this->degree < fe_q_other->degree)
153  return FiniteElementDomination::this_element_dominates;
154  else if (this->degree == fe_q_other->degree)
155  return FiniteElementDomination::either_element_can_dominate;
156  else
157  return FiniteElementDomination::other_element_dominates;
158  }
159  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
160  {
161  if (fe_nothing->is_dominating())
162  {
163  return FiniteElementDomination::other_element_dominates;
164  }
165  else
166  {
167  // the FE_Nothing has no degrees of freedom and it is typically used in
168  // a context where we don't require any continuity along the interface
169  return FiniteElementDomination::no_requirements;
170  }
171  }
172 
173  Assert (false, ExcNotImplemented());
174  return FiniteElementDomination::neither_element_dominates;
175 }
176 
177 
178 
179 template <int dim, int spacedim>
180 void
183  FullMatrix<double> &interpolation_matrix) const
184 {
186  interpolation_matrix);
187 }
188 
189 
190 
191 template <int dim, int spacedim>
192 void
195  const unsigned int subface,
196  FullMatrix<double> &interpolation_matrix) const
197 {
198  // this is the code from FE_FaceQ
199  Assert (interpolation_matrix.n() == this->dofs_per_face,
200  ExcDimensionMismatch (interpolation_matrix.n(),
201  this->dofs_per_face));
202  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
203  ExcDimensionMismatch (interpolation_matrix.m(),
204  x_source_fe.dofs_per_face));
205 
206  // see if source is a FaceQ element
207  if (const FE_TraceQ<dim,spacedim> *source_fe
208  = dynamic_cast<const FE_TraceQ<dim,spacedim> *>(&x_source_fe))
209  {
210  fe_q.get_subface_interpolation_matrix (source_fe->fe_q, subface, interpolation_matrix);
211  }
212  else if (dynamic_cast<const FE_Nothing<dim> *>(&x_source_fe) != 0)
213  {
214  // nothing to do here, the FE_Nothing has no degrees of freedom anyway
215  }
216  else
217  AssertThrow (false,(typename FiniteElement<dim,spacedim>::
218  ExcInterpolationNotImplemented()));
219 }
220 
221 
222 
223 template <int spacedim>
225  :
226  FE_FaceQ<1,spacedim> (degree)
227 {}
228 
229 
230 
231 template <int spacedim>
232 std::string
234 {
235  // note that the FETools::get_fe_from_name function depends on the
236  // particular format of the string this function returns, so they have to be
237  // kept in synch
238  std::ostringstream namebuf;
239  namebuf << "FE_TraceQ<"
240  << Utilities::dim_string(1,spacedim)
241  << ">(" << this->degree << ")";
242 
243  return namebuf.str();
244 }
245 
246 
247 
248 // explicit instantiations
249 #include "fe_trace.inst"
250 
251 
252 DEAL_II_NAMESPACE_CLOSE
size_type m() const
static const unsigned int invalid_unsigned_int
Definition: types.h:164
FullMatrix< double > interface_constraints
Definition: fe.h:2132
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_trace.cc:194
FE_TraceQ(unsigned int p)
Definition: fe_trace.cc:35
::ExceptionBase & ExcMessage(std::string arg1)
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
void set_numbering(const std::vector< unsigned int > &renumber)
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_trace.cc:182
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:294
std::vector< Point< dim-1 > > unit_face_support_points
Definition: fe.h:2151
FE_Q< dim, spacedim > fe_q
Definition: fe_trace.h:132
virtual bool hp_constraints_are_implemented() const
Definition: fe_trace.cc:138
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_trace.cc:87
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:158
const unsigned int dofs_per_cell
Definition: fe_base.h:283
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_trace.cc:147
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_trace.cc:107
virtual FiniteElement< dim, spacedim > * clone() const
Definition: fe_trace.cc:60
const unsigned int dofs_per_face
Definition: fe_base.h:276
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
Definition: fe_trace.cc:120
virtual std::string get_name() const
Definition: fe_trace.cc:69