Reference documentation for deal.II version 9.0.0
fe_dgp.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/fe/fe_dgp.h>
18 #include <deal.II/fe/fe_tools.h>
19 
20 #include <sstream>
21 #include <deal.II/base/std_cxx14/memory.h>
22 
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 template <int dim, int spacedim>
27 FE_DGP<dim,spacedim>::FE_DGP (const unsigned int degree)
28  :
29  FE_Poly<PolynomialSpace<dim>, dim, spacedim> (
30  PolynomialSpace<dim>(Polynomials::Legendre::generate_complete_basis(degree)),
31  FiniteElementData<dim>(get_dpo_vector(degree), 1, degree, FiniteElementData<dim>::L2),
32  std::vector<bool>(FiniteElementData<dim>(get_dpo_vector(degree), 1, degree).dofs_per_cell,true),
33  std::vector<ComponentMask>(FiniteElementData<dim>(
34  get_dpo_vector(degree), 1, degree).dofs_per_cell, std::vector<bool>(1,true)))
35 {
36  // Reinit the vectors of restriction and prolongation matrices to the right
37  // sizes
39  // Fill prolongation matrices with embedding operators
40  if (dim == spacedim)
41  {
43  // Fill restriction matrices with L2-projection
45  }
46 }
47 
48 
49 template <int dim, int spacedim>
50 std::string
52 {
53  // note that the FETools::get_fe_by_name function depends on the
54  // particular format of the string this function returns, so they have to be
55  // kept in sync
56 
57  std::ostringstream namebuf;
58  namebuf << "FE_DGP<"
59  << Utilities::dim_string(dim,spacedim)
60  << ">(" << this->degree << ")";
61 
62  return namebuf.str();
63 }
64 
65 
66 
67 template <int dim, int spacedim>
68 std::unique_ptr<FiniteElement<dim,spacedim> >
70 {
71  return std_cxx14::make_unique<FE_DGP<dim,spacedim>>(*this);
72 }
73 
74 
75 
76 //---------------------------------------------------------------------------
77 // Auxiliary functions
78 //---------------------------------------------------------------------------
79 
80 
81 template <int dim, int spacedim>
82 std::vector<unsigned int>
83 FE_DGP<dim,spacedim>::get_dpo_vector (const unsigned int deg)
84 {
85  std::vector<unsigned int> dpo(dim+1, 0U);
86  dpo[dim] = deg+1;
87  for (unsigned int i=1; i<dim; ++i)
88  {
89  dpo[dim] *= deg+1+i;
90  dpo[dim] /= i+1;
91  }
92  return dpo;
93 }
94 
95 
96 
97 template <int dim, int spacedim>
98 void
101  FullMatrix<double> &interpolation_matrix) const
102 {
103  // this is only implemented, if the source FE is also a DGP element. in that
104  // case, both elements have no dofs on their faces and the face
105  // interpolation matrix is necessarily empty -- i.e. there isn't much we
106  // need to do here.
107  (void)interpolation_matrix;
108  typedef FiniteElement<dim,spacedim> FE;
109  typedef FE_DGP<dim,spacedim> FEDGP;
110  AssertThrow ((x_source_fe.get_name().find ("FE_DGP<") == 0)
111  ||
112  (dynamic_cast<const FEDGP *>(&x_source_fe) != nullptr),
113  typename FE::
114  ExcInterpolationNotImplemented());
115 
116  Assert (interpolation_matrix.m() == 0,
117  ExcDimensionMismatch (interpolation_matrix.m(),
118  0));
119  Assert (interpolation_matrix.n() == 0,
120  ExcDimensionMismatch (interpolation_matrix.n(),
121  0));
122 }
123 
124 
125 
126 template <int dim, int spacedim>
127 void
130  const unsigned int,
131  FullMatrix<double> &interpolation_matrix) const
132 {
133  // this is only implemented, if the source FE is also a DGP element. in that
134  // case, both elements have no dofs on their faces and the face
135  // interpolation matrix is necessarily empty -- i.e. there isn't much we
136  // need to do here.
137  (void)interpolation_matrix;
138  typedef FiniteElement<dim,spacedim> FE;
139  typedef FE_DGP<dim,spacedim> FEDGP;
140  AssertThrow ((x_source_fe.get_name().find ("FE_DGP<") == 0)
141  ||
142  (dynamic_cast<const FEDGP *>(&x_source_fe) != nullptr),
143  typename FE::
144  ExcInterpolationNotImplemented());
145 
146  Assert (interpolation_matrix.m() == 0,
147  ExcDimensionMismatch (interpolation_matrix.m(),
148  0));
149  Assert (interpolation_matrix.n() == 0,
150  ExcDimensionMismatch (interpolation_matrix.n(),
151  0));
152 }
153 
154 
155 
156 template <int dim, int spacedim>
157 bool
159 {
160  return true;
161 }
162 
163 
164 
165 template <int dim, int spacedim>
166 std::vector<std::pair<unsigned int, unsigned int> >
169 {
170  // there are no such constraints for DGP elements at all
171  if (dynamic_cast<const FE_DGP<dim,spacedim>*>(&fe_other) != nullptr)
172  return
173  std::vector<std::pair<unsigned int, unsigned int> > ();
174  else
175  {
176  Assert (false, ExcNotImplemented());
177  return std::vector<std::pair<unsigned int, unsigned int> > ();
178  }
179 }
180 
181 
182 
183 template <int dim, int spacedim>
184 std::vector<std::pair<unsigned int, unsigned int> >
187 {
188  // there are no such constraints for DGP elements at all
189  if (dynamic_cast<const FE_DGP<dim,spacedim>*>(&fe_other) != nullptr)
190  return
191  std::vector<std::pair<unsigned int, unsigned int> > ();
192  else
193  {
194  Assert (false, ExcNotImplemented());
195  return std::vector<std::pair<unsigned int, unsigned int> > ();
196  }
197 }
198 
199 
200 
201 template <int dim, int spacedim>
202 std::vector<std::pair<unsigned int, unsigned int> >
205 {
206  // there are no such constraints for DGP elements at all
207  if (dynamic_cast<const FE_DGP<dim,spacedim>*>(&fe_other) != nullptr)
208  return
209  std::vector<std::pair<unsigned int, unsigned int> > ();
210  else
211  {
212  Assert (false, ExcNotImplemented());
213  return std::vector<std::pair<unsigned int, unsigned int> > ();
214  }
215 }
216 
217 
218 
219 template <int dim, int spacedim>
222 {
223  // check whether both are discontinuous elements, see the description of
224  // FiniteElementDomination::Domination
225  if (dynamic_cast<const FE_DGP<dim,spacedim>*>(&fe_other) != nullptr)
227 
228  Assert (false, ExcNotImplemented());
230 }
231 
232 
233 
234 template <int dim, int spacedim>
235 bool
237  const unsigned int) const
238 {
239  // all shape functions have support on all faces
240  return true;
241 }
242 
243 
244 
245 template <int dim, int spacedim>
246 std::pair<Table<2,bool>, std::vector<unsigned int> >
248 {
249  Table<2,bool> constant_modes(1, this->dofs_per_cell);
250  constant_modes(0,0) = true;
251  return std::pair<Table<2,bool>, std::vector<unsigned int> >
252  (constant_modes, std::vector<unsigned int>(1, 0));
253 }
254 
255 
256 
257 template <int dim, int spacedim>
258 std::size_t
260 {
261  Assert (false, ExcNotImplemented ());
262  return 0;
263 }
264 
265 
266 
267 // explicit instantiations
268 #include "fe_dgp.inst"
269 
270 
271 DEAL_II_NAMESPACE_CLOSE
size_type m() const
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2375
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe_dgp.cc:100
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgp.cc:168
virtual std::string get_name() const
Definition: fe_dgp.cc:51
virtual std::size_t memory_consumption() const
Definition: fe_dgp.cc:259
FE_DGP(const unsigned int p)
Definition: fe_dgp.cc:27
STL namespace.
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)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgp.cc:221
Definition: fe_dgp.h:309
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe_dgp.cc:236
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgp.cc:186
size_type n() const
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe_dgp.cc:129
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe_dgp.cc:247
static std::vector< unsigned int > get_dpo_vector(const unsigned int degree)
Definition: fe_dgp.cc:83
virtual std::string get_name() const =0
std::string dim_string(const int dim, const int spacedim)
Definition: utilities.cc:176
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const
Definition: fe_dgp.cc:69
void compute_projection_matrices(const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false)
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:274
static ::ExceptionBase & ExcNotImplemented()
virtual bool hp_constraints_are_implemented() const
Definition: fe_dgp.cc:158
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe_dgp.cc:204