Reference documentation for deal.II version 8.5.1
fe_trace.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 #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::generate_complete_Lagrange_basis(QGaussLobatto<1>(degree+1).get_points())),
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 
53  // initialize unit support points (this makes it possible to assign initial
54  // values to FE_TraceQ). Note that we simply take the points of fe_q but
55  // skip the last ones which are associated with the interior of FE_Q.
56  this->unit_support_points.resize(this->dofs_per_cell);
57  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
58  this->unit_support_points[i] = fe_q.get_unit_support_points()[i];
59 
60  // Initialize constraint matrices
61  this->interface_constraints = fe_q.constraints();
62 }
63 
64 
65 
66 template <int dim, int spacedim>
69 {
70  return new FE_TraceQ<dim,spacedim>(this->degree);
71 }
72 
73 
74 
75 template <int dim, int spacedim>
76 std::string
78 {
79  // note that the FETools::get_fe_by_name function depends on the
80  // particular format of the string this function returns, so they have to be
81  // kept in synch
82 
83  std::ostringstream namebuf;
84  namebuf << "FE_TraceQ<"
85  << Utilities::dim_string(dim,spacedim)
86  << ">(" << this->degree << ")";
87 
88  return namebuf.str();
89 }
90 
91 
92 
93 template <int dim, int spacedim>
94 bool
95 FE_TraceQ<dim,spacedim>::has_support_on_face (const unsigned int shape_index,
96  const unsigned int face_index) const
97 {
98  Assert (shape_index < this->dofs_per_cell,
99  ExcIndexRange (shape_index, 0, this->dofs_per_cell));
102 
103  // FE_TraceQ shares the numbering of elemental degrees of freedom with FE_Q
104  // except for the missing interior ones (quad dofs in 2D and hex dofs in
105  // 3D). Therefore, it is safe to ask fe_q for the corresponding
106  // information. The assertion 'shape_index < this->dofs_per_cell' will make
107  // sure that we only access the trace dofs.
108  return fe_q.has_support_on_face (shape_index, face_index);
109 }
110 
111 
112 
113 template <int dim, int spacedim>
114 std::pair<Table<2,bool>, std::vector<unsigned int> >
116 {
117  Table<2,bool> constant_modes(1, this->dofs_per_cell);
118  for (unsigned int i=0; i<this->dofs_per_cell; ++i)
119  constant_modes(0,i) = true;
120  return std::pair<Table<2,bool>, std::vector<unsigned int> >
121  (constant_modes, std::vector<unsigned int>(1, 0));
122 }
123 
124 
125 
126 template <int dim, int spacedim>
127 std::vector<unsigned int>
129 {
130  // This constructs FE_TraceQ in exactly the same way as FE_Q except for the
131  // interior degrees of freedom that are not present here (line in 1D, quad
132  // in 2D, hex in 3D).
133  AssertThrow(deg>0,ExcMessage("FE_TraceQ needs to be of degree > 0."));
134  std::vector<unsigned int> dpo(dim+1, 1U);
135  dpo[dim]=0;
136  dpo[0]=1;
137  for (unsigned int i=1; i<dim; ++i)
138  dpo[i] = dpo[i-1]*(deg-1);
139  return dpo;
140 }
141 
142 
143 
144 template <int dim, int spacedim>
145 bool
147 {
148  return fe_q.hp_constraints_are_implemented ();
149 }
150 
151 
152 template <int dim, int spacedim>
156 {
157  if (const FE_TraceQ<dim,spacedim> *fe_q_other
158  = dynamic_cast<const FE_TraceQ<dim,spacedim>*>(&fe_other))
159  {
160  if (this->degree < fe_q_other->degree)
162  else if (this->degree == fe_q_other->degree)
164  else
166  }
167  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
168  {
169  if (fe_nothing->is_dominating())
170  {
172  }
173  else
174  {
175  // the FE_Nothing has no degrees of freedom and it is typically used in
176  // a context where we don't require any continuity along the interface
178  }
179  }
180 
181  Assert (false, ExcNotImplemented());
183 }
184 
185 
186 
187 template <int dim, int spacedim>
188 void
191  FullMatrix<double> &interpolation_matrix) const
192 {
193  get_subface_interpolation_matrix (source_fe, numbers::invalid_unsigned_int,
194  interpolation_matrix);
195 }
196 
197 
198 
199 template <int dim, int spacedim>
200 void
203  const unsigned int subface,
204  FullMatrix<double> &interpolation_matrix) const
205 {
206  // this is the code from FE_FaceQ
207  Assert (interpolation_matrix.n() == this->dofs_per_face,
208  ExcDimensionMismatch (interpolation_matrix.n(),
209  this->dofs_per_face));
210  Assert (interpolation_matrix.m() == x_source_fe.dofs_per_face,
211  ExcDimensionMismatch (interpolation_matrix.m(),
212  x_source_fe.dofs_per_face));
213 
214  // see if source is a FaceQ element
215  if (const FE_TraceQ<dim,spacedim> *source_fe
216  = dynamic_cast<const FE_TraceQ<dim,spacedim> *>(&x_source_fe))
217  {
218  fe_q.get_subface_interpolation_matrix (source_fe->fe_q, subface, interpolation_matrix);
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 spacedim>
232 FE_TraceQ<1,spacedim>::FE_TraceQ (const unsigned int degree)
233  :
234  FE_FaceQ<1,spacedim> (degree)
235 {}
236 
237 
238 
239 template <int spacedim>
240 std::string
242 {
243  // note that the FETools::get_fe_by_name function depends on the
244  // particular format of the string this function returns, so they have to be
245  // kept in synch
246  std::ostringstream namebuf;
247  namebuf << "FE_TraceQ<"
248  << Utilities::dim_string(1,spacedim)
249  << ">(" << this->degree << ")";
250 
251  return namebuf.str();
252 }
253 
254 
255 
256 // explicit instantiations
257 #include "fe_trace.inst"
258 
259 
260 DEAL_II_NAMESPACE_CLOSE
size_type m() const
static const unsigned int invalid_unsigned_int
Definition: types.h:170
FullMatrix< double > interface_constraints
Definition: fe.h:2211
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_trace.cc:202
FE_TraceQ(unsigned int p)
Definition: fe_trace.cc:35
void hierarchic_to_lexicographic_numbering(unsigned int degree, std::vector< unsigned int > &h2l)
Definition: fe_tools.cc:2928
const unsigned int degree
Definition: fe_base.h:311
void set_numbering(const std::vector< unsigned int > &renumber)
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2223
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_trace.cc:190
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< Point< dim-1 > > unit_face_support_points
Definition: fe.h:2230
FE_Q< dim, spacedim > fe_q
Definition: fe_trace.h:130
virtual bool hp_constraints_are_implemented() const
Definition: fe_trace.cc:146
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_trace.cc:95
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:161
const unsigned int dofs_per_cell
Definition: fe_base.h:295
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_trace.cc:155
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_trace.cc:115
virtual FiniteElement< dim, spacedim > * clone() const
Definition: fe_trace.cc:68
const unsigned int dofs_per_face
Definition: fe_base.h:288
static ::ExceptionBase & ExcNotImplemented()
static std::vector< unsigned int > get_dpo_vector(const unsigned int deg)
Definition: fe_trace.cc:128
virtual std::string get_name() const
Definition: fe_trace.cc:77