Reference documentation for deal.II version 8.5.1
polynomials_raviart_thomas.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 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 
17 #include <deal.II/base/polynomials_raviart_thomas.h>
18 #include <deal.II/base/quadrature_lib.h>
19 #include <deal.II/base/thread_management.h>
20 #include <iostream>
21 #include <iomanip>
22 
23 //TODO[WB]: This class is not thread-safe: it uses mutable member variables that contain temporary state. this is not what one would want when one uses a finite element object in a number of different contexts on different threads: finite element objects should be stateless
24 //TODO:[GK] This can be achieved by writing a function in Polynomial space which does the rotated fill performed below and writes the data into the right data structures. The same function would be used
25 //by Nedelec polynomials.
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 
30 template <int dim>
32  :
33  my_degree(k),
34  polynomial_space (create_polynomials (k)),
35  n_pols(compute_n_pols(k))
36 {}
37 
38 
39 
40 template <int dim>
41 std::vector<std::vector< Polynomials::Polynomial< double > > >
43 {
44  std::vector<std::vector< Polynomials::Polynomial< double > > > pols(dim);
46  if (k == 0)
47  for (unsigned int d=1; d<dim; ++d)
49  else
50  for (unsigned int d=1; d<dim; ++d)
52 
53  return pols;
54 }
55 
56 
57 template <int dim>
58 void
60  std::vector<Tensor<1,dim> > &values,
61  std::vector<Tensor<2,dim> > &grads,
62  std::vector<Tensor<3, dim> > &grad_grads,
63  std::vector<Tensor<4, dim> > &third_derivatives,
64  std::vector<Tensor<5, dim> > &fourth_derivatives) const
65 {
66  Assert(values.size()==n_pols || values.size()==0,
67  ExcDimensionMismatch(values.size(), n_pols));
68  Assert(grads.size()==n_pols|| grads.size()==0,
69  ExcDimensionMismatch(grads.size(), n_pols));
70  Assert(grad_grads.size()==n_pols|| grad_grads.size()==0,
71  ExcDimensionMismatch(grad_grads.size(), n_pols));
72  Assert(third_derivatives.size()==n_pols|| third_derivatives.size()==0,
73  ExcDimensionMismatch(third_derivatives.size(), n_pols));
74  Assert(fourth_derivatives.size()==n_pols|| fourth_derivatives.size()==0,
75  ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
76 
77  // have a few scratch
78  // arrays. because we don't want to
79  // re-allocate them every time this
80  // function is called, we make them
81  // static. however, in return we
82  // have to ensure that the calls to
83  // the use of these variables is
84  // locked with a mutex. if the
85  // mutex is removed, several tests
86  // (notably
87  // deal.II/create_mass_matrix_05)
88  // will start to produce random
89  // results in multithread mode
90  static Threads::Mutex mutex;
91  Threads::Mutex::ScopedLock lock(mutex);
92 
93  static std::vector<double> p_values;
94  static std::vector<Tensor<1,dim> > p_grads;
95  static std::vector<Tensor<2,dim> > p_grad_grads;
96  static std::vector<Tensor<3,dim> > p_third_derivatives;
97  static std::vector<Tensor<4,dim> > p_fourth_derivatives;
98 
99  const unsigned int n_sub = polynomial_space.n();
100  p_values.resize((values.size() == 0) ? 0 : n_sub);
101  p_grads.resize((grads.size() == 0) ? 0 : n_sub);
102  p_grad_grads.resize((grad_grads.size() == 0) ? 0 : n_sub);
103  p_third_derivatives.resize((third_derivatives.size() == 0) ? 0 : n_sub);
104  p_fourth_derivatives.resize((fourth_derivatives.size() == 0) ? 0 : n_sub);
105 
106  for (unsigned int d=0; d<dim; ++d)
107  {
108  // First we copy the point. The
109  // polynomial space for
110  // component d consists of
111  // polynomials of degree k+1 in
112  // x_d and degree k in the
113  // other variables. in order to
114  // simplify this, we use the
115  // same AnisotropicPolynomial
116  // space and simply rotate the
117  // coordinates through all
118  // directions.
119  Point<dim> p;
120  for (unsigned int c=0; c<dim; ++c)
121  p(c) = unit_point((c+d)%dim);
122 
123  polynomial_space.compute (p, p_values, p_grads, p_grad_grads, p_third_derivatives, p_fourth_derivatives);
124 
125  for (unsigned int i=0; i<p_values.size(); ++i)
126  values[i+d*n_sub][d] = p_values[i];
127 
128  for (unsigned int i=0; i<p_grads.size(); ++i)
129  for (unsigned int d1=0; d1<dim; ++d1)
130  grads[i+d*n_sub][d][(d1+d)%dim] = p_grads[i][d1];
131 
132  for (unsigned int i=0; i<p_grad_grads.size(); ++i)
133  for (unsigned int d1=0; d1<dim; ++d1)
134  for (unsigned int d2=0; d2<dim; ++d2)
135  grad_grads[i+d*n_sub][d][(d1+d)%dim][(d2+d)%dim]
136  = p_grad_grads[i][d1][d2];
137 
138  for (unsigned int i=0; i<p_third_derivatives.size(); ++i)
139  for (unsigned int d1=0; d1<dim; ++d1)
140  for (unsigned int d2=0; d2<dim; ++d2)
141  for (unsigned int d3=0; d3<dim; ++d3)
142  third_derivatives[i+d*n_sub][d][(d1+d)%dim][(d2+d)%dim][(d3+d)%dim]
143  = p_third_derivatives[i][d1][d2][d3];
144 
145  for (unsigned int i=0; i<p_fourth_derivatives.size(); ++i)
146  for (unsigned int d1=0; d1<dim; ++d1)
147  for (unsigned int d2=0; d2<dim; ++d2)
148  for (unsigned int d3=0; d3<dim; ++d3)
149  for (unsigned int d4=0; d4<dim; ++d4)
150  fourth_derivatives[i+d*n_sub][d][(d1+d)%dim][(d2+d)%dim][(d3+d)%dim][(d4+d)%dim]
151  = p_fourth_derivatives[i][d1][d2][d3][d4];
152  }
153 }
154 
155 
156 template <int dim>
157 unsigned int
159 {
160  if (dim == 1) return k+1;
161  if (dim == 2) return 2*(k+1)*(k+2);
162  if (dim == 3) return 3*(k+1)*(k+1)*(k+2);
163 
164  Assert(false, ExcNotImplemented());
165  return 0;
166 }
167 
168 
169 template class PolynomialsRaviartThomas<1>;
170 template class PolynomialsRaviartThomas<2>;
171 template class PolynomialsRaviartThomas<3>;
172 
173 
174 DEAL_II_NAMESPACE_CLOSE
PolynomialsRaviartThomas(const unsigned int k)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:873
static unsigned int compute_n_pols(unsigned int degree)
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:806
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void compute(const Point< dim > &unit_point, std::vector< Tensor< 1, dim > > &values, std::vector< Tensor< 2, dim > > &grads, std::vector< Tensor< 3, dim > > &grad_grads, std::vector< Tensor< 4, dim > > &third_derivatives, std::vector< Tensor< 5, dim > > &fourth_derivatives) const
static std::vector< std::vector< Polynomials::Polynomial< double > > > create_polynomials(const unsigned int k)
static ::ExceptionBase & ExcNotImplemented()