Reference documentation for deal.II version 9.0.0
polynomials_abf.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 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/base/polynomials_abf.h>
18 #include <deal.II/base/quadrature_lib.h>
19 #include <iostream>
20 #include <iomanip>
21 
22 
23 DEAL_II_NAMESPACE_OPEN
24 
25 
26 
27 namespace
28 {
29  template <int dim>
30  std::vector<std::vector< Polynomials::Polynomial< double > > >
31  get_abf_polynomials (const unsigned int k)
32  {
33  std::vector<std::vector< Polynomials::Polynomial< double > > > pols(dim);
35 
36  if (k == 0)
37  for (unsigned int d=1; d<dim; ++d)
39  else
40  for (unsigned int d=1; d<dim; ++d)
42 
43  return pols;
44  }
45 }
46 
47 template <int dim>
49  :
50  my_degree(k),
51  polynomial_space(get_abf_polynomials<dim>(k)),
52  n_pols(compute_n_pols(k))
53 {
54  // check that the dimensions match. we only store one of the 'dim'
55  // anisotropic polynomials that make up the vector-valued space, so
56  // multiply by 'dim'
57  Assert (dim * polynomial_space.n() == compute_n_pols(k),
59 }
60 
61 
62 
63 template <int dim>
64 void
66  std::vector<Tensor<1,dim> > &values,
67  std::vector<Tensor<2,dim> > &grads,
68  std::vector<Tensor<3,dim> > &grad_grads,
69  std::vector<Tensor<4,dim> > &third_derivatives,
70  std::vector<Tensor<5,dim> > &fourth_derivatives) const
71 {
72  Assert(values.size()==n_pols || values.size()==0,
73  ExcDimensionMismatch(values.size(), n_pols));
74  Assert(grads.size()==n_pols|| grads.size()==0,
75  ExcDimensionMismatch(grads.size(), n_pols));
76  Assert(grad_grads.size()==n_pols|| grad_grads.size()==0,
77  ExcDimensionMismatch(grad_grads.size(), n_pols));
78  Assert(third_derivatives.size()==n_pols|| third_derivatives.size()==0,
79  ExcDimensionMismatch(third_derivatives.size(), n_pols));
80  Assert(fourth_derivatives.size()==n_pols|| fourth_derivatives.size()==0,
81  ExcDimensionMismatch(fourth_derivatives.size(), n_pols));
82 
83  const unsigned int n_sub = polynomial_space.n();
84  // guard access to the scratch
85  // arrays in the following block
86  // using a mutex to make sure they
87  // are not used by multiple threads
88  // at once
89  Threads::Mutex::ScopedLock lock(mutex);
90 
91  p_values.resize((values.size() == 0) ? 0 : n_sub);
92  p_grads.resize((grads.size() == 0) ? 0 : n_sub);
93  p_grad_grads.resize((grad_grads.size() == 0) ? 0 : n_sub);
94  p_third_derivatives.resize((third_derivatives.size() == 0) ? 0 : n_sub);
95  p_fourth_derivatives.resize((fourth_derivatives.size() == 0) ? 0 : n_sub);
96 
97  for (unsigned int d=0; d<dim; ++d)
98  {
99  // First we copy the point. The
100  // polynomial space for
101  // component d consists of
102  // polynomials of degree k+1 in
103  // x_d and degree k in the
104  // other variables. in order to
105  // simplify this, we use the
106  // same AnisotropicPolynomial
107  // space and simply rotate the
108  // coordinates through all
109  // directions.
110  Point<dim> p;
111  for (unsigned int c=0; c<dim; ++c)
112  p(c) = unit_point((c+d)%dim);
113 
114  polynomial_space.compute (p, p_values, p_grads, p_grad_grads,
115  p_third_derivatives, p_fourth_derivatives);
116 
117  for (unsigned int i=0; i<p_values.size(); ++i)
118  values[i+d*n_sub][d] = p_values[i];
119 
120  for (unsigned int i=0; i<p_grads.size(); ++i)
121  for (unsigned int d1=0; d1<dim; ++d1)
122  grads[i+d*n_sub][d][(d1+d)%dim] = p_grads[i][d1];
123 
124  for (unsigned int i=0; i<p_grad_grads.size(); ++i)
125  for (unsigned int d1=0; d1<dim; ++d1)
126  for (unsigned int d2=0; d2<dim; ++d2)
127  grad_grads[i+d*n_sub][d][(d1+d)%dim][(d2+d)%dim]
128  = p_grad_grads[i][d1][d2];
129 
130  for (unsigned int i=0; i<p_third_derivatives.size(); ++i)
131  for (unsigned int d1=0; d1<dim; ++d1)
132  for (unsigned int d2=0; d2<dim; ++d2)
133  for (unsigned int d3=0; d3<dim; ++d3)
134  third_derivatives[i+d*n_sub][d][(d1+d)%dim][(d2+d)%dim][(d3+d)%dim]
135  = p_third_derivatives[i][d1][d2][d3];
136 
137  for (unsigned int i=0; i<p_fourth_derivatives.size(); ++i)
138  for (unsigned int d1=0; d1<dim; ++d1)
139  for (unsigned int d2=0; d2<dim; ++d2)
140  for (unsigned int d3=0; d3<dim; ++d3)
141  for (unsigned int d4=0; d4<dim; ++d4)
142  fourth_derivatives[i+d*n_sub][d][(d1+d)%dim][(d2+d)%dim][(d3+d)%dim][(d4+d)%dim]
143  = p_fourth_derivatives[i][d1][d2][d3][d4];
144  }
145 }
146 
147 
148 template <int dim>
149 unsigned int
151 {
152  switch (dim)
153  {
154  case 1:
155  // in 1d, we simply have Q_{k+2}, which has dimension k+3
156  return k+3;
157 
158  case 2:
159  // the polynomial space is Q_{k+2,k} \times Q_{k,k+2}, which has
160  // 2(k+3)(k+1) DoFs
161  return 2*(k+3)*(k+1);
162 
163  case 3:
164  // the polynomial space is Q_{k+2,k,k} \times Q_{k,k+2,k} \times Q_{k,k,k+2},
165  // which has 3(k+3)(k+1)(k+1) DoFs
166  return 3*(k+3)*(k+1)*(k+1);
167 
168  default:
169  Assert(false, ExcNotImplemented());
170  }
171 
172  return 0;
173 }
174 
175 
176 template class PolynomialsABF<1>;
177 template class PolynomialsABF<2>;
178 template class PolynomialsABF<3>;
179 
180 
181 DEAL_II_NAMESPACE_CLOSE
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:874
static std::vector< Polynomial< double > > generate_complete_basis(const unsigned int degree)
Definition: polynomial.cc:807
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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 unsigned int compute_n_pols(unsigned int degree)
static ::ExceptionBase & ExcNotImplemented()
const AnisotropicPolynomials< dim > polynomial_space
PolynomialsABF(const unsigned int k)
static ::ExceptionBase & ExcInternalError()