Reference documentation for deal.II version 9.0.0
fe_series.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 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 
17 
18 #include <deal.II/fe/fe_series.h>
19 #include <deal.II/base/numbers.h>
20 
21 #include <cctype>
22 #include <iostream>
23 
24 #include <deal.II/base/config.h>
25 #ifdef DEAL_II_WITH_GSL
26 #include <gsl/gsl_sf_legendre.h>
27 #endif
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 namespace FESeries
33 {
34 
35  /*-------------- Fourier -------------------------------*/
36 
37  void set_k_vectors(Table<1, Tensor<1,1> > &k_vectors,
38  const unsigned int N)
39  {
40  k_vectors.reinit(TableIndices<1>(N));
41  for (unsigned int i=0; i<N; ++i)
42  k_vectors(i)[0] = 2. * numbers::PI * i;
43 
44  }
45 
46  void set_k_vectors(Table<2, Tensor<1,2> > &k_vectors,
47  const unsigned int N)
48  {
49  k_vectors.reinit(TableIndices<2>(N,N));
50  for (unsigned int i=0; i<N; ++i)
51  for (unsigned int j=0; j<N; ++j)
52  {
53  k_vectors(i,j)[0] = 2. * numbers::PI * i;
54  k_vectors(i,j)[1] = 2. * numbers::PI * j;
55  }
56  }
57 
58  void set_k_vectors(Table<3, Tensor<1,3> > &k_vectors,
59  const unsigned int N)
60  {
61  k_vectors.reinit(TableIndices<3>(N,N,N));
62  for (unsigned int i=0; i<N; ++i)
63  for (unsigned int j=0; j<N; ++j)
64  for (unsigned int k=0; k<N; ++k)
65  {
66  k_vectors(i,j,k)[0] = 2. * numbers::PI * i;
67  k_vectors(i,j,k)[0] = 2. * numbers::PI * j;
68  k_vectors(i,j,k)[0] = 2. * numbers::PI * k;
69  }
70  }
71 
72 
73  template <int dim>
74  Fourier<dim>::Fourier(const unsigned int N,
75  const hp::FECollection<dim> &fe_collection,
76  const hp::QCollection<dim> &q_collection)
77  :
78  fe_collection(&fe_collection),
79  q_collection(&q_collection),
80  fourier_transform_matrices(fe_collection.size())
81  {
82  set_k_vectors(k_vectors,N);
84  }
85 
86  template <int dim>
87  void Fourier<dim>::calculate(const Vector<double> &local_dof_values,
88  const unsigned int cell_active_fe_index,
89  Table<dim,std::complex<double> > &fourier_coefficients)
90  {
91  ensure_existence(cell_active_fe_index);
92  const FullMatrix<std::complex<double> > &matrix = fourier_transform_matrices[cell_active_fe_index];
93 
94  std::fill(unrolled_coefficients.begin(),
95  unrolled_coefficients.end(),
96  std::complex<double>(0.));
97 
98  Assert (unrolled_coefficients.size() == matrix.m(),
100 
101  Assert (local_dof_values.size() == matrix.n(),
102  ExcDimensionMismatch(local_dof_values.size(),matrix.n()));
103 
104  for (unsigned int i = 0; i < unrolled_coefficients.size(); i++)
105  for (unsigned int j = 0; j < local_dof_values.size(); j++)
106  unrolled_coefficients[i] += matrix[i][j] *
107  local_dof_values[j];
108 
109  fourier_coefficients.fill(unrolled_coefficients.begin());
110  }
111 
112  template <int dim>
113  std::complex<double> integrate(const FiniteElement<dim> &fe,
114  const Quadrature<dim> &quadrature,
115  const Tensor<1,dim> &k_vector,
116  const unsigned int j)
117  {
118  std::complex<double> sum = 0;
119  for (unsigned int q=0; q<quadrature.size(); ++q)
120  {
121  const Point<dim> &x_q = quadrature.point(q);
122  sum += std::exp(std::complex<double>(0,1) *
123  (k_vector * x_q)) *
124  fe.shape_value(j,x_q) *
125  quadrature.weight(q);
126  }
127  return sum;
128  }
129 
130 
131  template <>
132  void Fourier<2>::ensure_existence(const unsigned int fe)
133  {
134  Assert (fe < fe_collection->size(),
135  ExcIndexRange(fe,0,fe_collection->size()))
136 
137  if (fourier_transform_matrices[fe].m() == 0)
138  {
139  fourier_transform_matrices[fe].reinit(k_vectors.n_elements(),
140  (*fe_collection)[fe].dofs_per_cell);
141 
142  unsigned int k = 0;
143  for (unsigned int k1=0; k1<k_vectors.size(0); ++k1)
144  for (unsigned int k2=0; k2<k_vectors.size(1); ++k2,k++)
145  for (unsigned int j=0; j<(*fe_collection)[fe].dofs_per_cell; ++j)
146  fourier_transform_matrices[fe](k,j) = integrate((*fe_collection)[fe],
147  (*q_collection)[fe],
148  k_vectors(k1,k2),
149  j);
150  }
151  }
152 
153  template <>
154  void Fourier<3>::ensure_existence(const unsigned int fe)
155  {
156  Assert (fe < fe_collection->size(),
157  ExcIndexRange(fe,0,fe_collection->size()))
158 
159  if (fourier_transform_matrices[fe].m() == 0)
160  {
161  fourier_transform_matrices[fe].reinit(k_vectors.n_elements(),
162  (*fe_collection)[fe].dofs_per_cell);
163 
164  unsigned int k = 0;
165  for (unsigned int k1=0; k1<k_vectors.size(0); ++k1)
166  for (unsigned int k2=0; k2<k_vectors.size(1); ++k2)
167  for (unsigned int k3=0; k3<k_vectors.size(2); ++k3, k++)
168  for (unsigned int j=0; j<(*fe_collection)[fe].dofs_per_cell; ++j)
169  fourier_transform_matrices[fe](k,j) = integrate((*fe_collection)[fe],
170  (*q_collection)[fe],
171  k_vectors(k1,k2,k3),
172  j);
173  }
174  }
175 
176  template <>
177  void Fourier<1>::ensure_existence(const unsigned int fe)
178  {
179  Assert (fe < fe_collection->size(),
180  ExcIndexRange(fe,0,fe_collection->size()))
181 
182  if (fourier_transform_matrices[fe].m() == 0)
183  {
184  fourier_transform_matrices[fe].reinit(k_vectors.n_elements(),
185  (*fe_collection)[fe].dofs_per_cell);
186 
187  for (unsigned int k=0; k<k_vectors.size(0); ++k)
188  for (unsigned int j=0; j<(*fe_collection)[fe].dofs_per_cell; ++j)
189  fourier_transform_matrices[fe](k,j) = integrate((*fe_collection)[fe],
190  (*q_collection)[fe],
191  k_vectors(k),
192  j);
193  }
194  }
195 
196 
197  /*-------------- Legendre -------------------------------*/
198  DeclException2 (ExcLegendre, int, double,
199  << "x["<<arg1 << "] = "<<arg2 << " is not in [0,1]");
200 
201  /* dim dimensional Legendre function with indices @p indices
202  * evaluated at @p x_q in [0,1]^dim.
203  */
204  template <int dim>
205  double Lh(const Point<dim> &x_q,
206  const TableIndices<dim> &indices)
207  {
208 #ifdef DEAL_II_WITH_GSL
209  double res = 1.0;
210  for (unsigned int d = 0; d < dim; d++)
211  {
212  const double x = 2.0*(x_q[d]-0.5);
213  Assert ( (x_q[d] <= 1.0) && (x_q[d] >= 0.),
214  ExcLegendre(d,x_q[d]));
215  const int ind = indices[d];
216  res *= sqrt(2.0) * gsl_sf_legendre_Pl (ind, x);
217  }
218  return res;
219 
220 #else
221 
222  (void)x_q;
223  (void)indices;
224  AssertThrow(false, ExcMessage("deal.II has to be configured with GSL"
225  "in order to use Legendre transformation."));
226  return 0;
227 #endif
228  }
229 
230  /*
231  * Multiplier in Legendre coefficients
232  */
233  template <int dim>
234  double multiplier(const TableIndices<dim> &indices)
235  {
236  double res = 1.0;
237  for (unsigned int d = 0; d < dim; d++)
238  res *= (0.5+indices[d]);
239 
240  return res;
241  }
242 
243 
244  template <int dim>
245  Legendre<dim>::Legendre(const unsigned int size_in_each_direction,
246  const hp::FECollection<dim> &fe_collection,
247  const hp::QCollection<dim> &q_collection)
248  :
249  N(size_in_each_direction),
250  fe_collection(&fe_collection),
251  q_collection(&q_collection),
252  legendre_transform_matrices(fe_collection.size()),
253  unrolled_coefficients(Utilities::fixed_power<dim>(N),
254  0.)
255  {
256  }
257 
258  template <int dim>
259  void Legendre<dim>::calculate(const ::Vector<double> &local_dof_values,
260  const unsigned int cell_active_fe_index,
261  Table<dim,double> &legendre_coefficients)
262  {
263  ensure_existence(cell_active_fe_index);
264  const FullMatrix<double> &matrix = legendre_transform_matrices[cell_active_fe_index];
265 
266  std::fill(unrolled_coefficients.begin(),
267  unrolled_coefficients.end(),
268  0.);
269 
270  Assert (unrolled_coefficients.size() == matrix.m(),
271  ExcInternalError());
272 
273  Assert (local_dof_values.size() == matrix.n(),
274  ExcDimensionMismatch(local_dof_values.size(),matrix.n()));
275 
276  for (unsigned int i = 0; i < unrolled_coefficients.size(); i++)
277  for (unsigned int j = 0; j < local_dof_values.size(); j++)
278  unrolled_coefficients[i] += matrix[i][j] *
279  local_dof_values[j];
280 
281  legendre_coefficients.fill(unrolled_coefficients.begin());
282  }
283 
284 
285  template <int dim>
286  double integrate_Legendre(const FiniteElement<dim> &fe,
287  const Quadrature<dim> &quadrature,
288  const TableIndices<dim> &indices,
289  const unsigned int dof)
290  {
291  double sum = 0;
292  for (unsigned int q=0; q<quadrature.size(); ++q)
293  {
294  const Point<dim> &x_q = quadrature.point(q);
295  sum += Lh(x_q, indices) *
296  fe.shape_value(dof,x_q) *
297  quadrature.weight(q);
298  }
299  return sum * multiplier(indices);
300  }
301 
302  template <>
303  void Legendre<1>::ensure_existence(const unsigned int fe)
304  {
305  Assert (fe < fe_collection->size(),
306  ExcIndexRange(fe,0,fe_collection->size()))
307  if (legendre_transform_matrices[fe].m() == 0)
308  {
309  legendre_transform_matrices[fe].reinit(N,
310  (*fe_collection)[fe].dofs_per_cell);
311 
312  for (unsigned int k=0; k<N; ++k)
313  for (unsigned int j=0; j<(*fe_collection)[fe].dofs_per_cell; ++j)
314  legendre_transform_matrices[fe](k,j) =
315  integrate_Legendre((*fe_collection)[fe],
316  (*q_collection)[fe],
317  TableIndices<1>(k),
318  j);
319  }
320  }
321 
322 
323 
324  template <>
325  void Legendre<2>::ensure_existence(const unsigned int fe)
326  {
327  Assert (fe < fe_collection->size(),
328  ExcIndexRange(fe,0,fe_collection->size()))
329 
330  if (legendre_transform_matrices[fe].m() == 0)
331  {
332  legendre_transform_matrices[fe].reinit(N*N,
333  (*fe_collection)[fe].dofs_per_cell);
334 
335  unsigned int k = 0;
336  for (unsigned int k1=0; k1<N; ++k1)
337  for (unsigned int k2=0; k2<N; ++k2,k++)
338  for (unsigned int j=0; j<(*fe_collection)[fe].dofs_per_cell; ++j)
339  legendre_transform_matrices[fe](k,j) =
340  integrate_Legendre((*fe_collection)[fe],
341  (*q_collection)[fe],
342  TableIndices<2>(k1,k2),
343  j);
344  }
345  }
346 
347  template <>
348  void Legendre<3>::ensure_existence(const unsigned int fe)
349  {
350  Assert (fe < fe_collection->size(),
351  ExcIndexRange(fe,0,fe_collection->size()))
352 
353  if (legendre_transform_matrices[fe].m() == 0)
354  {
355  legendre_transform_matrices[fe].reinit(N*N*N,
356  (*fe_collection)[fe].dofs_per_cell);
357 
358  unsigned int k = 0;
359  for (unsigned int k1=0; k1<N; ++k1)
360  for (unsigned int k2=0; k2<N; ++k2)
361  for (unsigned int k3=0; k3<N; ++k3, k++)
362  for (unsigned int j=0; j<(*fe_collection)[fe].dofs_per_cell; ++j)
363  legendre_transform_matrices[fe](k,j) =
364  integrate_Legendre((*fe_collection)[fe],
365  (*q_collection)[fe],
366  TableIndices<3>(k1,k2,k3),
367  j);
368  }
369  }
370 
371  /*-------------- linear_regression -------------------------------*/
372  std::pair<double,double> linear_regression(const std::vector<double> &x,
373  const std::vector<double> &y)
374  {
375  FullMatrix<double> K(2,2), invK(2,2);
376  Vector<double> X(2), B(2);
377 
378  Assert (x.size() == y.size(),
379  ExcMessage("x and y are expected to have the same size"));
380 
381  Assert (x.size() >= 2,
382  ::ExcMessage("at least two points are required for linear regression fit"));
383 
384  double sum_1 = 0.0,
385  sum_x = 0.0,
386  sum_x2= 0.0,
387  sum_y = 0.0,
388  sum_xy= 0.0;
389 
390  for (unsigned int i = 0; i < x.size(); i++)
391  {
392  sum_1 += 1.0;
393  sum_x += x[i];
394  sum_x2 += x[i]*x[i];
395  sum_y += y[i];
396  sum_xy += x[i]*y[i];
397  }
398 
399  K(0,0) = sum_1;
400  K(0,1) = sum_x;
401  K(1,0) = sum_x;
402  K(1,1) = sum_x2;
403 
404  B(0) = sum_y;
405  B(1) = sum_xy;
406 
407  invK.invert(K);
408  invK.vmult(X,B,false);
409 
410  return std::make_pair(X(1),X(0));
411  }
412 
413 
414 } // end of namespace FESeries
415 
416 
417 
418 /*-------------- Explicit Instantiations -------------------------------*/
419 template class FESeries::Fourier<1>;
420 template class FESeries::Fourier<2>;
421 template class FESeries::Fourier<3>;
422 template class FESeries::Legendre<1>;
423 template class FESeries::Legendre<2>;
424 template class FESeries::Legendre<3>;
425 
426 
427 
428 DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcLegendre(int arg1, double arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:358
void ensure_existence(const unsigned int fe_index)
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
const Point< dim > & point(const unsigned int i) const
Legendre(const unsigned int size_in_each_direction, const hp::FECollection< dim > &fe_collection, const hp::QCollection< dim > &q_collection)
Definition: fe_series.cc:245
std::pair< double, double > linear_regression(const std::vector< double > &x, const std::vector< double > &y)
Definition: fe_series.cc:372
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
void calculate(const ::Vector< double > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, double > &legendre_coefficients)
Definition: fe_series.cc:259
Fourier(const unsigned int size_in_each_direction, const hp::FECollection< dim > &fe_collection, const hp::QCollection< dim > &q_collection)
Definition: fe_series.cc:74
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
size_type n_elements() const
std::vector< std::complex< double > > unrolled_coefficients
Definition: fe_series.h:140
static const double PI
Definition: numbers.h:127
static ::ExceptionBase & ExcMessage(std::string arg1)
void invert(const FullMatrix< number2 > &M)
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void calculate(const ::Vector< double > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, std::complex< double > > &fourier_coefficients)
Definition: fe_series.cc:87
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void ensure_existence(const unsigned int fe_index)
unsigned int size() const
Definition: cuda.h:31
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:160
Definition: mpi.h:53
size_type size() const
Definition: table.h:34
Table< dim, Tensor< 1, dim > > k_vectors
Definition: fe_series.h:130
void fill(InputIterator entries, const bool C_style_indexing=true)
double weight(const unsigned int i) const
static ::ExceptionBase & ExcInternalError()