Reference documentation for deal.II version 8.5.1
mapping_q_generic.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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/derivative_form.h>
18 #include <deal.II/base/quadrature.h>
19 #include <deal.II/base/qprojector.h>
20 #include <deal.II/base/quadrature_lib.h>
21 #include <deal.II/base/tensor_product_polynomials.h>
22 #include <deal.II/base/memory_consumption.h>
23 #include <deal.II/base/std_cxx11/array.h>
24 #include <deal.II/base/std_cxx11/unique_ptr.h>
25 #include <deal.II/lac/full_matrix.h>
26 #include <deal.II/grid/tria.h>
27 #include <deal.II/grid/tria_iterator.h>
28 #include <deal.II/grid/tria_boundary.h>
29 #include <deal.II/dofs/dof_accessor.h>
30 #include <deal.II/fe/fe_tools.h>
31 #include <deal.II/fe/fe.h>
32 #include <deal.II/fe/fe_values.h>
33 #include <deal.II/fe/mapping_q_generic.h>
34 #include <deal.II/fe/mapping_q1.h>
35 
36 #include <cmath>
37 #include <algorithm>
38 #include <numeric>
39 
40 
41 DEAL_II_NAMESPACE_OPEN
42 
43 
44 namespace
45 {
46  template <int dim>
47  std::vector<unsigned int>
48  get_dpo_vector (const unsigned int degree)
49  {
50  std::vector<unsigned int> dpo(dim+1, 1U);
51  for (unsigned int i=1; i<dpo.size(); ++i)
52  dpo[i]=dpo[i-1]*(degree-1);
53  return dpo;
54  }
55 }
56 
57 namespace internal
58 {
59  namespace MappingQ1
60  {
61  namespace
62  {
63 
64  // These are left as templates on the spatial dimension (even though dim
65  // == spacedim must be true for them to make sense) because templates are
66  // expanded before the compiler eliminates code due to the 'if (dim ==
67  // spacedim)' statement (see the body of the general
68  // transform_real_to_unit_cell).
69  template<int spacedim>
70  Point<1>
71  transform_real_to_unit_cell
72  (const std_cxx11::array<Point<spacedim>, GeometryInfo<1>::vertices_per_cell> &vertices,
73  const Point<spacedim> &p)
74  {
75  Assert(spacedim == 1, ExcInternalError());
76  return Point<1>((p[0] - vertices[0](0))/(vertices[1](0) - vertices[0](0)));
77  }
78 
79 
80 
81  template<int spacedim>
82  Point<2>
83  transform_real_to_unit_cell
84  (const std_cxx11::array<Point<spacedim>, GeometryInfo<2>::vertices_per_cell> &vertices,
85  const Point<spacedim> &p)
86  {
87  Assert(spacedim == 2, ExcInternalError());
88  const double x = p(0);
89  const double y = p(1);
90 
91  const double x0 = vertices[0](0);
92  const double x1 = vertices[1](0);
93  const double x2 = vertices[2](0);
94  const double x3 = vertices[3](0);
95 
96  const double y0 = vertices[0](1);
97  const double y1 = vertices[1](1);
98  const double y2 = vertices[2](1);
99  const double y3 = vertices[3](1);
100 
101  const double a = (x1 - x3)*(y0 - y2) - (x0 - x2)*(y1 - y3);
102  const double b = -(x0 - x1 - x2 + x3)*y + (x - 2*x1 + x3)*y0 - (x - 2*x0 + x2)*y1
103  - (x - x1)*y2 + (x - x0)*y3;
104  const double c = (x0 - x1)*y - (x - x1)*y0 + (x - x0)*y1;
105 
106  const double discriminant = b*b - 4*a*c;
107  // exit if the point is not in the cell (this is the only case where the
108  // discriminant is negative)
109  if (discriminant < 0.0)
110  {
111  AssertThrow (false,
113  }
114 
115  double eta1;
116  double eta2;
117  // special case #1: if a is zero, then use the linear formula
118  if (a == 0.0 && b != 0.0)
119  {
120  eta1 = -c/b;
121  eta2 = -c/b;
122  }
123  // special case #2: if c is very small or the square root of the
124  // discriminant is nearly b.
125  else if (std::abs(c) < 1e-12*std::abs(b)
126  || std::abs(std::sqrt(discriminant) - b) <= 1e-14*std::abs(b))
127  {
128  eta1 = (-b - std::sqrt(discriminant)) / (2*a);
129  eta2 = (-b + std::sqrt(discriminant)) / (2*a);
130  }
131  // finally, use the numerically stable version of the quadratic formula:
132  else
133  {
134  eta1 = 2*c / (-b - std::sqrt(discriminant));
135  eta2 = 2*c / (-b + std::sqrt(discriminant));
136  }
137  // pick the one closer to the center of the cell.
138  const double eta = (std::abs(eta1 - 0.5) < std::abs(eta2 - 0.5)) ? eta1 : eta2;
139 
140  /*
141  * There are two ways to compute xi from eta, but either one may have a
142  * zero denominator.
143  */
144  const double subexpr0 = -eta*x2 + x0*(eta - 1);
145  const double xi_denominator0 = eta*x3 - x1*(eta - 1) + subexpr0;
146  const double max_x = std::max(std::max(std::abs(x0), std::abs(x1)),
147  std::max(std::abs(x2), std::abs(x3)));
148 
149  if (std::abs(xi_denominator0) > 1e-10*max_x)
150  {
151  const double xi = (x + subexpr0)/xi_denominator0;
152  return Point<2>(xi, eta);
153  }
154  else
155  {
156  const double max_y = std::max(std::max(std::abs(y0), std::abs(y1)),
157  std::max(std::abs(y2), std::abs(y3)));
158  const double subexpr1 = -eta*y2 + y0*(eta - 1);
159  const double xi_denominator1 = eta*y3 - y1*(eta - 1) + subexpr1;
160  if (std::abs(xi_denominator1) > 1e-10*max_y)
161  {
162  const double xi = (subexpr1 + y)/xi_denominator1;
163  return Point<2>(xi, eta);
164  }
165  else // give up and try Newton iteration
166  {
167  AssertThrow (false,
169  }
170  }
171  // bogus return to placate compiler. It should not be possible to get
172  // here.
173  Assert(false, ExcInternalError());
174  return Point<2>(std::numeric_limits<double>::quiet_NaN(),
175  std::numeric_limits<double>::quiet_NaN());
176  }
177 
178 
179 
180  template<int spacedim>
181  Point<3>
182  transform_real_to_unit_cell
183  (const std_cxx11::array<Point<spacedim>, GeometryInfo<3>::vertices_per_cell> &/*vertices*/,
184  const Point<spacedim> &/*p*/)
185  {
186  // It should not be possible to get here
187  Assert(false, ExcInternalError());
188  return Point<3>();
189  }
190 
191 
192 
220  template <int dim>
221  struct TransformR2UInitialGuess
222  {
223  static const double KA[GeometryInfo<dim>::vertices_per_cell][dim];
224  static const double Kb[GeometryInfo<dim>::vertices_per_cell];
225  };
226 
227 
228  /*
229  Octave code:
230  M=[0 1; 1 1];
231  K1 = transpose(M) * inverse (M*transpose(M));
232  printf ("{%f, %f},\n", K1' );
233  */
234  template <>
235  const double
236  TransformR2UInitialGuess<1>::
238  {
239  {-1.000000},
240  {1.000000}
241  };
242 
243  template <>
244  const double
245  TransformR2UInitialGuess<1>::
246  Kb[GeometryInfo<1>::vertices_per_cell] = {1.000000, 0.000000};
247 
248 
249  /*
250  Octave code:
251  M=[0 1 0 1;0 0 1 1;1 1 1 1];
252  K2 = transpose(M) * inverse (M*transpose(M));
253  printf ("{%f, %f, %f},\n", K2' );
254  */
255  template <>
256  const double
257  TransformR2UInitialGuess<2>::
259  {
260  {-0.500000, -0.500000},
261  { 0.500000, -0.500000},
262  {-0.500000, 0.500000},
263  { 0.500000, 0.500000}
264  };
265 
266  /*
267  Octave code:
268  M=[0 1 0 1 0 1 0 1;0 0 1 1 0 0 1 1; 0 0 0 0 1 1 1 1; 1 1 1 1 1 1 1 1];
269  K3 = transpose(M) * inverse (M*transpose(M))
270  printf ("{%f, %f, %f, %f},\n", K3' );
271  */
272  template <>
273  const double
274  TransformR2UInitialGuess<2>::
276  {0.750000,0.250000,0.250000,-0.250000 };
277 
278 
279  template <>
280  const double
281  TransformR2UInitialGuess<3>::
283  {
284  {-0.250000, -0.250000, -0.250000},
285  { 0.250000, -0.250000, -0.250000},
286  {-0.250000, 0.250000, -0.250000},
287  { 0.250000, 0.250000, -0.250000},
288  {-0.250000, -0.250000, 0.250000},
289  { 0.250000, -0.250000, 0.250000},
290  {-0.250000, 0.250000, 0.250000},
291  { 0.250000, 0.250000, 0.250000}
292 
293  };
294 
295 
296  template <>
297  const double
298  TransformR2UInitialGuess<3>::
300  {0.500000,0.250000,0.250000,0.000000,0.250000,0.000000,0.000000,-0.250000};
301 
302  template<int dim, int spacedim>
303  Point<dim>
304  transform_real_to_unit_cell_initial_guess (const std::vector<Point<spacedim> > &vertex,
305  const Point<spacedim> &p)
306  {
307  Point<dim> p_unit;
308 
311 
312  KA.fill( (double *)(TransformR2UInitialGuess<dim>::KA) );
313  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
314  Kb(i) = TransformR2UInitialGuess<dim>::Kb[i];
315 
317  for (unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; v++)
318  for (unsigned int i=0; i<spacedim; ++i)
319  Y(i,v) = vertex[v][i];
320 
321  FullMatrix<double> A(spacedim,dim);
322  Y.mmult(A,KA); // A = Y*KA
323  ::Vector<double> b(spacedim);
324  Y.vmult(b,Kb); // b = Y*Kb
325 
326  for (unsigned int i=0; i<spacedim; ++i)
327  b(i) -= p[i];
328  b*=-1;
329 
330  ::Vector<double> dest(dim);
331 
332  FullMatrix<double> A_1(dim,spacedim);
333  if (dim<spacedim)
334  A_1.left_invert(A);
335  else
336  A_1.invert(A);
337 
338  A_1.vmult(dest,b); //A^{-1}*b
339 
340  for (unsigned int i=0; i<dim; ++i)
341  p_unit[i]=dest(i);
342 
343  return p_unit;
344  }
345 
346 
347  template <int dim, int spacedim>
348  void compute_shape_function_values_general (const unsigned int n_shape_functions,
349  const std::vector<Point<dim> > &unit_points,
350  typename ::MappingQGeneric<dim,spacedim>::InternalData &data)
351  {
352  const unsigned int n_points=unit_points.size();
353 
354  // Construct the tensor product polynomials used as shape functions for the
355  // Qp mapping of cells at the boundary.
357  tensor_pols (Polynomials::generate_complete_Lagrange_basis(data.line_support_points.get_points()));
358  Assert (n_shape_functions==tensor_pols.n(),
359  ExcInternalError());
360 
361  // then also construct the mapping from lexicographic to the Qp shape function numbering
362  const std::vector<unsigned int>
363  renumber (FETools::
364  lexicographic_to_hierarchic_numbering (
365  FiniteElementData<dim> (get_dpo_vector<dim>(data.polynomial_degree), 1,
366  data.polynomial_degree)));
367 
368  std::vector<double> values;
369  std::vector<Tensor<1,dim> > grads;
370  if (data.shape_values.size()!=0)
371  {
372  Assert(data.shape_values.size()==n_shape_functions*n_points,
373  ExcInternalError());
374  values.resize(n_shape_functions);
375  }
376  if (data.shape_derivatives.size()!=0)
377  {
378  Assert(data.shape_derivatives.size()==n_shape_functions*n_points,
379  ExcInternalError());
380  grads.resize(n_shape_functions);
381  }
382 
383  std::vector<Tensor<2,dim> > grad2;
384  if (data.shape_second_derivatives.size()!=0)
385  {
386  Assert(data.shape_second_derivatives.size()==n_shape_functions*n_points,
387  ExcInternalError());
388  grad2.resize(n_shape_functions);
389  }
390 
391  std::vector<Tensor<3,dim> > grad3;
392  if (data.shape_third_derivatives.size()!=0)
393  {
394  Assert(data.shape_third_derivatives.size()==n_shape_functions*n_points,
395  ExcInternalError());
396  grad3.resize(n_shape_functions);
397  }
398 
399  std::vector<Tensor<4,dim> > grad4;
400  if (data.shape_fourth_derivatives.size()!=0)
401  {
402  Assert(data.shape_fourth_derivatives.size()==n_shape_functions*n_points,
403  ExcInternalError());
404  grad4.resize(n_shape_functions);
405  }
406 
407 
408  if (data.shape_values.size()!=0 ||
409  data.shape_derivatives.size()!=0 ||
410  data.shape_second_derivatives.size()!=0 ||
411  data.shape_third_derivatives.size()!=0 ||
412  data.shape_fourth_derivatives.size()!=0 )
413  for (unsigned int point=0; point<n_points; ++point)
414  {
415  tensor_pols.compute(unit_points[point], values, grads, grad2, grad3, grad4);
416 
417  if (data.shape_values.size()!=0)
418  for (unsigned int i=0; i<n_shape_functions; ++i)
419  data.shape(point,renumber[i]) = values[i];
420 
421  if (data.shape_derivatives.size()!=0)
422  for (unsigned int i=0; i<n_shape_functions; ++i)
423  data.derivative(point,renumber[i]) = grads[i];
424 
425  if (data.shape_second_derivatives.size()!=0)
426  for (unsigned int i=0; i<n_shape_functions; ++i)
427  data.second_derivative(point,renumber[i]) = grad2[i];
428 
429  if (data.shape_third_derivatives.size()!=0)
430  for (unsigned int i=0; i<n_shape_functions; ++i)
431  data.third_derivative(point,renumber[i]) = grad3[i];
432 
433  if (data.shape_fourth_derivatives.size()!=0)
434  for (unsigned int i=0; i<n_shape_functions; ++i)
435  data.fourth_derivative(point,renumber[i]) = grad4[i];
436  }
437  }
438 
439 
440  void
441  compute_shape_function_values_hardcode (const unsigned int n_shape_functions,
442  const std::vector<Point<1> > &unit_points,
444  {
445  (void)n_shape_functions;
446  const unsigned int n_points=unit_points.size();
447  for (unsigned int k = 0 ; k < n_points ; ++k)
448  {
449  double x = unit_points[k](0);
450 
451  if (data.shape_values.size()!=0)
452  {
453  Assert(data.shape_values.size()==n_shape_functions*n_points,
454  ExcInternalError());
455  data.shape(k,0) = 1.-x;
456  data.shape(k,1) = x;
457  }
458  if (data.shape_derivatives.size()!=0)
459  {
460  Assert(data.shape_derivatives.size()==n_shape_functions*n_points,
461  ExcInternalError());
462  data.derivative(k,0)[0] = -1.;
463  data.derivative(k,1)[0] = 1.;
464  }
465  if (data.shape_second_derivatives.size()!=0)
466  {
467  Assert(data.shape_second_derivatives.size()==n_shape_functions*n_points,
468  ExcInternalError());
469  data.second_derivative(k,0)[0][0] = 0;
470  data.second_derivative(k,1)[0][0] = 0;
471  }
472  if (data.shape_third_derivatives.size()!=0)
473  {
474  Assert(data.shape_third_derivatives.size()==n_shape_functions*n_points,
475  ExcInternalError());
476 
477  Tensor<3,1> zero;
478  data.third_derivative(k,0) = zero;
479  data.third_derivative(k,1) = zero;
480  }
481  if (data.shape_fourth_derivatives.size()!=0)
482  {
483  Assert(data.shape_fourth_derivatives.size()==n_shape_functions*n_points,
484  ExcInternalError());
485 
486  Tensor<4,1> zero;
487  data.fourth_derivative(k,0) = zero;
488  data.fourth_derivative(k,1) = zero;
489  }
490  }
491  }
492 
493 
494  void
495  compute_shape_function_values_hardcode (const unsigned int n_shape_functions,
496  const std::vector<Point<2> > &unit_points,
498  {
499 
500  (void)n_shape_functions;
501  const unsigned int n_points=unit_points.size();
502  for (unsigned int k = 0 ; k < n_points ; ++k)
503  {
504  double x = unit_points[k](0);
505  double y = unit_points[k](1);
506 
507  if (data.shape_values.size()!=0)
508  {
509  Assert(data.shape_values.size()==n_shape_functions*n_points,
510  ExcInternalError());
511  data.shape(k,0) = (1.-x)*(1.-y);
512  data.shape(k,1) = x*(1.-y);
513  data.shape(k,2) = (1.-x)*y;
514  data.shape(k,3) = x*y;
515  }
516  if (data.shape_derivatives.size()!=0)
517  {
518  Assert(data.shape_derivatives.size()==n_shape_functions*n_points,
519  ExcInternalError());
520  data.derivative(k,0)[0] = (y-1.);
521  data.derivative(k,1)[0] = (1.-y);
522  data.derivative(k,2)[0] = -y;
523  data.derivative(k,3)[0] = y;
524  data.derivative(k,0)[1] = (x-1.);
525  data.derivative(k,1)[1] = -x;
526  data.derivative(k,2)[1] = (1.-x);
527  data.derivative(k,3)[1] = x;
528  }
529  if (data.shape_second_derivatives.size()!=0)
530  {
531  Assert(data.shape_second_derivatives.size()==n_shape_functions*n_points,
532  ExcInternalError());
533  data.second_derivative(k,0)[0][0] = 0;
534  data.second_derivative(k,1)[0][0] = 0;
535  data.second_derivative(k,2)[0][0] = 0;
536  data.second_derivative(k,3)[0][0] = 0;
537  data.second_derivative(k,0)[0][1] = 1.;
538  data.second_derivative(k,1)[0][1] = -1.;
539  data.second_derivative(k,2)[0][1] = -1.;
540  data.second_derivative(k,3)[0][1] = 1.;
541  data.second_derivative(k,0)[1][0] = 1.;
542  data.second_derivative(k,1)[1][0] = -1.;
543  data.second_derivative(k,2)[1][0] = -1.;
544  data.second_derivative(k,3)[1][0] = 1.;
545  data.second_derivative(k,0)[1][1] = 0;
546  data.second_derivative(k,1)[1][1] = 0;
547  data.second_derivative(k,2)[1][1] = 0;
548  data.second_derivative(k,3)[1][1] = 0;
549  }
550  if (data.shape_third_derivatives.size()!=0)
551  {
552  Assert(data.shape_third_derivatives.size()==n_shape_functions*n_points,
553  ExcInternalError());
554 
555  Tensor<3,2> zero;
556  for (unsigned int i=0; i<4; ++i)
557  data.third_derivative(k,i) = zero;
558  }
559  if (data.shape_fourth_derivatives.size()!=0)
560  {
561  Assert(data.shape_fourth_derivatives.size()==n_shape_functions*n_points,
562  ExcInternalError());
563  Tensor<4,2> zero;
564  for (unsigned int i=0; i<4; ++i)
565  data.fourth_derivative(k,i) = zero;
566  }
567  }
568  }
569 
570 
571 
572  void
573  compute_shape_function_values_hardcode (const unsigned int n_shape_functions,
574  const std::vector<Point<3> > &unit_points,
576  {
577  (void)n_shape_functions;
578  const unsigned int n_points=unit_points.size();
579  for (unsigned int k = 0 ; k < n_points ; ++k)
580  {
581  double x = unit_points[k](0);
582  double y = unit_points[k](1);
583  double z = unit_points[k](2);
584 
585  if (data.shape_values.size()!=0)
586  {
587  Assert(data.shape_values.size()==n_shape_functions*n_points,
588  ExcInternalError());
589  data.shape(k,0) = (1.-x)*(1.-y)*(1.-z);
590  data.shape(k,1) = x*(1.-y)*(1.-z);
591  data.shape(k,2) = (1.-x)*y*(1.-z);
592  data.shape(k,3) = x*y*(1.-z);
593  data.shape(k,4) = (1.-x)*(1.-y)*z;
594  data.shape(k,5) = x*(1.-y)*z;
595  data.shape(k,6) = (1.-x)*y*z;
596  data.shape(k,7) = x*y*z;
597  }
598  if (data.shape_derivatives.size()!=0)
599  {
600  Assert(data.shape_derivatives.size()==n_shape_functions*n_points,
601  ExcInternalError());
602  data.derivative(k,0)[0] = (y-1.)*(1.-z);
603  data.derivative(k,1)[0] = (1.-y)*(1.-z);
604  data.derivative(k,2)[0] = -y*(1.-z);
605  data.derivative(k,3)[0] = y*(1.-z);
606  data.derivative(k,4)[0] = (y-1.)*z;
607  data.derivative(k,5)[0] = (1.-y)*z;
608  data.derivative(k,6)[0] = -y*z;
609  data.derivative(k,7)[0] = y*z;
610  data.derivative(k,0)[1] = (x-1.)*(1.-z);
611  data.derivative(k,1)[1] = -x*(1.-z);
612  data.derivative(k,2)[1] = (1.-x)*(1.-z);
613  data.derivative(k,3)[1] = x*(1.-z);
614  data.derivative(k,4)[1] = (x-1.)*z;
615  data.derivative(k,5)[1] = -x*z;
616  data.derivative(k,6)[1] = (1.-x)*z;
617  data.derivative(k,7)[1] = x*z;
618  data.derivative(k,0)[2] = (x-1)*(1.-y);
619  data.derivative(k,1)[2] = x*(y-1.);
620  data.derivative(k,2)[2] = (x-1.)*y;
621  data.derivative(k,3)[2] = -x*y;
622  data.derivative(k,4)[2] = (1.-x)*(1.-y);
623  data.derivative(k,5)[2] = x*(1.-y);
624  data.derivative(k,6)[2] = (1.-x)*y;
625  data.derivative(k,7)[2] = x*y;
626  }
627  if (data.shape_second_derivatives.size()!=0)
628  {
629  Assert(data.shape_second_derivatives.size()==n_shape_functions*n_points,
630  ExcInternalError());
631  data.second_derivative(k,0)[0][0] = 0;
632  data.second_derivative(k,1)[0][0] = 0;
633  data.second_derivative(k,2)[0][0] = 0;
634  data.second_derivative(k,3)[0][0] = 0;
635  data.second_derivative(k,4)[0][0] = 0;
636  data.second_derivative(k,5)[0][0] = 0;
637  data.second_derivative(k,6)[0][0] = 0;
638  data.second_derivative(k,7)[0][0] = 0;
639  data.second_derivative(k,0)[1][1] = 0;
640  data.second_derivative(k,1)[1][1] = 0;
641  data.second_derivative(k,2)[1][1] = 0;
642  data.second_derivative(k,3)[1][1] = 0;
643  data.second_derivative(k,4)[1][1] = 0;
644  data.second_derivative(k,5)[1][1] = 0;
645  data.second_derivative(k,6)[1][1] = 0;
646  data.second_derivative(k,7)[1][1] = 0;
647  data.second_derivative(k,0)[2][2] = 0;
648  data.second_derivative(k,1)[2][2] = 0;
649  data.second_derivative(k,2)[2][2] = 0;
650  data.second_derivative(k,3)[2][2] = 0;
651  data.second_derivative(k,4)[2][2] = 0;
652  data.second_derivative(k,5)[2][2] = 0;
653  data.second_derivative(k,6)[2][2] = 0;
654  data.second_derivative(k,7)[2][2] = 0;
655 
656  data.second_derivative(k,0)[0][1] = (1.-z);
657  data.second_derivative(k,1)[0][1] = -(1.-z);
658  data.second_derivative(k,2)[0][1] = -(1.-z);
659  data.second_derivative(k,3)[0][1] = (1.-z);
660  data.second_derivative(k,4)[0][1] = z;
661  data.second_derivative(k,5)[0][1] = -z;
662  data.second_derivative(k,6)[0][1] = -z;
663  data.second_derivative(k,7)[0][1] = z;
664  data.second_derivative(k,0)[1][0] = (1.-z);
665  data.second_derivative(k,1)[1][0] = -(1.-z);
666  data.second_derivative(k,2)[1][0] = -(1.-z);
667  data.second_derivative(k,3)[1][0] = (1.-z);
668  data.second_derivative(k,4)[1][0] = z;
669  data.second_derivative(k,5)[1][0] = -z;
670  data.second_derivative(k,6)[1][0] = -z;
671  data.second_derivative(k,7)[1][0] = z;
672 
673  data.second_derivative(k,0)[0][2] = (1.-y);
674  data.second_derivative(k,1)[0][2] = -(1.-y);
675  data.second_derivative(k,2)[0][2] = y;
676  data.second_derivative(k,3)[0][2] = -y;
677  data.second_derivative(k,4)[0][2] = -(1.-y);
678  data.second_derivative(k,5)[0][2] = (1.-y);
679  data.second_derivative(k,6)[0][2] = -y;
680  data.second_derivative(k,7)[0][2] = y;
681  data.second_derivative(k,0)[2][0] = (1.-y);
682  data.second_derivative(k,1)[2][0] = -(1.-y);
683  data.second_derivative(k,2)[2][0] = y;
684  data.second_derivative(k,3)[2][0] = -y;
685  data.second_derivative(k,4)[2][0] = -(1.-y);
686  data.second_derivative(k,5)[2][0] = (1.-y);
687  data.second_derivative(k,6)[2][0] = -y;
688  data.second_derivative(k,7)[2][0] = y;
689 
690  data.second_derivative(k,0)[1][2] = (1.-x);
691  data.second_derivative(k,1)[1][2] = x;
692  data.second_derivative(k,2)[1][2] = -(1.-x);
693  data.second_derivative(k,3)[1][2] = -x;
694  data.second_derivative(k,4)[1][2] = -(1.-x);
695  data.second_derivative(k,5)[1][2] = -x;
696  data.second_derivative(k,6)[1][2] = (1.-x);
697  data.second_derivative(k,7)[1][2] = x;
698  data.second_derivative(k,0)[2][1] = (1.-x);
699  data.second_derivative(k,1)[2][1] = x;
700  data.second_derivative(k,2)[2][1] = -(1.-x);
701  data.second_derivative(k,3)[2][1] = -x;
702  data.second_derivative(k,4)[2][1] = -(1.-x);
703  data.second_derivative(k,5)[2][1] = -x;
704  data.second_derivative(k,6)[2][1] = (1.-x);
705  data.second_derivative(k,7)[2][1] = x;
706  }
707  if (data.shape_third_derivatives.size()!=0)
708  {
709  Assert(data.shape_third_derivatives.size()==n_shape_functions*n_points,
710  ExcInternalError());
711 
712  for (unsigned int i=0; i<3; ++i)
713  for (unsigned int j=0; j<3; ++j)
714  for (unsigned int l=0; l<3; ++l)
715  if ((i==j)||(j==l)||(l==i))
716  {
717  for (unsigned int m=0; m<8; ++m)
718  data.third_derivative(k,m)[i][j][l] = 0;
719  }
720  else
721  {
722  data.third_derivative(k,0)[i][j][l] = -1.;
723  data.third_derivative(k,1)[i][j][l] = 1.;
724  data.third_derivative(k,2)[i][j][l] = 1.;
725  data.third_derivative(k,3)[i][j][l] = -1.;
726  data.third_derivative(k,4)[i][j][l] = 1.;
727  data.third_derivative(k,5)[i][j][l] = -1.;
728  data.third_derivative(k,6)[i][j][l] = -1.;
729  data.third_derivative(k,7)[i][j][l] = 1.;
730  }
731 
732  }
733  if (data.shape_fourth_derivatives.size()!=0)
734  {
735  Assert(data.shape_fourth_derivatives.size()==n_shape_functions*n_points,
736  ExcInternalError());
737  Tensor<4,3> zero;
738  for (unsigned int i=0; i<8; ++i)
739  data.fourth_derivative(k,i) = zero;
740  }
741  }
742  }
743  }
744  }
745 }
746 
747 
748 
749 template<int dim, int spacedim>
750 MappingQGeneric<dim,spacedim>::InternalData::InternalData (const unsigned int polynomial_degree)
751  :
752  polynomial_degree (polynomial_degree),
753  n_shape_functions (Utilities::fixed_power<dim>(polynomial_degree+1)),
754  line_support_points(QGaussLobatto<1>(polynomial_degree+1))
755 {}
756 
757 
758 
759 template<int dim, int spacedim>
760 std::size_t
762 {
765  MemoryConsumption::memory_consumption (shape_derivatives) +
767  MemoryConsumption::memory_consumption (contravariant) +
768  MemoryConsumption::memory_consumption (unit_tangentials) +
770  MemoryConsumption::memory_consumption (mapping_support_points) +
771  MemoryConsumption::memory_consumption (cell_of_current_support_points) +
772  MemoryConsumption::memory_consumption (volume_elements) +
774  MemoryConsumption::memory_consumption (n_shape_functions));
775 }
776 
777 
778 template <int dim, int spacedim>
779 void
781 initialize (const UpdateFlags update_flags,
782  const Quadrature<dim> &q,
783  const unsigned int n_original_q_points)
784 {
785  // store the flags in the internal data object so we can access them
786  // in fill_fe_*_values()
787  this->update_each = update_flags;
788 
789  const unsigned int n_q_points = q.size();
790 
791  // see if we need the (transformation) shape function values
792  // and/or gradients and resize the necessary arrays
793  if (this->update_each & update_quadrature_points)
794  shape_values.resize(n_shape_functions * n_q_points);
795 
796  if (this->update_each & (update_covariant_transformation
809  shape_derivatives.resize(n_shape_functions * n_q_points);
810 
811  if (this->update_each & update_covariant_transformation)
812  covariant.resize(n_original_q_points);
813 
814  if (this->update_each & update_contravariant_transformation)
815  contravariant.resize(n_original_q_points);
816 
817  if (this->update_each & update_volume_elements)
818  volume_elements.resize(n_original_q_points);
819 
820  if (this->update_each &
822  shape_second_derivatives.resize(n_shape_functions * n_q_points);
823 
824  if (this->update_each &
826  shape_third_derivatives.resize(n_shape_functions * n_q_points);
827 
828  if (this->update_each &
830  shape_fourth_derivatives.resize(n_shape_functions * n_q_points);
831 
832  // now also fill the various fields with their correct values
833  compute_shape_function_values (q.get_points());
834 }
835 
836 
837 
838 template <int dim, int spacedim>
839 void
841 initialize_face (const UpdateFlags update_flags,
842  const Quadrature<dim> &q,
843  const unsigned int n_original_q_points)
844 {
845  initialize (update_flags, q, n_original_q_points);
846 
847  if (dim > 1)
848  {
849  if (this->update_each & (update_boundary_forms |
854  {
855  aux.resize (dim-1, std::vector<Tensor<1,spacedim> > (n_original_q_points));
856 
857  // Compute tangentials to the unit cell.
858  for (unsigned int i=0; i<unit_tangentials.size(); ++i)
859  unit_tangentials[i].resize (n_original_q_points);
860  switch (dim)
861  {
862  case 2:
863  {
864  // ensure a counterclockwise orientation of tangentials
865  static const int tangential_orientation[4]= {-1,1,1,-1};
866  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
867  {
868  Tensor<1,dim> tang;
869  tang[1-i/2] = tangential_orientation[i];
870  std::fill (unit_tangentials[i].begin(),
871  unit_tangentials[i].end(),
872  tang);
873  }
874 
875  break;
876  }
877 
878  case 3:
879  {
880  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
881  {
882  Tensor<1,dim> tang1, tang2;
883 
884  const unsigned int nd=
886 
887  // first tangential
888  // vector in direction
889  // of the (nd+1)%3 axis
890  // and inverted in case
891  // of unit inward normal
892  tang1[(nd+1)%dim]=GeometryInfo<dim>::unit_normal_orientation[i];
893  // second tangential
894  // vector in direction
895  // of the (nd+2)%3 axis
896  tang2[(nd+2)%dim]=1.;
897 
898  // same unit tangents
899  // for all quadrature
900  // points on this face
901  std::fill (unit_tangentials[i].begin(),
902  unit_tangentials[i].end(),
903  tang1);
904  std::fill (unit_tangentials[GeometryInfo<dim>::faces_per_cell+i].begin(),
905  unit_tangentials[GeometryInfo<dim>::faces_per_cell+i].end(),
906  tang2);
907  }
908 
909  break;
910  }
911 
912  default:
913  Assert (false, ExcNotImplemented());
914  }
915  }
916  }
917 }
918 
919 
920 
921 
922 
923 template<>
924 void
926 compute_shape_function_values (const std::vector<Point<1> > &unit_points)
927 {
928  // if the polynomial degree is one, then we can simplify code a bit
929  // by using hard-coded shape functions.
930  if (polynomial_degree == 1)
931  internal::MappingQ1::compute_shape_function_values_hardcode (n_shape_functions,
932  unit_points, *this);
933  else
934  {
935  // otherwise ask an object that describes the polynomial space
936  internal::MappingQ1::compute_shape_function_values_general<1,1>(n_shape_functions,
937  unit_points,*this);
938  }
939 }
940 
941 template<>
942 void
944 compute_shape_function_values (const std::vector<Point<2> > &unit_points)
945 {
946  // if the polynomial degree is one, then we can simplify code a bit
947  // by using hard-coded shape functions.
948  if (polynomial_degree == 1)
949  internal::MappingQ1::compute_shape_function_values_hardcode (n_shape_functions,
950  unit_points, *this);
951  else
952  {
953  // otherwise ask an object that describes the polynomial space
954  internal::MappingQ1::compute_shape_function_values_general<2,2>(n_shape_functions,
955  unit_points,*this);
956  }
957 }
958 
959 template<>
960 void
962 compute_shape_function_values (const std::vector<Point<3> > &unit_points)
963 {
964  // if the polynomial degree is one, then we can simplify code a bit
965  // by using hard-coded shape functions.
966  if (polynomial_degree == 1)
967  internal::MappingQ1::compute_shape_function_values_hardcode (n_shape_functions,
968  unit_points, *this);
969  else
970  {
971  // otherwise ask an object that describes the polynomial space
972  internal::MappingQ1::compute_shape_function_values_general<3,3>(n_shape_functions,
973  unit_points,*this);
974  }
975 }
976 
977 template<int dim, int spacedim>
978 void
980 compute_shape_function_values (const std::vector<Point<dim> > &unit_points)
981 {
982  // for non-matching combinations of dim and spacedim, just run the general
983  // case
984  internal::MappingQ1::compute_shape_function_values_general<dim,spacedim>(n_shape_functions,
985  unit_points,*this);
986 }
987 
988 
989 namespace
990 {
1000  template<int dim>
1002  compute_laplace_vector(const unsigned int polynomial_degree)
1003  {
1004  Table<2,double> lvs;
1005 
1006  Assert(lvs.n_rows()==0, ExcInternalError());
1007  Assert(dim==2 || dim==3, ExcNotImplemented());
1008 
1009  // for degree==1, we shouldn't have to compute any support points, since all
1010  // of them are on the vertices
1012 
1013  const unsigned int n_inner = Utilities::fixed_power<dim>(polynomial_degree-1);
1014  const unsigned int n_outer = (dim==1) ? 2 :
1015  ((dim==2) ?
1016  4+4*(polynomial_degree-1) :
1018 
1019 
1020  // compute the shape gradients at the quadrature points on the unit cell
1021  const QGauss<dim> quadrature(polynomial_degree+1);
1022  const unsigned int n_q_points=quadrature.size();
1023 
1024  typename MappingQGeneric<dim>::InternalData quadrature_data(polynomial_degree);
1025  quadrature_data.shape_derivatives.resize(quadrature_data.n_shape_functions *
1026  n_q_points);
1027  quadrature_data.compute_shape_function_values(quadrature.get_points());
1028 
1029  // Compute the stiffness matrix of the inner dofs
1030  FullMatrix<long double> S(n_inner);
1031  for (unsigned int point=0; point<n_q_points; ++point)
1032  for (unsigned int i=0; i<n_inner; ++i)
1033  for (unsigned int j=0; j<n_inner; ++j)
1034  {
1035  long double res = 0.;
1036  for (unsigned int l=0; l<dim; ++l)
1037  res += (long double)quadrature_data.derivative(point, n_outer+i)[l] *
1038  (long double)quadrature_data.derivative(point, n_outer+j)[l];
1039 
1040  S(i,j) += res * (long double)quadrature.weight(point);
1041  }
1042 
1043  // Compute the components of T to be the product of gradients of inner and
1044  // outer shape functions.
1045  FullMatrix<long double> T(n_inner, n_outer);
1046  for (unsigned int point=0; point<n_q_points; ++point)
1047  for (unsigned int i=0; i<n_inner; ++i)
1048  for (unsigned int k=0; k<n_outer; ++k)
1049  {
1050  long double res = 0.;
1051  for (unsigned int l=0; l<dim; ++l)
1052  res += (long double)quadrature_data.derivative(point, n_outer+i)[l] *
1053  (long double)quadrature_data.derivative(point, k)[l];
1054 
1055  T(i,k) += res *(long double)quadrature.weight(point);
1056  }
1057 
1058  FullMatrix<long double> S_1(n_inner);
1059  S_1.invert(S);
1060 
1061  FullMatrix<long double> S_1_T(n_inner, n_outer);
1062 
1063  // S:=S_1*T
1064  S_1.mmult(S_1_T,T);
1065 
1066  // Resize and initialize the lvs
1067  lvs.reinit (n_inner, n_outer);
1068  for (unsigned int i=0; i<n_inner; ++i)
1069  for (unsigned int k=0; k<n_outer; ++k)
1070  lvs(i,k) = -S_1_T(i,k);
1071 
1072  return lvs;
1073  }
1074 
1075 
1088  compute_support_point_weights_on_quad(const unsigned int polynomial_degree)
1089  {
1090  Table<2,double> loqvs;
1091 
1092  // we are asked to compute weights for interior support points, but
1093  // there are no interior points if degree==1
1094  if (polynomial_degree == 1)
1095  return loqvs;
1096 
1097  const unsigned int n_inner_2d=(polynomial_degree-1)*(polynomial_degree-1);
1098  const unsigned int n_outer_2d=4+4*(polynomial_degree-1);
1099 
1100  // first check whether we have precomputed the values for some polynomial
1101  // degree; the sizes of arrays is n_inner_2d*n_outer_2d
1102  if (polynomial_degree == 2)
1103  {
1104  // (checked these values against the output of compute_laplace_vector
1105  // again, and found they're indeed right -- just in case someone wonders
1106  // where they come from -- WB)
1107  static const double loqv2[1*8]
1108  = {1/16., 1/16., 1/16., 1/16., 3/16., 3/16., 3/16., 3/16.};
1109  Assert (sizeof(loqv2)/sizeof(loqv2[0]) ==
1110  n_inner_2d * n_outer_2d,
1111  ExcInternalError());
1112 
1113  // copy and return
1114  loqvs.reinit(n_inner_2d, n_outer_2d);
1115  for (unsigned int unit_point=0; unit_point<n_inner_2d; ++unit_point)
1116  for (unsigned int k=0; k<n_outer_2d; ++k)
1117  loqvs[unit_point][k] = loqv2[unit_point*n_outer_2d+k];
1118  }
1119  else
1120  {
1121  // not precomputed, then do so now
1122  loqvs = compute_laplace_vector<2>(polynomial_degree);
1123  }
1124 
1125  // the sum of weights of the points at the outer rim should be one. check
1126  // this
1127  for (unsigned int unit_point=0; unit_point<loqvs.n_rows(); ++unit_point)
1128  Assert(std::fabs(std::accumulate(loqvs[unit_point].begin(),
1129  loqvs[unit_point].end(),0.)-1)<1e-13*polynomial_degree,
1130  ExcInternalError());
1131 
1132  return loqvs;
1133  }
1134 
1135 
1136 
1147  compute_support_point_weights_on_hex(const unsigned int polynomial_degree)
1148  {
1149  Table<2,double> lohvs;
1150 
1151  // we are asked to compute weights for interior support points, but
1152  // there are no interior points if degree==1
1153  if (polynomial_degree == 1)
1154  return lohvs;
1155 
1156  const unsigned int n_inner = Utilities::fixed_power<3>(polynomial_degree-1);
1157  const unsigned int n_outer = 8+12*(polynomial_degree-1)+6*(polynomial_degree-1)*(polynomial_degree-1);
1158 
1159  // first check whether we have precomputed the values for some polynomial
1160  // degree; the sizes of arrays is n_inner_2d*n_outer_2d
1161  if (polynomial_degree == 2)
1162  {
1163  static const double lohv2[26]
1164  = {1/128., 1/128., 1/128., 1/128., 1/128., 1/128., 1/128., 1/128.,
1165  7/192., 7/192., 7/192., 7/192., 7/192., 7/192., 7/192., 7/192.,
1166  7/192., 7/192., 7/192., 7/192.,
1167  1/12., 1/12., 1/12., 1/12., 1/12., 1/12.
1168  };
1169 
1170  // copy and return
1171  lohvs.reinit(n_inner, n_outer);
1172  for (unsigned int unit_point=0; unit_point<n_inner; ++unit_point)
1173  for (unsigned int k=0; k<n_outer; ++k)
1174  lohvs[unit_point][k] = lohv2[unit_point*n_outer+k];
1175  }
1176  else
1177  {
1178  // not precomputed, then do so now
1179  lohvs = compute_laplace_vector<3>(polynomial_degree);
1180  }
1181 
1182  // the sum of weights of the points at the outer rim should be one. check
1183  // this
1184  for (unsigned int unit_point=0; unit_point<n_inner; ++unit_point)
1185  Assert(std::fabs(std::accumulate(lohvs[unit_point].begin(),
1186  lohvs[unit_point].end(),0.) - 1)<1e-13*polynomial_degree,
1187  ExcInternalError());
1188 
1189  return lohvs;
1190  }
1191 
1196  std::vector<Table<2,double> >
1197  compute_support_point_weights_perimeter_to_interior(const unsigned int polynomial_degree,
1198  const unsigned int dim)
1199  {
1200  Assert(dim > 0 && dim <= 3, ExcImpossibleInDim(dim));
1201  std::vector<Table<2,double> > output(dim);
1202  if (polynomial_degree <= 1)
1203  return output;
1204 
1205  // fill the 1D interior weights
1206  QGaussLobatto<1> quadrature(polynomial_degree+1);
1208  for (unsigned int q=0; q<polynomial_degree-1; ++q)
1209  for (unsigned int i=0; i<GeometryInfo<1>::vertices_per_cell; ++i)
1210  output[0](q,i) = GeometryInfo<1>::d_linear_shape_function(quadrature.point(q+1),
1211  i);
1212 
1213  if (dim > 1)
1214  output[1] = compute_support_point_weights_on_quad(polynomial_degree);
1215 
1216  if (dim > 2)
1217  output[2] = compute_support_point_weights_on_hex(polynomial_degree);
1218 
1219  return output;
1220  }
1221 
1225  template <int dim>
1227  compute_support_point_weights_cell(const unsigned int polynomial_degree)
1228  {
1229  Assert(dim > 0 && dim <= 3, ExcImpossibleInDim(dim));
1230  if (polynomial_degree <= 1)
1231  return Table<2,double>();
1232 
1233  QGaussLobatto<dim> quadrature(polynomial_degree+1);
1234  std::vector<unsigned int> h2l(quadrature.size());
1235  FETools::hierarchic_to_lexicographic_numbering<dim>(polynomial_degree, h2l);
1236 
1237  Table<2,double> output(quadrature.size() - GeometryInfo<dim>::vertices_per_cell,
1239  for (unsigned int q=0; q<output.size(0); ++q)
1240  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
1242  i);
1243 
1244  return output;
1245  }
1246 }
1247 
1248 
1249 
1250 
1251 template<int dim, int spacedim>
1253  :
1254  polynomial_degree(p),
1255  line_support_points(this->polynomial_degree+1),
1256  fe_q(dim == 3 ? new FE_Q<dim>(this->polynomial_degree) : 0),
1257  support_point_weights_perimeter_to_interior (compute_support_point_weights_perimeter_to_interior(this->polynomial_degree, dim)),
1258  support_point_weights_cell (compute_support_point_weights_cell<dim>(this->polynomial_degree))
1259 {
1260  Assert (p >= 1, ExcMessage ("It only makes sense to create polynomial mappings "
1261  "with a polynomial degree greater or equal to one."));
1262 }
1263 
1264 
1265 
1266 template<int dim, int spacedim>
1268  :
1269  polynomial_degree(mapping.polynomial_degree),
1270  line_support_points(mapping.line_support_points),
1271  fe_q(dim == 3 ? new FE_Q<dim>(*mapping.fe_q) : 0),
1272  support_point_weights_perimeter_to_interior (mapping.support_point_weights_perimeter_to_interior),
1273  support_point_weights_cell (mapping.support_point_weights_cell)
1274 {}
1275 
1276 
1277 
1278 
1279 template<int dim, int spacedim>
1282 {
1283  return new MappingQGeneric<dim,spacedim>(*this);
1284 }
1285 
1286 
1287 
1288 
1289 template<int dim, int spacedim>
1290 unsigned int
1292 {
1293  return polynomial_degree;
1294 }
1295 
1296 
1297 
1298 template<int dim, int spacedim>
1302  const Point<dim> &p) const
1303 {
1304  // set up the polynomial space
1306  tensor_pols (Polynomials::generate_complete_Lagrange_basis(line_support_points.get_points()));
1307  Assert (tensor_pols.n() == Utilities::fixed_power<dim>(polynomial_degree+1),
1308  ExcInternalError());
1309 
1310  // then also construct the mapping from lexicographic to the Qp shape function numbering
1311  const std::vector<unsigned int>
1312  renumber (FETools::
1313  lexicographic_to_hierarchic_numbering (
1314  FiniteElementData<dim> (get_dpo_vector<dim>(polynomial_degree), 1,
1315  polynomial_degree)));
1316 
1317  const std::vector<Point<spacedim> > support_points
1318  = this->compute_mapping_support_points(cell);
1319 
1320  Point<spacedim> mapped_point;
1321  for (unsigned int i=0; i<tensor_pols.n(); ++i)
1322  mapped_point += support_points[renumber[i]] * tensor_pols.compute_value (i, p);
1323 
1324  return mapped_point;
1325 }
1326 
1327 
1328 // In the code below, GCC tries to instantiate MappingQGeneric<3,4> when
1329 // seeing which of the overloaded versions of
1330 // do_transform_real_to_unit_cell_internal() to call. This leads to bad
1331 // error messages and, generally, nothing very good. Avoid this by ensuring
1332 // that this class exists, but does not have an inner InternalData
1333 // type, thereby ruling out the codim-1 version of the function
1334 // below when doing overload resolution.
1335 template <>
1336 class MappingQGeneric<3,4>
1337 {};
1338 
1339 namespace
1340 {
1348  template<int dim, int spacedim>
1350  compute_mapped_location_of_point (const typename MappingQGeneric<dim,spacedim>::InternalData &data)
1351  {
1352  AssertDimension (data.shape_values.size(),
1353  data.mapping_support_points.size());
1354 
1355  // use now the InternalData to compute the point in real space.
1356  Point<spacedim> p_real;
1357  for (unsigned int i=0; i<data.mapping_support_points.size(); ++i)
1358  p_real += data.mapping_support_points[i] * data.shape(0,i);
1359 
1360  return p_real;
1361  }
1362 
1363 
1367  template <int dim>
1368  Point<dim>
1369  do_transform_real_to_unit_cell_internal
1370  (const typename Triangulation<dim,dim>::cell_iterator &cell,
1371  const Point<dim> &p,
1372  const Point<dim> &initial_p_unit,
1374  {
1375  const unsigned int spacedim = dim;
1376 
1377  const unsigned int n_shapes=mdata.shape_values.size();
1378  (void)n_shapes;
1379  Assert(n_shapes!=0, ExcInternalError());
1380  AssertDimension (mdata.shape_derivatives.size(), n_shapes);
1381 
1382  std::vector<Point<spacedim> > &points=mdata.mapping_support_points;
1383  AssertDimension (points.size(), n_shapes);
1384 
1385 
1386  // Newton iteration to solve
1387  // f(x)=p(x)-p=0
1388  // where we are looking for 'x' and p(x) is the forward transformation
1389  // from unit to real cell. We solve this using a Newton iteration
1390  // x_{n+1}=x_n-[f'(x)]^{-1}f(x)
1391  // The start value is set to be the linear approximation to the cell
1392 
1393  // The shape values and derivatives of the mapping at this point are
1394  // previously computed.
1395 
1396  Point<dim> p_unit = initial_p_unit;
1397 
1398  mdata.compute_shape_function_values(std::vector<Point<dim> > (1, p_unit));
1399 
1400  Point<spacedim> p_real = compute_mapped_location_of_point<dim,spacedim>(mdata);
1401  Tensor<1,spacedim> f = p_real-p;
1402 
1403  // early out if we already have our point
1404  if (f.norm_square() < 1e-24 * cell->diameter() * cell->diameter())
1405  return p_unit;
1406 
1407  // we need to compare the position of the computed p(x) against the given
1408  // point 'p'. We will terminate the iteration and return 'x' if they are
1409  // less than eps apart. The question is how to choose eps -- or, put maybe
1410  // more generally: in which norm we want these 'p' and 'p(x)' to be eps
1411  // apart.
1412  //
1413  // the question is difficult since we may have to deal with very elongated
1414  // cells where we may achieve 1e-12*h for the distance of these two points
1415  // in the 'long' direction, but achieving this tolerance in the 'short'
1416  // direction of the cell may not be possible
1417  //
1418  // what we do instead is then to terminate iterations if
1419  // \| p(x) - p \|_A < eps
1420  // where the A-norm is somehow induced by the transformation of the cell.
1421  // in particular, we want to measure distances relative to the sizes of
1422  // the cell in its principal directions.
1423  //
1424  // to define what exactly A should be, note that to first order we have
1425  // the following (assuming that x* is the solution of the problem, i.e.,
1426  // p(x*)=p):
1427  // p(x) - p = p(x) - p(x*)
1428  // = -grad p(x) * (x*-x) + higher order terms
1429  // This suggest to measure with a norm that corresponds to
1430  // A = {[grad p(x]^T [grad p(x)]}^{-1}
1431  // because then
1432  // \| p(x) - p \|_A \approx \| x - x* \|
1433  // Consequently, we will try to enforce that
1434  // \| p(x) - p \|_A = \| f \| <= eps
1435  //
1436  // Note that using this norm is a bit dangerous since the norm changes
1437  // in every iteration (A isn't fixed by depends on xk). However, if the
1438  // cell is not too deformed (it may be stretched, but not twisted) then
1439  // the mapping is almost linear and A is indeed constant or nearly so.
1440  const double eps = 1.e-11;
1441  const unsigned int newton_iteration_limit = 20;
1442 
1443  unsigned int newton_iteration = 0;
1444  double last_f_weighted_norm;
1445  do
1446  {
1447 #ifdef DEBUG_TRANSFORM_REAL_TO_UNIT_CELL
1448  std::cout << "Newton iteration " << newton_iteration << std::endl;
1449 #endif
1450 
1451  // f'(x)
1452  Tensor<2,spacedim> df;
1453  for (unsigned int k=0; k<mdata.n_shape_functions; ++k)
1454  {
1455  const Tensor<1,dim> &grad_transform=mdata.derivative(0,k);
1456  const Point<spacedim> &point=points[k];
1457 
1458  for (unsigned int i=0; i<spacedim; ++i)
1459  for (unsigned int j=0; j<dim; ++j)
1460  df[i][j]+=point[i]*grad_transform[j];
1461  }
1462 
1463  // Solve [f'(x)]d=f(x)
1464  Tensor<2,spacedim> df_inverse = invert(df);
1465  const Tensor<1,spacedim> delta = df_inverse * static_cast<const Tensor<1,spacedim>&>(f);
1466 
1467 #ifdef DEBUG_TRANSFORM_REAL_TO_UNIT_CELL
1468  std::cout << " delta=" << delta << std::endl;
1469 #endif
1470 
1471  // do a line search
1472  double step_length = 1;
1473  do
1474  {
1475  // update of p_unit. The spacedim-th component of transformed point
1476  // is simply ignored in codimension one case. When this component is
1477  // not zero, then we are projecting the point to the surface or
1478  // curve identified by the cell.
1479  Point<dim> p_unit_trial = p_unit;
1480  for (unsigned int i=0; i<dim; ++i)
1481  p_unit_trial[i] -= step_length * delta[i];
1482 
1483  // shape values and derivatives
1484  // at new p_unit point
1485  mdata.compute_shape_function_values(std::vector<Point<dim> > (1, p_unit_trial));
1486 
1487  // f(x)
1488  Point<spacedim> p_real_trial = compute_mapped_location_of_point<dim,spacedim>(mdata);
1489  const Tensor<1,spacedim> f_trial = p_real_trial-p;
1490 
1491 #ifdef DEBUG_TRANSFORM_REAL_TO_UNIT_CELL
1492  std::cout << " step_length=" << step_length << std::endl
1493  << " ||f || =" << f.norm() << std::endl
1494  << " ||f*|| =" << f_trial.norm() << std::endl
1495  << " ||f*||_A =" << (df_inverse * f_trial).norm() << std::endl;
1496 #endif
1497 
1498  // see if we are making progress with the current step length
1499  // and if not, reduce it by a factor of two and try again
1500  //
1501  // strictly speaking, we should probably use the same norm as we use
1502  // for the outer algorithm. in practice, line search is just a
1503  // crutch to find a "reasonable" step length, and so using the l2
1504  // norm is probably just fine
1505  if (f_trial.norm() < f.norm())
1506  {
1507  p_real = p_real_trial;
1508  p_unit = p_unit_trial;
1509  f = f_trial;
1510  break;
1511  }
1512  else if (step_length > 0.05)
1513  step_length /= 2;
1514  else
1515  AssertThrow (false,
1517  }
1518  while (true);
1519 
1520  ++newton_iteration;
1521  if (newton_iteration > newton_iteration_limit)
1522  AssertThrow (false,
1524  last_f_weighted_norm = (df_inverse * f).norm();
1525  }
1526  while (last_f_weighted_norm > eps);
1527 
1528  return p_unit;
1529  }
1530 
1531 
1532 
1536  template <int dim>
1537  Point<dim>
1538  do_transform_real_to_unit_cell_internal_codim1
1539  (const typename Triangulation<dim,dim+1>::cell_iterator &cell,
1540  const Point<dim+1> &p,
1541  const Point<dim> &initial_p_unit,
1543  {
1544  const unsigned int spacedim = dim+1;
1545 
1546  const unsigned int n_shapes=mdata.shape_values.size();
1547  (void)n_shapes;
1548  Assert(n_shapes!=0, ExcInternalError());
1549  Assert(mdata.shape_derivatives.size()==n_shapes, ExcInternalError());
1550  Assert(mdata.shape_second_derivatives.size()==n_shapes, ExcInternalError());
1551 
1552  std::vector<Point<spacedim> > &points=mdata.mapping_support_points;
1553  Assert(points.size()==n_shapes, ExcInternalError());
1554 
1555  Point<spacedim> p_minus_F;
1556 
1557  Tensor<1,spacedim> DF[dim];
1558  Tensor<1,spacedim> D2F[dim][dim];
1559 
1560  Point<dim> p_unit = initial_p_unit;
1561  Point<dim> f;
1562  Tensor<2,dim> df;
1563 
1564  // Evaluate first and second derivatives
1565  mdata.compute_shape_function_values(std::vector<Point<dim> > (1, p_unit));
1566 
1567  for (unsigned int k=0; k<mdata.n_shape_functions; ++k)
1568  {
1569  const Tensor<1,dim> &grad_phi_k = mdata.derivative(0,k);
1570  const Tensor<2,dim> &hessian_k = mdata.second_derivative(0,k);
1571  const Point<spacedim> &point_k = points[k];
1572 
1573  for (unsigned int j=0; j<dim; ++j)
1574  {
1575  DF[j] += grad_phi_k[j] * point_k;
1576  for (unsigned int l=0; l<dim; ++l)
1577  D2F[j][l] += hessian_k[j][l] * point_k;
1578  }
1579  }
1580 
1581  p_minus_F = p;
1582  p_minus_F -= compute_mapped_location_of_point<dim,spacedim>(mdata);
1583 
1584 
1585  for (unsigned int j=0; j<dim; ++j)
1586  f[j] = DF[j] * p_minus_F;
1587 
1588  for (unsigned int j=0; j<dim; ++j)
1589  {
1590  f[j] = DF[j] * p_minus_F;
1591  for (unsigned int l=0; l<dim; ++l)
1592  df[j][l] = -DF[j]*DF[l] + D2F[j][l] * p_minus_F;
1593  }
1594 
1595 
1596  const double eps = 1.e-12*cell->diameter();
1597  const unsigned int loop_limit = 10;
1598 
1599  unsigned int loop=0;
1600 
1601  while (f.norm()>eps && loop++<loop_limit)
1602  {
1603  // Solve [df(x)]d=f(x)
1604  const Tensor<1,dim> d = invert(df) * static_cast<const Tensor<1,dim>&>(f);
1605  p_unit -= d;
1606 
1607  for (unsigned int j=0; j<dim; ++j)
1608  {
1609  DF[j].clear();
1610  for (unsigned int l=0; l<dim; ++l)
1611  D2F[j][l].clear();
1612  }
1613 
1614  mdata.compute_shape_function_values(std::vector<Point<dim> > (1, p_unit));
1615 
1616  for (unsigned int k=0; k<mdata.n_shape_functions; ++k)
1617  {
1618  const Tensor<1,dim> &grad_phi_k = mdata.derivative(0,k);
1619  const Tensor<2,dim> &hessian_k = mdata.second_derivative(0,k);
1620  const Point<spacedim> &point_k = points[k];
1621 
1622  for (unsigned int j=0; j<dim; ++j)
1623  {
1624  DF[j] += grad_phi_k[j] * point_k;
1625  for (unsigned int l=0; l<dim; ++l)
1626  D2F[j][l] += hessian_k[j][l] * point_k;
1627  }
1628  }
1629 
1630  //TODO: implement a line search here in much the same way as for
1631  // the corresponding function above that does so for dim==spacedim
1632  p_minus_F = p;
1633  p_minus_F -= compute_mapped_location_of_point<dim,spacedim>(mdata);
1634 
1635  for (unsigned int j=0; j<dim; ++j)
1636  {
1637  f[j] = DF[j] * p_minus_F;
1638  for (unsigned int l=0; l<dim; ++l)
1639  df[j][l] = -DF[j]*DF[l] + D2F[j][l] * p_minus_F;
1640  }
1641 
1642  }
1643 
1644 
1645  // Here we check that in the last execution of while the first
1646  // condition was already wrong, meaning the residual was below
1647  // eps. Only if the first condition failed, loop will have been
1648  // increased and tested, and thus have reached the limit.
1650 
1651  return p_unit;
1652  }
1653 
1654 
1655 }
1656 
1657 
1658 
1659 // visual studio freaks out when trying to determine if
1660 // do_transform_real_to_unit_cell_internal with dim=3 and spacedim=4 is a good
1661 // candidate. So instead of letting the compiler pick the correct overload, we
1662 // use template specialization to make sure we pick up the right function to
1663 // call:
1664 
1665 template<int dim, int spacedim>
1666 Point<dim>
1670  const Point<spacedim> &,
1671  const Point<dim> &) const
1672 {
1673  // default implementation (should never be called)
1674  Assert(false, ExcInternalError());
1675  return Point<dim>();
1676 }
1677 
1678 template<>
1679 Point<1>
1683  const Point<1> &p,
1684  const Point<1> &initial_p_unit) const
1685 {
1686  const int dim = 1;
1687  const int spacedim = 1;
1688 
1689  const Quadrature<dim> point_quadrature(initial_p_unit);
1690 
1692  if (spacedim>dim)
1693  update_flags |= update_jacobian_grads;
1694  std_cxx11::unique_ptr<InternalData> mdata (get_data(update_flags,
1695  point_quadrature));
1696 
1697  mdata->mapping_support_points = this->compute_mapping_support_points (cell);
1698 
1699  // dispatch to the various specializations for spacedim=dim,
1700  // spacedim=dim+1, etc
1701  return do_transform_real_to_unit_cell_internal<1>(cell, p, initial_p_unit, *mdata);
1702 }
1703 
1704 template<>
1705 Point<2>
1709  const Point<2> &p,
1710  const Point<2> &initial_p_unit) const
1711 {
1712  const int dim = 2;
1713  const int spacedim = 2;
1714 
1715  const Quadrature<dim> point_quadrature(initial_p_unit);
1716 
1718  if (spacedim>dim)
1719  update_flags |= update_jacobian_grads;
1720  std_cxx11::unique_ptr<InternalData> mdata (get_data(update_flags,
1721  point_quadrature));
1722 
1723  mdata->mapping_support_points = this->compute_mapping_support_points (cell);
1724 
1725  // dispatch to the various specializations for spacedim=dim,
1726  // spacedim=dim+1, etc
1727  return do_transform_real_to_unit_cell_internal<2>(cell, p, initial_p_unit, *mdata);
1728 }
1729 
1730 template<>
1731 Point<3>
1735  const Point<3> &p,
1736  const Point<3> &initial_p_unit) const
1737 {
1738  const int dim = 3;
1739  const int spacedim = 3;
1740 
1741  const Quadrature<dim> point_quadrature(initial_p_unit);
1742 
1744  if (spacedim>dim)
1745  update_flags |= update_jacobian_grads;
1746  std_cxx11::unique_ptr<InternalData> mdata (get_data(update_flags,
1747  point_quadrature));
1748 
1749  mdata->mapping_support_points = this->compute_mapping_support_points (cell);
1750 
1751  // dispatch to the various specializations for spacedim=dim,
1752  // spacedim=dim+1, etc
1753  return do_transform_real_to_unit_cell_internal<3>(cell, p, initial_p_unit, *mdata);
1754 }
1755 
1756 template<>
1757 Point<1>
1761  const Point<2> &p,
1762  const Point<1> &initial_p_unit) const
1763 {
1764  const int dim = 1;
1765  const int spacedim = 2;
1766 
1767  const Quadrature<dim> point_quadrature(initial_p_unit);
1768 
1770  if (spacedim>dim)
1771  update_flags |= update_jacobian_grads;
1772  std_cxx11::unique_ptr<InternalData> mdata (get_data(update_flags,
1773  point_quadrature));
1774 
1775  mdata->mapping_support_points = this->compute_mapping_support_points (cell);
1776 
1777  // dispatch to the various specializations for spacedim=dim,
1778  // spacedim=dim+1, etc
1779  return do_transform_real_to_unit_cell_internal_codim1<1>(cell, p, initial_p_unit, *mdata);
1780 }
1781 
1782 template<>
1783 Point<2>
1787  const Point<3> &p,
1788  const Point<2> &initial_p_unit) const
1789 {
1790  const int dim = 2;
1791  const int spacedim = 3;
1792 
1793  const Quadrature<dim> point_quadrature(initial_p_unit);
1794 
1796  if (spacedim>dim)
1797  update_flags |= update_jacobian_grads;
1798  std_cxx11::unique_ptr<InternalData> mdata (get_data(update_flags,
1799  point_quadrature));
1800 
1801  mdata->mapping_support_points = this->compute_mapping_support_points (cell);
1802 
1803  // dispatch to the various specializations for spacedim=dim,
1804  // spacedim=dim+1, etc
1805  return do_transform_real_to_unit_cell_internal_codim1<2>(cell, p, initial_p_unit, *mdata);
1806 }
1807 
1808 template<>
1809 Point<1>
1813  const Point<3> &,
1814  const Point<1> &) const
1815 {
1816  Assert (false, ExcNotImplemented());
1817  return Point<1>();
1818 }
1819 
1820 
1821 
1822 template<int dim, int spacedim>
1823 Point<dim>
1826  const Point<spacedim> &p) const
1827 {
1828  // Use an exact formula if one is available. this is only the case
1829  // for Q1 mappings in 1d, and in 2d if dim==spacedim
1830  if ((polynomial_degree == 1) &&
1831  ((dim == 1)
1832  ||
1833  ((dim == 2) && (dim == spacedim))))
1834  {
1835  // The dimension-dependent algorithms are much faster (about 25-45x in
1836  // 2D) but fail most of the time when the given point (p) is not in the
1837  // cell. The dimension-independent Newton algorithm given below is
1838  // slower, but more robust (though it still sometimes fails). Therefore
1839  // this function implements the following strategy based on the
1840  // p's dimension:
1841  //
1842  // * In 1D this mapping is linear, so the mapping is always invertible
1843  // (and the exact formula is known) as long as the cell has non-zero
1844  // length.
1845  // * In 2D the exact (quadratic) formula is called first. If either the
1846  // exact formula does not succeed (negative discriminant in the
1847  // quadratic formula) or succeeds but finds a solution outside of the
1848  // unit cell, then the Newton solver is called. The rationale for the
1849  // second choice is that the exact formula may provide two different
1850  // answers when mapping a point outside of the real cell, but the
1851  // Newton solver (if it converges) will only return one answer.
1852  // Otherwise the exact formula successfully found a point in the unit
1853  // cell and that value is returned.
1854  // * In 3D there is no (known to the authors) exact formula, so the Newton
1855  // algorithm is used.
1856  const std_cxx11::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
1857  vertices = this->get_vertices(cell);
1858  try
1859  {
1860  switch (dim)
1861  {
1862  case 1:
1863  {
1864  // formula not subject to any issues in 1d
1865  if (spacedim == 1)
1866  return internal::MappingQ1::transform_real_to_unit_cell(vertices, p);
1867  else
1868  {
1869  const std::vector<Point<spacedim> > a (vertices.begin(),
1870  vertices.end());
1871  return internal::MappingQ1::transform_real_to_unit_cell_initial_guess<dim,spacedim>(a,p);
1872  }
1873  }
1874 
1875  case 2:
1876  {
1877  const Point<dim> point
1878  = internal::MappingQ1::transform_real_to_unit_cell(vertices, p);
1879 
1880  // formula not guaranteed to work for points outside of
1881  // the cell. only take the computed point if it lies
1882  // inside the reference cell
1883  const double eps = 1e-15;
1884  if (-eps <= point(1) && point(1) <= 1 + eps &&
1885  -eps <= point(0) && point(0) <= 1 + eps)
1886  {
1887  return point;
1888  }
1889  else
1890  break;
1891  }
1892 
1893  default:
1894  {
1895  // we should get here, based on the if-condition at the top
1896  Assert(false, ExcInternalError());
1897  }
1898  }
1899  }
1901  {
1902  // simply fall through and continue on to the standard Newton code
1903  }
1904  }
1905  else
1906  {
1907  // we can't use an explicit formula,
1908  }
1909 
1910 
1911  Point<dim> initial_p_unit;
1912  if (polynomial_degree == 1)
1913  {
1914  // Find the initial value for the Newton iteration by a normal
1915  // projection to the least square plane determined by the vertices
1916  // of the cell
1917  const std::vector<Point<spacedim> > a
1918  = this->compute_mapping_support_points (cell);
1920  ExcInternalError());
1921  initial_p_unit = internal::MappingQ1::transform_real_to_unit_cell_initial_guess<dim,spacedim>(a,p);
1922  }
1923  else
1924  {
1925  try
1926  {
1927  // Find the initial value for the Newton iteration by a normal
1928  // projection to the least square plane determined by the vertices
1929  // of the cell
1930  //
1931  // we do this by first getting all support points, then
1932  // throwing away all but the vertices, and finally calling
1933  // the same function as above
1934  std::vector<Point<spacedim> > a
1935  = this->compute_mapping_support_points (cell);
1937  initial_p_unit = internal::MappingQ1::transform_real_to_unit_cell_initial_guess<dim,spacedim>(a,p);
1938  }
1939  catch (const typename Mapping<dim,spacedim>::ExcTransformationFailed &)
1940  {
1941  for (unsigned int d=0; d<dim; ++d)
1942  initial_p_unit[d] = 0.5;
1943  }
1944 
1945  // in case the function above should have given us something
1946  // back that lies outside the unit cell (that might happen
1947  // because we may have given a point 'p' that lies inside the
1948  // cell with the higher order mapping, but outside the Q1-mapped
1949  // reference cell), then project it back into the reference cell
1950  // in hopes that this gives a better starting point to the
1951  // following iteration
1952  initial_p_unit = GeometryInfo<dim>::project_to_unit_cell(initial_p_unit);
1953  }
1954 
1955  // perform the Newton iteration and return the result. note that
1956  // this statement may throw an exception, which we simply pass up to
1957  // the caller
1958  return this->transform_real_to_unit_cell_internal(cell, p, initial_p_unit);
1959 }
1960 
1961 
1962 
1963 template<int dim, int spacedim>
1966 {
1967  // add flags if the respective quantities are necessary to compute
1968  // what we need. note that some flags appear in both the conditions
1969  // and in subsequent set operations. this leads to some circular
1970  // logic. the only way to treat this is to iterate. since there are
1971  // 5 if-clauses in the loop, it will take at most 5 iterations to
1972  // converge. do them:
1973  UpdateFlags out = in;
1974  for (unsigned int i=0; i<5; ++i)
1975  {
1976  // The following is a little incorrect:
1977  // If not applied on a face,
1978  // update_boundary_forms does not
1979  // make sense. On the other hand,
1980  // it is necessary on a
1981  // face. Currently,
1982  // update_boundary_forms is simply
1983  // ignored for the interior of a
1984  // cell.
1985  if (out & (update_JxW_values
1987  out |= update_boundary_forms;
1988 
1996 
1997  if (out & (update_inverse_jacobians
2002 
2003  // The contravariant transformation is used in the Piola
2004  // transformation, which requires the determinant of the Jacobi
2005  // matrix of the transformation. Because we have no way of
2006  // knowing here whether the finite element wants to use the
2007  // contravariant or the Piola transforms, we add the JxW values
2008  // to the list of flags to be updated for each cell.
2010  out |= update_volume_elements;
2011 
2012  // the same is true when computing normal vectors: they require
2013  // the determinant of the Jacobian
2014  if (out & update_normal_vectors)
2015  out |= update_volume_elements;
2016  }
2017 
2018  return out;
2019 }
2020 
2021 
2022 
2023 template<int dim, int spacedim>
2026  const Quadrature<dim> &q) const
2027 {
2028  InternalData *data = new InternalData(polynomial_degree);
2029  data->initialize (this->requires_update_flags(update_flags), q, q.size());
2030 
2031  return data;
2032 }
2033 
2034 
2035 
2036 template<int dim, int spacedim>
2039  const Quadrature<dim-1> &quadrature) const
2040 {
2041  InternalData *data = new InternalData(polynomial_degree);
2042  data->initialize_face (this->requires_update_flags(update_flags),
2044  quadrature.size());
2045 
2046  return data;
2047 }
2048 
2049 
2050 
2051 template<int dim, int spacedim>
2054  const Quadrature<dim-1>& quadrature) const
2055 {
2056  InternalData *data = new InternalData(polynomial_degree);
2057  data->initialize_face (this->requires_update_flags(update_flags),
2059  quadrature.size());
2060 
2061  return data;
2062 }
2063 
2064 
2065 
2066 namespace internal
2067 {
2068  namespace
2069  {
2076  template <int dim, int spacedim>
2077  void
2078  maybe_compute_q_points (const typename QProjector<dim>::DataSetDescriptor data_set,
2079  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2080  std::vector<Point<spacedim> > &quadrature_points)
2081  {
2082  const UpdateFlags update_flags = data.update_each;
2083 
2084  if (update_flags & update_quadrature_points)
2085  {
2086  for (unsigned int point=0; point<quadrature_points.size(); ++point)
2087  {
2088  const double *shape = &data.shape(point+data_set,0);
2089  Point<spacedim> result = (shape[0] *
2090  data.mapping_support_points[0]);
2091  for (unsigned int k=1; k<data.n_shape_functions; ++k)
2092  for (unsigned int i=0; i<spacedim; ++i)
2093  result[i] += shape[k] * data.mapping_support_points[k][i];
2094  quadrature_points[point] = result;
2095  }
2096  }
2097  }
2098 
2099 
2107  template <int dim, int spacedim>
2108  void
2109  maybe_update_Jacobians (const CellSimilarity::Similarity cell_similarity,
2110  const typename ::QProjector<dim>::DataSetDescriptor data_set,
2111  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data)
2112  {
2113  const UpdateFlags update_flags = data.update_each;
2114 
2115  if (update_flags & update_contravariant_transformation)
2116  // if the current cell is just a
2117  // translation of the previous one, no
2118  // need to recompute jacobians...
2119  if (cell_similarity != CellSimilarity::translation)
2120  {
2121  const unsigned int n_q_points = data.contravariant.size();
2122 
2123  std::fill(data.contravariant.begin(), data.contravariant.end(),
2125 
2126  Assert (data.n_shape_functions > 0, ExcInternalError());
2127  const Tensor<1,spacedim> *supp_pts =
2128  &data.mapping_support_points[0];
2129 
2130  for (unsigned int point=0; point<n_q_points; ++point)
2131  {
2132  const Tensor<1,dim> *data_derv =
2133  &data.derivative(point+data_set, 0);
2134 
2135  double result [spacedim][dim];
2136 
2137  // peel away part of sum to avoid zeroing the
2138  // entries and adding for the first time
2139  for (unsigned int i=0; i<spacedim; ++i)
2140  for (unsigned int j=0; j<dim; ++j)
2141  result[i][j] = data_derv[0][j] * supp_pts[0][i];
2142  for (unsigned int k=1; k<data.n_shape_functions; ++k)
2143  for (unsigned int i=0; i<spacedim; ++i)
2144  for (unsigned int j=0; j<dim; ++j)
2145  result[i][j] += data_derv[k][j] * supp_pts[k][i];
2146 
2147  // write result into contravariant data. for
2148  // j=dim in the case dim<spacedim, there will
2149  // never be any nonzero data that arrives in
2150  // here, so it is ok anyway because it was
2151  // initialized to zero at the initialization
2152  for (unsigned int i=0; i<spacedim; ++i)
2153  for (unsigned int j=0; j<dim; ++j)
2154  data.contravariant[point][i][j] = result[i][j];
2155  }
2156  }
2157 
2158  if (update_flags & update_covariant_transformation)
2159  if (cell_similarity != CellSimilarity::translation)
2160  {
2161  const unsigned int n_q_points = data.contravariant.size();
2162  for (unsigned int point=0; point<n_q_points; ++point)
2163  {
2164  data.covariant[point] = (data.contravariant[point]).covariant_form();
2165  }
2166  }
2167 
2168  if (update_flags & update_volume_elements)
2169  if (cell_similarity != CellSimilarity::translation)
2170  {
2171  const unsigned int n_q_points = data.contravariant.size();
2172  for (unsigned int point=0; point<n_q_points; ++point)
2173  data.volume_elements[point] = data.contravariant[point].determinant();
2174  }
2175 
2176  }
2177 
2184  template <int dim, int spacedim>
2185  void
2186  maybe_update_jacobian_grads (const CellSimilarity::Similarity cell_similarity,
2187  const typename QProjector<dim>::DataSetDescriptor data_set,
2188  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2189  std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads)
2190  {
2191  const UpdateFlags update_flags = data.update_each;
2192  if (update_flags & update_jacobian_grads)
2193  {
2194  const unsigned int n_q_points = jacobian_grads.size();
2195 
2196  if (cell_similarity != CellSimilarity::translation)
2197  {
2198  for (unsigned int point=0; point<n_q_points; ++point)
2199  {
2200  const Tensor<2,dim> *second =
2201  &data.second_derivative(point+data_set, 0);
2202  double result [spacedim][dim][dim];
2203  for (unsigned int i=0; i<spacedim; ++i)
2204  for (unsigned int j=0; j<dim; ++j)
2205  for (unsigned int l=0; l<dim; ++l)
2206  result[i][j][l] = (second[0][j][l] *
2207  data.mapping_support_points[0][i]);
2208  for (unsigned int k=1; k<data.n_shape_functions; ++k)
2209  for (unsigned int i=0; i<spacedim; ++i)
2210  for (unsigned int j=0; j<dim; ++j)
2211  for (unsigned int l=0; l<dim; ++l)
2212  result[i][j][l]
2213  += (second[k][j][l]
2214  *
2215  data.mapping_support_points[k][i]);
2216 
2217  for (unsigned int i=0; i<spacedim; ++i)
2218  for (unsigned int j=0; j<dim; ++j)
2219  for (unsigned int l=0; l<dim; ++l)
2220  jacobian_grads[point][i][j][l] = result[i][j][l];
2221  }
2222  }
2223  }
2224  }
2225 
2232  template <int dim, int spacedim>
2233  void
2234  maybe_update_jacobian_pushed_forward_grads (const CellSimilarity::Similarity cell_similarity,
2235  const typename QProjector<dim>::DataSetDescriptor data_set,
2236  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2237  std::vector<Tensor<3,spacedim> > &jacobian_pushed_forward_grads)
2238  {
2239  const UpdateFlags update_flags = data.update_each;
2240  if (update_flags & update_jacobian_pushed_forward_grads)
2241  {
2242  const unsigned int n_q_points = jacobian_pushed_forward_grads.size();
2243 
2244  if (cell_similarity != CellSimilarity::translation)
2245  {
2246  double tmp[spacedim][spacedim][spacedim];
2247  for (unsigned int point=0; point<n_q_points; ++point)
2248  {
2249  const Tensor<2,dim> *second =
2250  &data.second_derivative(point+data_set, 0);
2251  double result [spacedim][dim][dim];
2252  for (unsigned int i=0; i<spacedim; ++i)
2253  for (unsigned int j=0; j<dim; ++j)
2254  for (unsigned int l=0; l<dim; ++l)
2255  result[i][j][l] = (second[0][j][l] *
2256  data.mapping_support_points[0][i]);
2257  for (unsigned int k=1; k<data.n_shape_functions; ++k)
2258  for (unsigned int i=0; i<spacedim; ++i)
2259  for (unsigned int j=0; j<dim; ++j)
2260  for (unsigned int l=0; l<dim; ++l)
2261  result[i][j][l]
2262  += (second[k][j][l]
2263  *
2264  data.mapping_support_points[k][i]);
2265 
2266  // first push forward the j-components
2267  for (unsigned int i=0; i<spacedim; ++i)
2268  for (unsigned int j=0; j<spacedim; ++j)
2269  for (unsigned int l=0; l<dim; ++l)
2270  {
2271  tmp[i][j][l] = result[i][0][l] *
2272  data.covariant[point][j][0];
2273  for (unsigned int jr=1; jr<dim; ++jr)
2274  {
2275  tmp[i][j][l] += result[i][jr][l] *
2276  data.covariant[point][j][jr];
2277  }
2278  }
2279 
2280  // now, pushing forward the l-components
2281  for (unsigned int i=0; i<spacedim; ++i)
2282  for (unsigned int j=0; j<spacedim; ++j)
2283  for (unsigned int l=0; l<spacedim; ++l)
2284  {
2285  jacobian_pushed_forward_grads[point][i][j][l] = tmp[i][j][0] *
2286  data.covariant[point][l][0];
2287  for (unsigned int lr=1; lr<dim; ++lr)
2288  {
2289  jacobian_pushed_forward_grads[point][i][j][l] += tmp[i][j][lr] *
2290  data.covariant[point][l][lr];
2291  }
2292 
2293  }
2294  }
2295  }
2296  }
2297  }
2298 
2305  template <int dim, int spacedim>
2306  void
2307  maybe_update_jacobian_2nd_derivatives (const CellSimilarity::Similarity cell_similarity,
2308  const typename QProjector<dim>::DataSetDescriptor data_set,
2309  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2310  std::vector<DerivativeForm<3,dim,spacedim> > &jacobian_2nd_derivatives)
2311  {
2312  const UpdateFlags update_flags = data.update_each;
2313  if (update_flags & update_jacobian_2nd_derivatives)
2314  {
2315  const unsigned int n_q_points = jacobian_2nd_derivatives.size();
2316 
2317  if (cell_similarity != CellSimilarity::translation)
2318  {
2319  for (unsigned int point=0; point<n_q_points; ++point)
2320  {
2321  const Tensor<3,dim> *third =
2322  &data.third_derivative(point+data_set, 0);
2323  double result [spacedim][dim][dim][dim];
2324  for (unsigned int i=0; i<spacedim; ++i)
2325  for (unsigned int j=0; j<dim; ++j)
2326  for (unsigned int l=0; l<dim; ++l)
2327  for (unsigned int m=0; m<dim; ++m)
2328  result[i][j][l][m] = (third[0][j][l][m] *
2329  data.mapping_support_points[0][i]);
2330  for (unsigned int k=1; k<data.n_shape_functions; ++k)
2331  for (unsigned int i=0; i<spacedim; ++i)
2332  for (unsigned int j=0; j<dim; ++j)
2333  for (unsigned int l=0; l<dim; ++l)
2334  for (unsigned int m=0; m<dim; ++m)
2335  result[i][j][l][m]
2336  += (third[k][j][l][m]
2337  *
2338  data.mapping_support_points[k][i]);
2339 
2340  for (unsigned int i=0; i<spacedim; ++i)
2341  for (unsigned int j=0; j<dim; ++j)
2342  for (unsigned int l=0; l<dim; ++l)
2343  for (unsigned int m=0; m<dim; ++m)
2344  jacobian_2nd_derivatives[point][i][j][l][m] = result[i][j][l][m];
2345  }
2346  }
2347  }
2348  }
2349 
2357  template <int dim, int spacedim>
2358  void
2359  maybe_update_jacobian_pushed_forward_2nd_derivatives (const CellSimilarity::Similarity cell_similarity,
2360  const typename QProjector<dim>::DataSetDescriptor data_set,
2361  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2362  std::vector<Tensor<4,spacedim> > &jacobian_pushed_forward_2nd_derivatives)
2363  {
2364  const UpdateFlags update_flags = data.update_each;
2366  {
2367  const unsigned int n_q_points = jacobian_pushed_forward_2nd_derivatives.size();
2368 
2369  if (cell_similarity != CellSimilarity::translation)
2370  {
2371  double tmp[spacedim][spacedim][spacedim][spacedim];
2372  for (unsigned int point=0; point<n_q_points; ++point)
2373  {
2374  const Tensor<3,dim> *third =
2375  &data.third_derivative(point+data_set, 0);
2376  double result [spacedim][dim][dim][dim];
2377  for (unsigned int i=0; i<spacedim; ++i)
2378  for (unsigned int j=0; j<dim; ++j)
2379  for (unsigned int l=0; l<dim; ++l)
2380  for (unsigned int m=0; m<dim; ++m)
2381  result[i][j][l][m] = (third[0][j][l][m] *
2382  data.mapping_support_points[0][i]);
2383  for (unsigned int k=1; k<data.n_shape_functions; ++k)
2384  for (unsigned int i=0; i<spacedim; ++i)
2385  for (unsigned int j=0; j<dim; ++j)
2386  for (unsigned int l=0; l<dim; ++l)
2387  for (unsigned int m=0; m<dim; ++m)
2388  result[i][j][l][m]
2389  += (third[k][j][l][m]
2390  *
2391  data.mapping_support_points[k][i]);
2392 
2393  // push forward the j-coordinate
2394  for (unsigned int i=0; i<spacedim; ++i)
2395  for (unsigned int j=0; j<spacedim; ++j)
2396  for (unsigned int l=0; l<dim; ++l)
2397  for (unsigned int m=0; m<dim; ++m)
2398  {
2399  jacobian_pushed_forward_2nd_derivatives[point][i][j][l][m]
2400  = result[i][0][l][m]*
2401  data.covariant[point][j][0];
2402  for (unsigned int jr=1; jr<dim; ++jr)
2403  jacobian_pushed_forward_2nd_derivatives[point][i][j][l][m]
2404  += result[i][jr][l][m]*
2405  data.covariant[point][j][jr];
2406  }
2407 
2408  // push forward the l-coordinate
2409  for (unsigned int i=0; i<spacedim; ++i)
2410  for (unsigned int j=0; j<spacedim; ++j)
2411  for (unsigned int l=0; l<spacedim; ++l)
2412  for (unsigned int m=0; m<dim; ++m)
2413  {
2414  tmp[i][j][l][m]
2415  = jacobian_pushed_forward_2nd_derivatives[point][i][j][0][m]*
2416  data.covariant[point][l][0];
2417  for (unsigned int lr=1; lr<dim; ++lr)
2418  tmp[i][j][l][m]
2419  += jacobian_pushed_forward_2nd_derivatives[point][i][j][lr][m]*
2420  data.covariant[point][l][lr];
2421  }
2422 
2423  // push forward the m-coordinate
2424  for (unsigned int i=0; i<spacedim; ++i)
2425  for (unsigned int j=0; j<spacedim; ++j)
2426  for (unsigned int l=0; l<spacedim; ++l)
2427  for (unsigned int m=0; m<spacedim; ++m)
2428  {
2429  jacobian_pushed_forward_2nd_derivatives[point][i][j][l][m]
2430  = tmp[i][j][l][0]*
2431  data.covariant[point][m][0];
2432  for (unsigned int mr=1; mr<dim; ++mr)
2433  jacobian_pushed_forward_2nd_derivatives[point][i][j][l][m]
2434  += tmp[i][j][l][mr]*
2435  data.covariant[point][m][mr];
2436  }
2437  }
2438  }
2439  }
2440  }
2441 
2448  template <int dim, int spacedim>
2449  void
2451  const typename QProjector<dim>::DataSetDescriptor data_set,
2452  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2453  std::vector<DerivativeForm<4,dim,spacedim> > &jacobian_3rd_derivatives)
2454  {
2455  const UpdateFlags update_flags = data.update_each;
2456  if (update_flags & update_jacobian_3rd_derivatives)
2457  {
2458  const unsigned int n_q_points = jacobian_3rd_derivatives.size();
2459 
2460  if (cell_similarity != CellSimilarity::translation)
2461  {
2462  for (unsigned int point=0; point<n_q_points; ++point)
2463  {
2464  const Tensor<4,dim> *fourth =
2465  &data.fourth_derivative(point+data_set, 0);
2466  double result [spacedim][dim][dim][dim][dim];
2467  for (unsigned int i=0; i<spacedim; ++i)
2468  for (unsigned int j=0; j<dim; ++j)
2469  for (unsigned int l=0; l<dim; ++l)
2470  for (unsigned int m=0; m<dim; ++m)
2471  for (unsigned int n=0; n<dim; ++n)
2472  result[i][j][l][m][n] = (fourth[0][j][l][m][n] *
2473  data.mapping_support_points[0][i]);
2474  for (unsigned int k=1; k<data.n_shape_functions; ++k)
2475  for (unsigned int i=0; i<spacedim; ++i)
2476  for (unsigned int j=0; j<dim; ++j)
2477  for (unsigned int l=0; l<dim; ++l)
2478  for (unsigned int m=0; m<dim; ++m)
2479  for (unsigned int n=0; n<dim; ++n)
2480  result[i][j][l][m][n]
2481  += (fourth[k][j][l][m][n]
2482  *
2483  data.mapping_support_points[k][i]);
2484 
2485  for (unsigned int i=0; i<spacedim; ++i)
2486  for (unsigned int j=0; j<dim; ++j)
2487  for (unsigned int l=0; l<dim; ++l)
2488  for (unsigned int m=0; m<dim; ++m)
2489  for (unsigned int n=0; n<dim; ++n)
2490  jacobian_3rd_derivatives[point][i][j][l][m][n] = result[i][j][l][m][n];
2491  }
2492  }
2493  }
2494  }
2495 
2502  template <int dim, int spacedim>
2503  void
2505  const typename QProjector<dim>::DataSetDescriptor data_set,
2506  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2507  std::vector<Tensor<5,spacedim> > &jacobian_pushed_forward_3rd_derivatives)
2508  {
2509  const UpdateFlags update_flags = data.update_each;
2511  {
2512  const unsigned int n_q_points = jacobian_pushed_forward_3rd_derivatives.size();
2513 
2514  if (cell_similarity != CellSimilarity::translation)
2515  {
2516  double tmp[spacedim][spacedim][spacedim][spacedim][spacedim];
2517  for (unsigned int point=0; point<n_q_points; ++point)
2518  {
2519  const Tensor<4,dim> *fourth =
2520  &data.fourth_derivative(point+data_set, 0);
2521  double result [spacedim][dim][dim][dim][dim];
2522  for (unsigned int i=0; i<spacedim; ++i)
2523  for (unsigned int j=0; j<dim; ++j)
2524  for (unsigned int l=0; l<dim; ++l)
2525  for (unsigned int m=0; m<dim; ++m)
2526  for (unsigned int n=0; n<dim; ++n)
2527  result[i][j][l][m][n] = (fourth[0][j][l][m][n] *
2528  data.mapping_support_points[0][i]);
2529  for (unsigned int k=1; k<data.n_shape_functions; ++k)
2530  for (unsigned int i=0; i<spacedim; ++i)
2531  for (unsigned int j=0; j<dim; ++j)
2532  for (unsigned int l=0; l<dim; ++l)
2533  for (unsigned int m=0; m<dim; ++m)
2534  for (unsigned int n=0; n<dim; ++n)
2535  result[i][j][l][m][n]
2536  += (fourth[k][j][l][m][n]
2537  *
2538  data.mapping_support_points[k][i]);
2539 
2540  // push-forward the j-coordinate
2541  for (unsigned int i=0; i<spacedim; ++i)
2542  for (unsigned int j=0; j<spacedim; ++j)
2543  for (unsigned int l=0; l<dim; ++l)
2544  for (unsigned int m=0; m<dim; ++m)
2545  for (unsigned int n=0; n<dim; ++n)
2546  {
2547  tmp[i][j][l][m][n] = result[i][0][l][m][n] *
2548  data.covariant[point][j][0];
2549  for (unsigned int jr=1; jr<dim; ++jr)
2550  tmp[i][j][l][m][n] += result[i][jr][l][m][n] *
2551  data.covariant[point][j][jr];
2552  }
2553 
2554  // push-forward the l-coordinate
2555  for (unsigned int i=0; i<spacedim; ++i)
2556  for (unsigned int j=0; j<spacedim; ++j)
2557  for (unsigned int l=0; l<spacedim; ++l)
2558  for (unsigned int m=0; m<dim; ++m)
2559  for (unsigned int n=0; n<dim; ++n)
2560  {
2561  jacobian_pushed_forward_3rd_derivatives[point][i][j][l][m][n]
2562  = tmp[i][j][0][m][n] *
2563  data.covariant[point][l][0];
2564  for (unsigned int lr=1; lr<dim; ++lr)
2565  jacobian_pushed_forward_3rd_derivatives[point][i][j][l][m][n]
2566  += tmp[i][j][lr][m][n] *
2567  data.covariant[point][l][lr];
2568  }
2569 
2570  // push-forward the m-coordinate
2571  for (unsigned int i=0; i<spacedim; ++i)
2572  for (unsigned int j=0; j<spacedim; ++j)
2573  for (unsigned int l=0; l<spacedim; ++l)
2574  for (unsigned int m=0; m<spacedim; ++m)
2575  for (unsigned int n=0; n<dim; ++n)
2576  {
2577  tmp[i][j][l][m][n]
2578  = jacobian_pushed_forward_3rd_derivatives[point][i][j][l][0][n] *
2579  data.covariant[point][m][0];
2580  for (unsigned int mr=1; mr<dim; ++mr)
2581  tmp[i][j][l][m][n]
2582  += jacobian_pushed_forward_3rd_derivatives[point][i][j][l][mr][n] *
2583  data.covariant[point][m][mr];
2584  }
2585 
2586  // push-forward the n-coordinate
2587  for (unsigned int i=0; i<spacedim; ++i)
2588  for (unsigned int j=0; j<spacedim; ++j)
2589  for (unsigned int l=0; l<spacedim; ++l)
2590  for (unsigned int m=0; m<spacedim; ++m)
2591  for (unsigned int n=0; n<spacedim; ++n)
2592  {
2593  jacobian_pushed_forward_3rd_derivatives[point][i][j][l][m][n]
2594  = tmp[i][j][l][m][0] *
2595  data.covariant[point][n][0];
2596  for (unsigned int nr=1; nr<dim; ++nr)
2597  jacobian_pushed_forward_3rd_derivatives[point][i][j][l][m][n]
2598  += tmp[i][j][l][m][nr] *
2599  data.covariant[point][n][nr];
2600  }
2601  }
2602  }
2603  }
2604  }
2605  }
2606 }
2607 
2608 
2609 
2610 
2611 template<int dim, int spacedim>
2615  const CellSimilarity::Similarity cell_similarity,
2616  const Quadrature<dim> &quadrature,
2617  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
2619 {
2620  // ensure that the following static_cast is really correct:
2621  Assert (dynamic_cast<const InternalData *>(&internal_data) != 0,
2622  ExcInternalError());
2623  const InternalData &data = static_cast<const InternalData &>(internal_data);
2624 
2625  const unsigned int n_q_points=quadrature.size();
2626 
2627  // recompute the support points of the transformation of this
2628  // cell. we tried to be clever here in an earlier version of the
2629  // library by checking whether the cell is the same as the one we
2630  // had visited last, but it turns out to be difficult to determine
2631  // that because a cell for the purposes of a mapping is
2632  // characterized not just by its (triangulation, level, index)
2633  // triple, but also by the locations of its vertices, the manifold
2634  // object attached to the cell and all of its bounding faces/edges,
2635  // etc. to reliably test that the "cell" we are on is, therefore,
2636  // not easily done
2637  data.mapping_support_points = this->compute_mapping_support_points(cell);
2638  data.cell_of_current_support_points = cell;
2639 
2640  // if the order of the mapping is greater than 1, then do not reuse any cell
2641  // similarity information. This is necessary because the cell similarity
2642  // value is computed with just cell vertices and does not take into account
2643  // cell curvature.
2644  const CellSimilarity::Similarity computed_cell_similarity =
2645  (polynomial_degree == 1 ? cell_similarity : CellSimilarity::none);
2646 
2647  internal::maybe_compute_q_points<dim,spacedim> (QProjector<dim>::DataSetDescriptor::cell (),
2648  data,
2649  output_data.quadrature_points);
2650  internal::maybe_update_Jacobians<dim,spacedim> (computed_cell_similarity,
2652  data);
2653 
2654  const UpdateFlags update_flags = data.update_each;
2655  const std::vector<double> &weights=quadrature.get_weights();
2656 
2657  // Multiply quadrature weights by absolute value of Jacobian determinants or
2658  // the area element g=sqrt(DX^t DX) in case of codim > 0
2659 
2660  if (update_flags & (update_normal_vectors
2661  | update_JxW_values))
2662  {
2663  AssertDimension (output_data.JxW_values.size(), n_q_points);
2664 
2665  Assert( !(update_flags & update_normal_vectors ) ||
2666  (output_data.normal_vectors.size() == n_q_points),
2667  ExcDimensionMismatch(output_data.normal_vectors.size(), n_q_points));
2668 
2669 
2670  if (computed_cell_similarity != CellSimilarity::translation)
2671  for (unsigned int point=0; point<n_q_points; ++point)
2672  {
2673 
2674  if (dim == spacedim)
2675  {
2676  const double det = data.contravariant[point].determinant();
2677 
2678  // check for distorted cells.
2679 
2680  // TODO: this allows for anisotropies of up to 1e6 in 3D and
2681  // 1e12 in 2D. might want to find a finer
2682  // (dimension-independent) criterion
2683  Assert (det > 1e-12*Utilities::fixed_power<dim>(cell->diameter()/
2684  std::sqrt(double(dim))),
2685  (typename Mapping<dim,spacedim>::ExcDistortedMappedCell(cell->center(), det, point)));
2686 
2687  output_data.JxW_values[point] = weights[point] * det;
2688  }
2689  // if dim==spacedim, then there is no cell normal to
2690  // compute. since this is for FEValues (and not FEFaceValues),
2691  // there are also no face normals to compute
2692  else //codim>0 case
2693  {
2694  Tensor<1, spacedim> DX_t [dim];
2695  for (unsigned int i=0; i<spacedim; ++i)
2696  for (unsigned int j=0; j<dim; ++j)
2697  DX_t[j][i] = data.contravariant[point][i][j];
2698 
2699  Tensor<2, dim> G; //First fundamental form
2700  for (unsigned int i=0; i<dim; ++i)
2701  for (unsigned int j=0; j<dim; ++j)
2702  G[i][j] = DX_t[i] * DX_t[j];
2703 
2704  output_data.JxW_values[point]
2705  = sqrt(determinant(G)) * weights[point];
2706 
2707  if (computed_cell_similarity == CellSimilarity::inverted_translation)
2708  {
2709  // we only need to flip the normal
2710  if (update_flags & update_normal_vectors)
2711  output_data.normal_vectors[point] *= -1.;
2712  }
2713  else
2714  {
2715  if (update_flags & update_normal_vectors)
2716  {
2717  Assert(spacedim == dim+1,
2718  ExcMessage("There is no (unique) cell normal for "
2719  + Utilities::int_to_string(dim) +
2720  "-dimensional cells in "
2721  + Utilities::int_to_string(spacedim) +
2722  "-dimensional space. This only works if the "
2723  "space dimension is one greater than the "
2724  "dimensionality of the mesh cells."));
2725 
2726  if (dim==1)
2727  output_data.normal_vectors[point] =
2728  cross_product_2d(-DX_t[0]);
2729  else //dim == 2
2730  output_data.normal_vectors[point] =
2731  cross_product_3d(DX_t[0], DX_t[1]);
2732 
2733  output_data.normal_vectors[point] /= output_data.normal_vectors[point].norm();
2734 
2735  if (cell->direction_flag() == false)
2736  output_data.normal_vectors[point] *= -1.;
2737  }
2738 
2739  }
2740  } //codim>0 case
2741 
2742  }
2743  }
2744 
2745 
2746 
2747  // copy values from InternalData to vector given by reference
2748  if (update_flags & update_jacobians)
2749  {
2750  AssertDimension (output_data.jacobians.size(), n_q_points);
2751  if (computed_cell_similarity != CellSimilarity::translation)
2752  for (unsigned int point=0; point<n_q_points; ++point)
2753  output_data.jacobians[point] = data.contravariant[point];
2754  }
2755 
2756  // copy values from InternalData to vector given by reference
2757  if (update_flags & update_inverse_jacobians)
2758  {
2759  AssertDimension (output_data.inverse_jacobians.size(), n_q_points);
2760  if (computed_cell_similarity != CellSimilarity::translation)
2761  for (unsigned int point=0; point<n_q_points; ++point)
2762  output_data.inverse_jacobians[point] = data.covariant[point].transpose();
2763  }
2764 
2765  internal::maybe_update_jacobian_grads<dim,spacedim> (computed_cell_similarity,
2767  data,
2768  output_data.jacobian_grads);
2769 
2770  internal::maybe_update_jacobian_pushed_forward_grads<dim,spacedim> (computed_cell_similarity,
2772  data,
2773  output_data.jacobian_pushed_forward_grads);
2774 
2775  internal::maybe_update_jacobian_2nd_derivatives<dim,spacedim> (computed_cell_similarity,
2777  data,
2778  output_data.jacobian_2nd_derivatives);
2779 
2780  internal::maybe_update_jacobian_pushed_forward_2nd_derivatives<dim,spacedim> (computed_cell_similarity,
2782  data,
2784 
2785  internal::maybe_update_jacobian_3rd_derivatives<dim,spacedim> (computed_cell_similarity,
2787  data,
2788  output_data.jacobian_3rd_derivatives);
2789 
2790  internal::maybe_update_jacobian_pushed_forward_3rd_derivatives<dim,spacedim> (computed_cell_similarity,
2792  data,
2794 
2795  return computed_cell_similarity;
2796 }
2797 
2798 
2799 
2800 
2801 
2802 
2803 namespace internal
2804 {
2805  namespace
2806  {
2816  template <int dim, int spacedim>
2817  void
2818  maybe_compute_face_data (const ::MappingQGeneric<dim,spacedim> &mapping,
2819  const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
2820  const unsigned int face_no,
2821  const unsigned int subface_no,
2822  const unsigned int n_q_points,
2823  const std::vector<double> &weights,
2824  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2826  {
2827  const UpdateFlags update_flags = data.update_each;
2828 
2829  if (update_flags & (update_boundary_forms |
2834  {
2835  if (update_flags & update_boundary_forms)
2836  AssertDimension (output_data.boundary_forms.size(), n_q_points);
2837  if (update_flags & update_normal_vectors)
2838  AssertDimension (output_data.normal_vectors.size(), n_q_points);
2839  if (update_flags & update_JxW_values)
2840  AssertDimension (output_data.JxW_values.size(), n_q_points);
2841 
2842  Assert (data.aux.size()+1 >= dim, ExcInternalError());
2843 
2844  // first compute some common data that is used for evaluating
2845  // all of the flags below
2846 
2847  // map the unit tangentials to the real cell. checking for d!=dim-1
2848  // eliminates compiler warnings regarding unsigned int expressions <
2849  // 0.
2850  for (unsigned int d=0; d!=dim-1; ++d)
2851  {
2853  data.unit_tangentials.size(),
2854  ExcInternalError());
2855  Assert (data.aux[d].size() <=
2856  data.unit_tangentials[face_no+GeometryInfo<dim>::faces_per_cell*d].size(),
2857  ExcInternalError());
2858 
2859  mapping.transform (make_array_view(data.unit_tangentials[face_no+GeometryInfo<dim>::faces_per_cell*d]),
2861  data,
2862  make_array_view(data.aux[d]));
2863  }
2864 
2865  if (update_flags & update_boundary_forms)
2866  {
2867  // if dim==spacedim, we can use the unit tangentials to compute the
2868  // boundary form by simply taking the cross product
2869  if (dim == spacedim)
2870  {
2871  for (unsigned int i=0; i<n_q_points; ++i)
2872  switch (dim)
2873  {
2874  case 1:
2875  // in 1d, we don't have access to any of the data.aux
2876  // fields (because it has only dim-1 components), but we
2877  // can still compute the boundary form by simply
2878  // looking at the number of the face
2879  output_data.boundary_forms[i][0] = (face_no == 0 ?
2880  -1 : +1);
2881  break;
2882  case 2:
2883  output_data.boundary_forms[i] =
2884  cross_product_2d(data.aux[0][i]);
2885  break;
2886  case 3:
2887  output_data.boundary_forms[i] =
2888  cross_product_3d(data.aux[0][i], data.aux[1][i]);
2889  break;
2890  default:
2891  Assert(false, ExcNotImplemented());
2892  }
2893  }
2894  else //(dim < spacedim)
2895  {
2896  // in the codim-one case, the boundary form results from the
2897  // cross product of all the face tangential vectors and the cell
2898  // normal vector
2899  //
2900  // to compute the cell normal, use the same method used in
2901  // fill_fe_values for cells above
2902  AssertDimension (data.contravariant.size(), n_q_points);
2903 
2904  for (unsigned int point=0; point<n_q_points; ++point)
2905  {
2906  if (dim==1)
2907  {
2908  // J is a tangent vector
2909  output_data.boundary_forms[point] = data.contravariant[point].transpose()[0];
2910  output_data.boundary_forms[point] /=
2911  (face_no == 0 ? -1. : +1.) * output_data.boundary_forms[point].norm();
2912  }
2913 
2914  if (dim==2)
2915  {
2916  const DerivativeForm<1,spacedim,dim> DX_t =
2917  data.contravariant[point].transpose();
2918 
2919  Tensor<1, spacedim> cell_normal =
2920  cross_product_3d(DX_t[0], DX_t[1]);
2921  cell_normal /= cell_normal.norm();
2922 
2923  // then compute the face normal from the face tangent
2924  // and the cell normal:
2925  output_data.boundary_forms[point] =
2926  cross_product_3d(data.aux[0][point], cell_normal);
2927  }
2928  }
2929  }
2930  }
2931 
2932  if (update_flags & update_JxW_values)
2933  for (unsigned int i=0; i<output_data.boundary_forms.size(); ++i)
2934  {
2935  output_data.JxW_values[i] = output_data.boundary_forms[i].norm() * weights[i];
2936 
2937  if (subface_no != numbers::invalid_unsigned_int)
2938  {
2939  const double area_ratio = GeometryInfo<dim>::subface_ratio(cell->subface_case(face_no),
2940  subface_no);
2941  output_data.JxW_values[i] *= area_ratio;
2942  }
2943  }
2944 
2945  if (update_flags & update_normal_vectors)
2946  for (unsigned int i=0; i<output_data.normal_vectors.size(); ++i)
2947  output_data.normal_vectors[i] = Point<spacedim>(output_data.boundary_forms[i] /
2948  output_data.boundary_forms[i].norm());
2949 
2950  if (update_flags & update_jacobians)
2951  for (unsigned int point=0; point<n_q_points; ++point)
2952  output_data.jacobians[point] = data.contravariant[point];
2953 
2954  if (update_flags & update_inverse_jacobians)
2955  for (unsigned int point=0; point<n_q_points; ++point)
2956  output_data.inverse_jacobians[point] = data.covariant[point].transpose();
2957  }
2958  }
2959 
2960 
2967  template<int dim, int spacedim>
2968  void
2969  do_fill_fe_face_values (const ::MappingQGeneric<dim,spacedim> &mapping,
2970  const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
2971  const unsigned int face_no,
2972  const unsigned int subface_no,
2973  const typename QProjector<dim>::DataSetDescriptor data_set,
2974  const Quadrature<dim-1> &quadrature,
2975  const typename ::MappingQGeneric<dim,spacedim>::InternalData &data,
2977  {
2978  maybe_compute_q_points<dim,spacedim> (data_set,
2979  data,
2980  output_data.quadrature_points);
2981  maybe_update_Jacobians<dim,spacedim> (CellSimilarity::none,
2982  data_set,
2983  data);
2984  maybe_update_jacobian_grads<dim,spacedim> (CellSimilarity::none,
2985  data_set,
2986  data,
2987  output_data.jacobian_grads);
2988  maybe_update_jacobian_pushed_forward_grads<dim,spacedim> (CellSimilarity::none,
2989  data_set,
2990  data,
2991  output_data.jacobian_pushed_forward_grads);
2992  maybe_update_jacobian_2nd_derivatives<dim,spacedim> (CellSimilarity::none,
2993  data_set,
2994  data,
2995  output_data.jacobian_2nd_derivatives);
2996  maybe_update_jacobian_pushed_forward_2nd_derivatives<dim,spacedim> (CellSimilarity::none,
2997  data_set,
2998  data,
3000  maybe_update_jacobian_3rd_derivatives<dim,spacedim> (CellSimilarity::none,
3001  data_set,
3002  data,
3003  output_data.jacobian_3rd_derivatives);
3004  maybe_update_jacobian_pushed_forward_3rd_derivatives<dim,spacedim> (CellSimilarity::none,
3005  data_set,
3006  data,
3008 
3009  maybe_compute_face_data (mapping,
3010  cell, face_no, subface_no, quadrature.size(),
3011  quadrature.get_weights(), data,
3012  output_data);
3013  }
3014  }
3015 }
3016 
3017 
3018 
3019 template<int dim, int spacedim>
3020 void
3023  const unsigned int face_no,
3024  const Quadrature<dim-1> &quadrature,
3025  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
3027 {
3028  // ensure that the following cast is really correct:
3029  Assert ((dynamic_cast<const InternalData *>(&internal_data) != 0),
3030  ExcInternalError());
3031  const InternalData &data
3032  = static_cast<const InternalData &>(internal_data);
3033 
3034  // if necessary, recompute the support points of the transformation of this cell
3035  // (note that we need to first check the triangulation pointer, since otherwise
3036  // the second test might trigger an exception if the triangulations are not the
3037  // same)
3038  if ((data.mapping_support_points.size() == 0)
3039  ||
3040  (&cell->get_triangulation() !=
3041  &data.cell_of_current_support_points->get_triangulation())
3042  ||
3043  (cell != data.cell_of_current_support_points))
3044  {
3045  data.mapping_support_points = this->compute_mapping_support_points(cell);
3046  data.cell_of_current_support_points = cell;
3047  }
3048 
3050  cell, face_no, numbers::invalid_unsigned_int,
3052  cell->face_orientation(face_no),
3053  cell->face_flip(face_no),
3054  cell->face_rotation(face_no),
3055  quadrature.size()),
3056  quadrature,
3057  data,
3058  output_data);
3059 }
3060 
3061 
3062 
3063 template<int dim, int spacedim>
3064 void
3067  const unsigned int face_no,
3068  const unsigned int subface_no,
3069  const Quadrature<dim-1> &quadrature,
3070  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
3072 {
3073  // ensure that the following cast is really correct:
3074  Assert ((dynamic_cast<const InternalData *>(&internal_data) != 0),
3075  ExcInternalError());
3076  const InternalData &data
3077  = static_cast<const InternalData &>(internal_data);
3078 
3079  // if necessary, recompute the support points of the transformation of this cell
3080  // (note that we need to first check the triangulation pointer, since otherwise
3081  // the second test might trigger an exception if the triangulations are not the
3082  // same)
3083  if ((data.mapping_support_points.size() == 0)
3084  ||
3085  (&cell->get_triangulation() !=
3086  &data.cell_of_current_support_points->get_triangulation())
3087  ||
3088  (cell != data.cell_of_current_support_points))
3089  {
3090  data.mapping_support_points = this->compute_mapping_support_points(cell);
3091  data.cell_of_current_support_points = cell;
3092  }
3093 
3095  cell, face_no, subface_no,
3096  QProjector<dim>::DataSetDescriptor::subface (face_no, subface_no,
3097  cell->face_orientation(face_no),
3098  cell->face_flip(face_no),
3099  cell->face_rotation(face_no),
3100  quadrature.size(),
3101  cell->subface_case(face_no)),
3102  quadrature,
3103  data,
3104  output_data);
3105 }
3106 
3107 
3108 
3109 namespace
3110 {
3111  template <int dim, int spacedim, int rank>
3112  void
3113  transform_fields(const ArrayView<const Tensor<rank,dim> > &input,
3114  const MappingType mapping_type,
3115  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3116  const ArrayView<Tensor<rank,spacedim> > &output)
3117  {
3118  AssertDimension (input.size(), output.size());
3119  Assert ((dynamic_cast<const typename MappingQGeneric<dim,spacedim>::InternalData *>(&mapping_data) != 0),
3120  ExcInternalError());
3122  &data = static_cast<const typename MappingQGeneric<dim,spacedim>::InternalData &>(mapping_data);
3123 
3124  switch (mapping_type)
3125  {
3126  case mapping_contravariant:
3127  {
3129  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
3130 
3131  for (unsigned int i=0; i<output.size(); ++i)
3132  output[i] = apply_transformation(data.contravariant[i], input[i]);
3133 
3134  return;
3135  }
3136 
3137  case mapping_piola:
3138  {
3140  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
3142  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
3143  Assert (rank==1, ExcMessage("Only for rank 1"));
3144  if (rank!=1)
3145  return;
3146 
3147  for (unsigned int i=0; i<output.size(); ++i)
3148  {
3149  output[i] = apply_transformation(data.contravariant[i], input[i]);
3150  output[i] /= data.volume_elements[i];
3151  }
3152  return;
3153  }
3154  //We still allow this operation as in the
3155  //reference cell Derivatives are Tensor
3156  //rather than DerivativeForm
3157  case mapping_covariant:
3158  {
3160  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3161 
3162  for (unsigned int i=0; i<output.size(); ++i)
3163  output[i] = apply_transformation(data.covariant[i], input[i]);
3164 
3165  return;
3166  }
3167 
3168  default:
3169  Assert(false, ExcNotImplemented());
3170  }
3171  }
3172 
3173 
3174  template <int dim, int spacedim, int rank>
3175  void
3176  transform_gradients(const ArrayView<const Tensor<rank,dim> > &input,
3177  const MappingType mapping_type,
3178  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3179  const ArrayView<Tensor<rank,spacedim> > &output)
3180  {
3181  AssertDimension (input.size(), output.size());
3182  Assert ((dynamic_cast<const typename MappingQGeneric<dim,spacedim>::InternalData *>(&mapping_data) != 0),
3183  ExcInternalError());
3185  &data = static_cast<const typename MappingQGeneric<dim,spacedim>::InternalData &>(mapping_data);
3186 
3187  switch (mapping_type)
3188  {
3190  {
3192  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3194  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
3195  Assert (rank==2, ExcMessage("Only for rank 2"));
3196 
3197  for (unsigned int i=0; i<output.size(); ++i)
3198  {
3200  apply_transformation(data.contravariant[i], transpose(input[i]) );
3201  output[i] = apply_transformation(data.covariant[i], A.transpose() );
3202  }
3203 
3204  return;
3205  }
3206 
3208  {
3210  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3211  Assert (rank==2, ExcMessage("Only for rank 2"));
3212 
3213  for (unsigned int i=0; i<output.size(); ++i)
3214  {
3216  apply_transformation(data.covariant[i], transpose(input[i]) );
3217  output[i] = apply_transformation(data.covariant[i], A.transpose() );
3218  }
3219 
3220  return;
3221  }
3222 
3224  {
3226  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3228  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
3230  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
3231  Assert (rank==2, ExcMessage("Only for rank 2"));
3232 
3233  for (unsigned int i=0; i<output.size(); ++i)
3234  {
3236  apply_transformation(data.covariant[i], input[i] );
3237  Tensor<2,spacedim> T =
3239 
3240  output[i] = transpose(T);
3241  output[i] /= data.volume_elements[i];
3242  }
3243 
3244  return;
3245  }
3246 
3247  default:
3248  Assert(false, ExcNotImplemented());
3249  }
3250  }
3251 
3252 
3253 
3254 
3255  template <int dim, int spacedim>
3256  void
3257  transform_hessians(const ArrayView<const Tensor<3,dim> > &input,
3258  const MappingType mapping_type,
3259  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3260  const ArrayView<Tensor<3,spacedim> > &output)
3261  {
3262  AssertDimension (input.size(), output.size());
3263  Assert ((dynamic_cast<const typename MappingQGeneric<dim,spacedim>::InternalData *>(&mapping_data) != 0),
3264  ExcInternalError());
3266  &data = static_cast<const typename MappingQGeneric<dim,spacedim>::InternalData &>(mapping_data);
3267 
3268  switch (mapping_type)
3269  {
3271  {
3273  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3275  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
3276 
3277  for (unsigned int q=0; q<output.size(); ++q)
3278  for (unsigned int i=0; i<spacedim; ++i)
3279  {
3280  double tmp1[dim][dim];
3281  for (unsigned int J=0; J<dim; ++J)
3282  for (unsigned int K=0; K<dim; ++K)
3283  {
3284  tmp1[J][K] = data.contravariant[q][i][0] * input[q][0][J][K];
3285  for (unsigned int I=1; I<dim; ++I)
3286  tmp1[J][K] += data.contravariant[q][i][I] * input[q][I][J][K];
3287  }
3288  for (unsigned int j=0; j<spacedim; ++j)
3289  {
3290  double tmp2[dim];
3291  for (unsigned int K=0; K<dim; ++K)
3292  {
3293  tmp2[K] = data.covariant[q][j][0] * tmp1[0][K];
3294  for (unsigned int J=1; J<dim; ++J)
3295  tmp2[K] += data.covariant[q][j][J] * tmp1[J][K];
3296  }
3297  for (unsigned int k=0; k<spacedim; ++k)
3298  {
3299  output[q][i][j][k] = data.covariant[q][k][0] * tmp2[0];
3300  for (unsigned int K=1; K<dim; ++K)
3301  output[q][i][j][k] += data.covariant[q][k][K] * tmp2[K];
3302  }
3303  }
3304  }
3305  return;
3306  }
3307 
3309  {
3311  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3312 
3313  for (unsigned int q=0; q<output.size(); ++q)
3314  for (unsigned int i=0; i<spacedim; ++i)
3315  {
3316  double tmp1[dim][dim];
3317  for (unsigned int J=0; J<dim; ++J)
3318  for (unsigned int K=0; K<dim; ++K)
3319  {
3320  tmp1[J][K] = data.covariant[q][i][0] * input[q][0][J][K];
3321  for (unsigned int I=1; I<dim; ++I)
3322  tmp1[J][K] += data.covariant[q][i][I] * input[q][I][J][K];
3323  }
3324  for (unsigned int j=0; j<spacedim; ++j)
3325  {
3326  double tmp2[dim];
3327  for (unsigned int K=0; K<dim; ++K)
3328  {
3329  tmp2[K] = data.covariant[q][j][0] * tmp1[0][K];
3330  for (unsigned int J=1; J<dim; ++J)
3331  tmp2[K] += data.covariant[q][j][J] * tmp1[J][K];
3332  }
3333  for (unsigned int k=0; k<spacedim; ++k)
3334  {
3335  output[q][i][j][k] = data.covariant[q][k][0] * tmp2[0];
3336  for (unsigned int K=1; K<dim; ++K)
3337  output[q][i][j][k] += data.covariant[q][k][K] * tmp2[K];
3338  }
3339  }
3340  }
3341 
3342  return;
3343  }
3344 
3345  case mapping_piola_hessian:
3346  {
3348  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3350  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
3352  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
3353 
3354  for (unsigned int q=0; q<output.size(); ++q)
3355  for (unsigned int i=0; i<spacedim; ++i)
3356  {
3357  double factor[dim];
3358  for (unsigned int I=0; I<dim; ++I)
3359  factor[I] = data.contravariant[q][i][I] / data.volume_elements[q];
3360  double tmp1[dim][dim];
3361  for (unsigned int J=0; J<dim; ++J)
3362  for (unsigned int K=0; K<dim; ++K)
3363  {
3364  tmp1[J][K] = factor[0] * input[q][0][J][K];
3365  for (unsigned int I=1; I<dim; ++I)
3366  tmp1[J][K] += factor[I] * input[q][I][J][K];
3367  }
3368  for (unsigned int j=0; j<spacedim; ++j)
3369  {
3370  double tmp2[dim];
3371  for (unsigned int K=0; K<dim; ++K)
3372  {
3373  tmp2[K] = data.covariant[q][j][0] * tmp1[0][K];
3374  for (unsigned int J=1; J<dim; ++J)
3375  tmp2[K] += data.covariant[q][j][J] * tmp1[J][K];
3376  }
3377  for (unsigned int k=0; k<spacedim; ++k)
3378  {
3379  output[q][i][j][k] = data.covariant[q][k][0] * tmp2[0];
3380  for (unsigned int K=1; K<dim; ++K)
3381  output[q][i][j][k] += data.covariant[q][k][K] * tmp2[K];
3382  }
3383  }
3384  }
3385 
3386  return;
3387  }
3388 
3389  default:
3390  Assert(false, ExcNotImplemented());
3391  }
3392  }
3393 
3394 
3395 
3396 
3397  template<int dim, int spacedim, int rank>
3398  void
3399  transform_differential_forms(const ArrayView<const DerivativeForm<rank, dim,spacedim> > &input,
3400  const MappingType mapping_type,
3401  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3402  const ArrayView<Tensor<rank+1, spacedim> > &output)
3403  {
3404  AssertDimension (input.size(), output.size());
3405  Assert ((dynamic_cast<const typename MappingQGeneric<dim,spacedim>::InternalData *>(&mapping_data) != 0),
3406  ExcInternalError());
3408  &data = static_cast<const typename MappingQGeneric<dim,spacedim>::InternalData &>(mapping_data);
3409 
3410  switch (mapping_type)
3411  {
3412  case mapping_covariant:
3413  {
3415  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3416 
3417  for (unsigned int i=0; i<output.size(); ++i)
3418  output[i] = apply_transformation(data.covariant[i], input[i]);
3419 
3420  return;
3421  }
3422  default:
3423  Assert(false, ExcNotImplemented());
3424  }
3425  }
3426 }
3427 
3428 
3429 
3430 template<int dim, int spacedim>
3431 void
3433 transform (const ArrayView<const Tensor<1, dim> > &input,
3434  const MappingType mapping_type,
3435  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3436  const ArrayView<Tensor<1, spacedim> > &output) const
3437 {
3438  transform_fields(input, mapping_type, mapping_data, output);
3439 }
3440 
3441 
3442 
3443 template<int dim, int spacedim>
3444 void
3447  const MappingType mapping_type,
3448  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3449  const ArrayView<Tensor<2, spacedim> > &output) const
3450 {
3451  transform_differential_forms(input, mapping_type, mapping_data, output);
3452 }
3453 
3454 
3455 
3456 template<int dim, int spacedim>
3457 void
3459 transform (const ArrayView<const Tensor<2, dim> > &input,
3460  const MappingType mapping_type,
3461  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3462  const ArrayView<Tensor<2, spacedim> > &output) const
3463 {
3464  switch (mapping_type)
3465  {
3466  case mapping_contravariant:
3467  transform_fields(input, mapping_type, mapping_data, output);
3468  return;
3469 
3473  transform_gradients(input, mapping_type, mapping_data, output);
3474  return;
3475  default:
3476  Assert(false, ExcNotImplemented());
3477  }
3478 }
3479 
3480 
3481 
3482 template<int dim, int spacedim>
3483 void
3486  const MappingType mapping_type,
3487  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3488  const ArrayView<Tensor<3,spacedim> > &output) const
3489 {
3490 
3491  AssertDimension (input.size(), output.size());
3492  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
3493  ExcInternalError());
3494  const InternalData &data = static_cast<const InternalData &>(mapping_data);
3495 
3496  switch (mapping_type)
3497  {
3499  {
3501  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
3502 
3503  for (unsigned int q=0; q<output.size(); ++q)
3504  for (unsigned int i=0; i<spacedim; ++i)
3505  for (unsigned int j=0; j<spacedim; ++j)
3506  {
3507  double tmp[dim];
3508  for (unsigned int K=0; K<dim; ++K)
3509  {
3510  tmp[K] = data.covariant[q][j][0] * input[q][i][0][K];
3511  for (unsigned int J=1; J<dim; ++J)
3512  tmp[K] += data.covariant[q][j][J] * input[q][i][J][K];
3513  }
3514  for (unsigned int k=0; k<spacedim; ++k)
3515  {
3516  output[q][i][j][k] = data.covariant[q][k][0] * tmp[0];
3517  for (unsigned int K=1; K<dim; ++K)
3518  output[q][i][j][k] += data.covariant[q][k][K] * tmp[K];
3519  }
3520  }
3521  return;
3522  }
3523 
3524  default:
3525  Assert(false, ExcNotImplemented());
3526  }
3527 }
3528 
3529 
3530 
3531 template<int dim, int spacedim>
3532 void
3534 transform (const ArrayView<const Tensor<3,dim> > &input,
3535  const MappingType mapping_type,
3536  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
3537  const ArrayView<Tensor<3,spacedim> > &output) const
3538 {
3539  switch (mapping_type)
3540  {
3541  case mapping_piola_hessian:
3544  transform_hessians(input, mapping_type, mapping_data, output);
3545  return;
3546  default:
3547  Assert(false, ExcNotImplemented());
3548  }
3549 }
3550 
3551 
3552 
3553 namespace
3554 {
3555  // We cannot query a manifold from the faces of a 1D elements (i.e.,
3556  // vertices), which is why we add a specialization for the 3D case here
3557  template <typename Iterator>
3558  bool check_identical_manifolds_of_quads(const Iterator &)
3559  {
3560  Assert(false, ExcNotImplemented());
3561  return true;
3562  }
3563 
3564  bool check_identical_manifolds_of_quads(const Triangulation<3,3>::cell_iterator &cell)
3565  {
3566  for (unsigned int f=0; f<GeometryInfo<3>::faces_per_cell; ++f)
3567  if (&cell->face(f)->get_manifold() != &cell->get_manifold())
3568  return false;
3569  return true;
3570  }
3571 }
3572 
3573 
3574 
3575 template <int dim, int spacedim>
3576 void
3579  std::vector<Point<spacedim> > &a) const
3580 {
3581  // if we only need the midpoint, then ask for it.
3582  if (this->polynomial_degree==2)
3583  {
3584  for (unsigned int line_no=0; line_no<GeometryInfo<dim>::lines_per_cell; ++line_no)
3585  {
3586  const typename Triangulation<dim,spacedim>::line_iterator line =
3587  (dim == 1 ?
3588  static_cast<typename Triangulation<dim,spacedim>::line_iterator>(cell) :
3589  cell->line(line_no));
3590 
3591  const Manifold<dim,spacedim> &manifold =
3592  ( ( line->manifold_id() == numbers::invalid_manifold_id ) &&
3593  ( dim < spacedim )
3594  ?
3595  cell->get_manifold()
3596  :
3597  line->get_manifold() );
3598  a.push_back(manifold.get_new_point_on_line(line));
3599  }
3600  }
3601  else
3602  // otherwise call the more complicated functions and ask for inner points
3603  // from the boundary description
3604  {
3605  std::vector<Point<spacedim> > tmp_points;
3606  // loop over each of the lines, and if it is at the boundary, then first
3607  // get the boundary description and second compute the points on it
3608  for (unsigned int line_no=0; line_no<GeometryInfo<dim>::lines_per_cell; ++line_no)
3609  {
3610  const typename Triangulation<dim,spacedim>::line_iterator
3611  line = (dim == 1
3612  ?
3613  static_cast<typename Triangulation<dim,spacedim>::line_iterator>(cell)
3614  :
3615  cell->line(line_no));
3616 
3617  const Manifold<dim,spacedim> &manifold =
3618  ( ( line->manifold_id() == numbers::invalid_manifold_id ) &&
3619  ( dim < spacedim )
3620  ?
3621  cell->get_manifold() :
3622  line->get_manifold() );
3623 
3624  if (const Boundary<dim,spacedim> *boundary
3625  = dynamic_cast<const Boundary<dim,spacedim> *>(&manifold))
3626  {
3627  tmp_points.resize(this->polynomial_degree-1);
3628  boundary->get_intermediate_points_on_line(line, tmp_points);
3629  if (dim != 3 || cell->line_orientation(line_no))
3630  a.insert (a.end(), tmp_points.begin(), tmp_points.end());
3631  else
3632  a.insert (a.end(), tmp_points.rbegin(), tmp_points.rend());
3633  }
3634  else
3635  {
3636  tmp_points.resize(2);
3637  tmp_points[0] = cell->vertex(GeometryInfo<dim>::line_to_cell_vertices(line_no, 0));
3638  tmp_points[1] = cell->vertex(GeometryInfo<dim>::line_to_cell_vertices(line_no, 1));
3639  manifold.add_new_points(tmp_points,
3640  support_point_weights_perimeter_to_interior[0], a);
3641  }
3642  }
3643  }
3644 }
3645 
3646 
3647 
3648 template <>
3649 void
3652  std::vector<Point<3> > &a) const
3653 {
3654  const unsigned int faces_per_cell = GeometryInfo<3>::faces_per_cell;
3655 
3656  // used if face quad at boundary or entirely in the interior of the domain
3657  std::vector<Point<3> > tmp_points;
3658 
3659  // loop over all faces and collect points on them
3660  for (unsigned int face_no=0; face_no<faces_per_cell; ++face_no)
3661  {
3662  const Triangulation<3>::face_iterator face = cell->face(face_no);
3663 
3664  // select the correct mappings for the present face
3665  const bool face_orientation = cell->face_orientation(face_no),
3666  face_flip = cell->face_flip (face_no),
3667  face_rotation = cell->face_rotation (face_no);
3668 
3669 #ifdef DEBUG
3670  const unsigned int vertices_per_face = GeometryInfo<3>::vertices_per_face,
3671  lines_per_face = GeometryInfo<3>::lines_per_face;
3672 
3673  // some sanity checks up front
3674  for (unsigned int i=0; i<vertices_per_face; ++i)
3675  Assert(face->vertex_index(i)==cell->vertex_index(
3677  face_orientation,
3678  face_flip,
3679  face_rotation)),
3680  ExcInternalError());
3681 
3682  // indices of the lines that bound a face are given by GeometryInfo<3>::
3683  // face_to_cell_lines
3684  for (unsigned int i=0; i<lines_per_face; ++i)
3685  Assert(face->line(i)==cell->line(GeometryInfo<3>::face_to_cell_lines(
3686  face_no, i, face_orientation, face_flip, face_rotation)),
3687  ExcInternalError());
3688 #endif
3689 
3690  // On a quad, we have to check whether the manifold should determine the
3691  // point distribution from all surrounding points (new manifold code) or
3692  // the old-style Boundary code should simply return the intermediate
3693  // points. The second check is to find out whether the Boundary object
3694  // is actually a StraightBoundary (the default flat manifold assigned to
3695  // the triangulation if no manifold is assigned).
3696  const Boundary<3,3> *boundary =
3697  dynamic_cast<const Boundary<3,3> *>(&face->get_manifold());
3698  if (boundary != NULL &&
3699  std::string(typeid(*boundary).name()).find("StraightBoundary") ==
3700  std::string::npos)
3701  {
3702  // ask the boundary/manifold object to return intermediate points on it
3703  tmp_points.resize((polynomial_degree-1)*(polynomial_degree-1));
3704  boundary->get_intermediate_points_on_quad(face, tmp_points);
3705  for (unsigned int i=0; i<tmp_points.size(); ++i)
3706  a.push_back(tmp_points[fe_q->adjust_quad_dof_index_for_face_orientation(i,
3707  face_orientation,
3708  face_flip,
3709  face_rotation)]);
3710  }
3711  else
3712  {
3713  // extract the points surrounding a quad from the points
3714  // already computed. First get the 4 vertices and then the points on
3715  // the four lines
3716  tmp_points.resize(4 + 4*(polynomial_degree-1));
3717  for (unsigned int v=0; v<GeometryInfo<2>::vertices_per_cell; ++v)
3718  tmp_points[v] = a[GeometryInfo<3>::face_to_cell_vertices(face_no,v)];
3719  if (polynomial_degree > 1)
3720  for (unsigned int line=0; line<GeometryInfo<2>::lines_per_cell; ++line)
3721  for (unsigned int i=0; i<polynomial_degree-1; ++i)
3722  tmp_points[4+line*(polynomial_degree-1)+i] =
3724  (polynomial_degree-1)*
3725  GeometryInfo<3>::face_to_cell_lines(face_no,line) + i];
3726  face->get_manifold().add_new_points (tmp_points,
3727  support_point_weights_perimeter_to_interior[1],
3728  a);
3729  }
3730  }
3731 }
3732 
3733 
3734 
3735 template <>
3736 void
3739  std::vector<Point<3> > &a) const
3740 {
3741  if (const Boundary<2,3> *boundary =
3742  dynamic_cast<const Boundary<2,3> *>(&cell->get_manifold()))
3743  {
3744  std::vector<Point<3> > points((polynomial_degree-1)*(polynomial_degree-1));
3745  boundary->get_intermediate_points_on_quad(cell, points);
3746  a.insert(a.end(), points.begin(), points.end());
3747  }
3748  else
3749  {
3750  std::vector<Point<3> > vertices;
3751  for (unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
3752  vertices.push_back(cell->vertex(i));
3753  Table<2,double> weights(Utilities::fixed_power<2>(polynomial_degree-1),
3755  for (unsigned int q=0, q2=0; q2<polynomial_degree-1; ++q2)
3756  for (unsigned int q1=0; q1<polynomial_degree-1; ++q1, ++q)
3757  {
3758  Point<2> point(line_support_points.point(q1+1)[0],
3759  line_support_points.point(q2+1)[0]);
3760  for (unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
3761  weights(q,i) = GeometryInfo<2>::d_linear_shape_function(point, i);
3762  }
3763  // TODO: use all surrounding points once Boundary path is removed
3764  cell->get_manifold().add_new_points(vertices, weights, a);
3765  }
3766 }
3767 
3768 
3769 
3770 template <int dim, int spacedim>
3771 void
3774  std::vector<Point<spacedim> > &) const
3775 {
3776  Assert (false, ExcInternalError());
3777 }
3778 
3779 
3780 
3781 template<int dim, int spacedim>
3782 std::vector<Point<spacedim> >
3785 {
3786  // get the vertices first
3787  std::vector<Point<spacedim> > a;
3788  a.reserve(Utilities::fixed_power<dim>(polynomial_degree+1));
3789  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
3790  a.push_back(cell->vertex(i));
3791 
3792  if (this->polynomial_degree > 1)
3793  {
3794  // check if all entities have the same manifold id which is when we can
3795  // simply ask the manifold for all points
3796  Assert(dim<=3, ExcImpossibleInDim(dim));
3797  bool all_manifold_ids_are_equal = (dim == spacedim);
3798  if (dim > 1)
3799  for (unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
3800  if (&cell->line(l)->get_manifold() != &cell->get_manifold())
3801  all_manifold_ids_are_equal = false;
3802  if (dim == 3)
3803  if (check_identical_manifolds_of_quads(cell) == false)
3804  all_manifold_ids_are_equal = false;
3805  if (all_manifold_ids_are_equal)
3806  {
3807  std::vector<Point<spacedim> > vertices(a);
3808  cell->get_manifold().add_new_points(vertices, support_point_weights_cell, a);
3809  }
3810  else
3811  switch (dim)
3812  {
3813  case 1:
3814  add_line_support_points(cell, a);
3815  break;
3816  case 2:
3817  // in 2d, add the points on the four bounding lines to the exterior
3818  // (outer) points
3819  add_line_support_points(cell, a);
3820 
3821  // then get the interior support points
3822  if (dim != spacedim)
3823  add_quad_support_points(cell, a);
3824  else
3825  {
3826  std::vector<Point<spacedim> > tmp_points(a);
3827  cell->get_manifold().add_new_points(tmp_points,
3828  support_point_weights_perimeter_to_interior[1],
3829  a);
3830  }
3831  break;
3832 
3833  case 3:
3834  // in 3d also add the points located on the boundary faces
3835  add_line_support_points (cell, a);
3836  add_quad_support_points (cell, a);
3837 
3838  // then compute the interior points
3839  {
3840  std::vector<Point<spacedim> > tmp_points(a);
3841  cell->get_manifold().add_new_points(tmp_points,
3842  support_point_weights_perimeter_to_interior[2],
3843  a);
3844  }
3845  break;
3846 
3847  default:
3848  Assert(false, ExcNotImplemented());
3849  break;
3850  }
3851  }
3852 
3853  return a;
3854 }
3855 
3856 
3857 
3858 //--------------------------- Explicit instantiations -----------------------
3859 #include "mapping_q_generic.inst"
3860 
3861 
3862 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
unsigned int n() const
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
Definition: tria.h:67
static const unsigned int invalid_unsigned_int
Definition: types.h:170
const std_cxx11::unique_ptr< FE_Q< dim > > fe_q
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
std::vector< Tensor< 2, dim > > shape_second_derivatives
const unsigned int polynomial_degree
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const
Contravariant transformation.
const std::vector< Point< dim > > & get_points() const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const
const std::vector< double > & get_weights() const
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
virtual InternalData * get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
MappingType
Definition: mapping.h:50
virtual InternalData * get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const
std::vector< Table< 2, double > > support_point_weights_perimeter_to_interior
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
virtual Mapping< dim, spacedim > * clone() const
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const
Volume element.
Outer normal vector, not normalized.
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
Definition: polynomial.cc:829
Table< 2, double > support_point_weights_cell
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
Determinant of the Jacobian.
Triangulation< dim, spacedim >::cell_iterator cell_of_current_support_points
Transformed quadrature points.
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
MappingQGeneric(const unsigned int polynomial_degree)
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
void maybe_compute_face_data(const ::Mapping< dim, spacedim > &mapping, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const std::vector< double > &weights, const typename ::MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data, internal::FEValues::MappingRelatedData< dim, spacedim > &output_data)
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:991
static DataSetDescriptor cell()
Definition: qprojector.h:348
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
Definition: point.h:89
std::vector< Tensor< 1, spacedim > > boundary_forms
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
InternalData(const unsigned int polynomial_degree)
void maybe_update_jacobian_pushed_forward_3rd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data, const FiniteElement< dim, spacedim > &fe, const ComponentMask &fe_mask, const std::vector< unsigned int > &fe_to_real, std::vector< Tensor< 5, spacedim > > &jacobian_pushed_forward_3rd_derivatives)
std::vector< DerivativeForm< 2, dim, spacedim > > jacobian_grads
std::vector< Tensor< 1, dim > > shape_derivatives
virtual void add_quad_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
void compute_shape_function_values(const std::vector< Point< dim > > &unit_points)
virtual void get_intermediate_points_on_quad(const typename Triangulation< dim, spacedim >::quad_iterator &quad, std::vector< Point< spacedim > > &points) const
Definition: tria.h:52
static ::ExceptionBase & ExcMessage(std::string arg1)
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
virtual Point< spacedim > get_new_point_on_line(const typename Triangulation< dim, spacedim >::line_iterator &line) const
Definition: manifold.cc:375
Definition: fe_q.h:552
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
Tensor< 1, spacedim, Number > apply_transformation(const DerivativeForm< 1, dim, spacedim, Number > &DF, const Tensor< 1, dim, Number > &T)
void reinit(const TableIndices< N > &new_size, const bool omit_default_initialization=false)
Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
Definition: tensor.h:1672
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
Abstract base class for mapping classes.
Definition: dof_tools.h:46
std::vector< Tensor< 1, spacedim > > normal_vectors
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const
void loop(ITERATOR begin, typename identity< ITERATOR >::type end, DOFINFO &dinfo, INFOBOX &info, const std_cxx11::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &cell_worker, const std_cxx11::function< void(DOFINFO &, typename INFOBOX::CellInfo &)> &boundary_worker, const std_cxx11::function< void(DOFINFO &, DOFINFO &, typename INFOBOX::CellInfo &, typename INFOBOX::CellInfo &)> &face_worker, ASSEMBLER &assembler, const LoopControl &lctrl=LoopControl())
Definition: loop.h:383
std::vector< Point< spacedim > > mapping_support_points
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
DerivativeForm< 1, spacedim, dim, Number > transpose() const
void do_fill_fe_face_values(const ::Mapping< dim, spacedim > &mapping, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const typename ::QProjector< dim >::DataSetDescriptor data_set, const Quadrature< dim-1 > &quadrature, const typename ::MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data, const FiniteElement< dim, spacedim > &fe, const ComponentMask &fe_mask, const std::vector< unsigned int > &fe_to_real, internal::FEValues::MappingRelatedData< dim, spacedim > &output_data)
void initialize(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
std::vector< double > volume_elements
Gradient of volume element.
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:85
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
const unsigned int n_shape_functions
unsigned int size() const
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Tensor< 1, dim, Number > cross_product_3d(const Tensor< 1, dim, Number > &src1, const Tensor< 1, dim, Number > &src2)
Definition: tensor.h:1698
Point< dim > transform_real_to_unit_cell_internal(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit) const
virtual void add_new_points(const std::vector< Point< spacedim > > &surrounding_points, const Table< 2, double > &weights, std::vector< Point< spacedim > > &new_points) const
Definition: manifold.cc:157
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
virtual std::size_t memory_consumption() const
unsigned int get_degree() const
double norm(const FEValuesBase< dim > &fe, const VectorSlice< const std::vector< std::vector< Tensor< 1, dim > > > > &Du)
Definition: divergence.h:532
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
Number determinant(const SymmetricTensor< 2, dim, Number > &t)
VectorizedArray< Number > sqrt(const ::VectorizedArray< Number > &x)
Definition: mpi.h:48
const types::manifold_id invalid_manifold_id
Definition: types.h:226
void maybe_update_jacobian_3rd_derivatives(const CellSimilarity::Similarity cell_similarity, const typename ::QProjector< dim >::DataSetDescriptor data_set, const typename ::MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data, const FiniteElement< dim, spacedim > &fe, const ComponentMask &fe_mask, const std::vector< unsigned int > &fe_to_real, std::vector< DerivativeForm< 4, dim, spacedim > > &jacobian_3rd_derivatives)
Number determinant(const Tensor< 2, dim, Number > &t)
Definition: tensor.h:1728
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const
Definition: mpi.h:41
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
std::vector< Point< spacedim > > quadrature_points
Normal vectors.
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual void add_line_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
virtual InternalData * get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
static Point< dim > project_to_unit_cell(const Point< dim > &p)
static ::ExceptionBase & ExcNotImplemented()
std::vector< Tensor< 3, dim > > shape_third_derivatives
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
SymmetricTensor< 2, dim, Number > invert(const SymmetricTensor< 2, dim, Number > &)
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
double compute_value(const unsigned int i, const Point< dim > &p) const
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
void clear()
Definition: tensor.h:1072
std::vector< double > shape_values
Point< 3 > point(const gp_Pnt &p)
Definition: utilities.cc:176
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
static double subface_ratio(const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
const Manifold< dim, spacedim > & get_manifold(const types::manifold_id number) const
Definition: tria.cc:9236
Covariant transformation.
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
UpdateFlags update_each
Definition: mapping.h:560
static ::ExceptionBase & ExcInternalError()