Reference documentation for deal.II version 8.5.1
mapping_fe_field.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 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 #include <deal.II/base/utilities.h>
17 #include <deal.II/base/polynomial.h>
18 #include <deal.II/base/quadrature.h>
19 #include <deal.II/base/quadrature_lib.h>
20 #include <deal.II/base/qprojector.h>
21 #include <deal.II/base/memory_consumption.h>
22 #include <deal.II/base/tensor_product_polynomials.h>
23 #include <deal.II/base/std_cxx11/unique_ptr.h>
24 #include <deal.II/lac/full_matrix.h>
25 #include <deal.II/lac/vector.h>
26 #include <deal.II/lac/block_vector.h>
27 #include <deal.II/lac/la_vector.h>
28 #include <deal.II/lac/parallel_vector.h>
29 #include <deal.II/lac/parallel_block_vector.h>
30 #include <deal.II/lac/petsc_vector.h>
31 #include <deal.II/lac/petsc_block_vector.h>
32 #include <deal.II/lac/trilinos_vector.h>
33 #include <deal.II/lac/trilinos_block_vector.h>
34 #include <deal.II/lac/trilinos_parallel_block_vector.h>
35 #include <deal.II/grid/tria_iterator.h>
36 #include <deal.II/grid/tria_boundary.h>
37 #include <deal.II/dofs/dof_accessor.h>
38 #include <deal.II/fe/fe_tools.h>
39 #include <deal.II/fe/fe_values.h>
40 #include <deal.II/fe/fe_system.h>
41 #include <deal.II/fe/mapping_fe_field.h>
42 #include <deal.II/fe/fe_q.h>
43 #include <deal.II/fe/mapping.h>
44 #include <deal.II/fe/mapping_q1.h>
45 #include <deal.II/base/qprojector.h>
46 #include <deal.II/base/thread_management.h>
47 #include <deal.II/numerics/vector_tools.h>
48 
49 #include <numeric>
50 #include <fstream>
51 
52 
53 
54 DEAL_II_NAMESPACE_OPEN
55 
56 
57 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
60  const ComponentMask mask)
61  :
62  n_shape_functions (fe.dofs_per_cell),
63  mask (mask),
64  local_dof_indices(fe.dofs_per_cell),
65  local_dof_values(fe.dofs_per_cell)
66 {}
67 
68 
69 
70 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
71 std::size_t
73 {
74  Assert (false, ExcNotImplemented());
75  return 0;
76 }
77 
78 
79 
80 template<int dim, int spacedim, typename DoFHandlerType, typename VectorType>
81 double &
83 (const unsigned int qpoint,
84  const unsigned int shape_nr)
85 {
86  Assert(qpoint*n_shape_functions + shape_nr < shape_values.size(),
87  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
88  shape_values.size()));
89  return shape_values [qpoint*n_shape_functions + shape_nr];
90 }
91 
92 
93 template<int dim, int spacedim, typename DoFHandlerType, typename VectorType>
94 const Tensor<1,dim> &
96 (const unsigned int qpoint,
97  const unsigned int shape_nr) const
98 {
99  Assert(qpoint*n_shape_functions + shape_nr < shape_derivatives.size(),
100  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
101  shape_derivatives.size()));
102  return shape_derivatives [qpoint*n_shape_functions + shape_nr];
103 }
104 
105 
106 
107 template<int dim, int spacedim, typename DoFHandlerType, typename VectorType>
110 (const unsigned int qpoint,
111  const unsigned int shape_nr)
112 {
113  Assert(qpoint*n_shape_functions + shape_nr < shape_derivatives.size(),
114  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
115  shape_derivatives.size()));
116  return shape_derivatives [qpoint*n_shape_functions + shape_nr];
117 }
118 
119 
120 template <int dim, int spacedim, typename DoFHandlerType, typename VectorType>
121 const Tensor<2,dim> &
123 (const unsigned int qpoint,
124  const unsigned int shape_nr) const
125 {
126  Assert(qpoint*n_shape_functions + shape_nr < shape_second_derivatives.size(),
127  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
128  shape_second_derivatives.size()));
129  return shape_second_derivatives [qpoint*n_shape_functions + shape_nr];
130 }
131 
132 
133 
134 template <int dim, int spacedim, typename DoFHandlerType, typename VectorType>
137 (const unsigned int qpoint,
138  const unsigned int shape_nr)
139 {
140  Assert(qpoint*n_shape_functions + shape_nr < shape_second_derivatives.size(),
141  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
142  shape_second_derivatives.size()));
143  return shape_second_derivatives [qpoint*n_shape_functions + shape_nr];
144 }
145 
146 
147 template <int dim, int spacedim, typename DoFHandlerType, typename VectorType>
148 const Tensor<3,dim> &
150 (const unsigned int qpoint,
151  const unsigned int shape_nr) const
152 {
153  Assert(qpoint*n_shape_functions + shape_nr < shape_third_derivatives.size(),
154  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
155  shape_third_derivatives.size()));
156  return shape_third_derivatives [qpoint*n_shape_functions + shape_nr];
157 }
158 
159 
160 
161 template <int dim, int spacedim, typename DoFHandlerType, typename VectorType>
164 (const unsigned int qpoint,
165  const unsigned int shape_nr)
166 {
167  Assert(qpoint*n_shape_functions + shape_nr < shape_third_derivatives.size(),
168  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
169  shape_third_derivatives.size()));
170  return shape_third_derivatives [qpoint*n_shape_functions + shape_nr];
171 }
172 
173 
174 template <int dim, int spacedim, typename DoFHandlerType, typename VectorType>
175 const Tensor<4,dim> &
177 (const unsigned int qpoint,
178  const unsigned int shape_nr) const
179 {
180  Assert(qpoint*n_shape_functions + shape_nr < shape_fourth_derivatives.size(),
181  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
182  shape_fourth_derivatives.size()));
183  return shape_fourth_derivatives [qpoint*n_shape_functions + shape_nr];
184 }
185 
186 
187 
188 template <int dim, int spacedim, typename DoFHandlerType, typename VectorType>
191 (const unsigned int qpoint,
192  const unsigned int shape_nr)
193 {
194  Assert(qpoint*n_shape_functions + shape_nr < shape_fourth_derivatives.size(),
195  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
196  shape_fourth_derivatives.size()));
197  return shape_fourth_derivatives [qpoint*n_shape_functions + shape_nr];
198 }
199 
200 
201 
202 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
204 (const DoFHandlerType &euler_dof_handler,
205  const VectorType &euler_vector,
206  const ComponentMask mask)
207  :
209  fe(&euler_dof_handler.get_fe()),
211  fe_mask(mask.size() ? mask :
212  ComponentMask(fe->get_nonzero_components(0).size(), true)),
214 {
215  unsigned int size = 0;
216  for (unsigned int i=0; i<fe_mask.size(); ++i)
217  {
218  if (fe_mask[i])
219  fe_to_real[i] = size++;
220  }
221  AssertDimension(size,spacedim);
222 }
223 
224 
225 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
228  :
229  euler_vector(mapping.euler_vector),
230  fe(mapping.fe),
232  fe_mask(mapping.fe_mask),
233  fe_to_real(mapping.fe_to_real)
234 {}
235 
236 
237 
238 template<int dim, int spacedim, typename DoFHandlerType, typename VectorType>
239 inline
240 const double &
242 (const unsigned int qpoint,
243  const unsigned int shape_nr) const
244 {
245  Assert(qpoint*n_shape_functions + shape_nr < shape_values.size(),
246  ExcIndexRange(qpoint*n_shape_functions + shape_nr, 0,
247  shape_values.size()));
248  return shape_values [qpoint*n_shape_functions + shape_nr];
249 }
250 
251 
252 
253 template <int dim, int spacedim, typename DoFHandlerType, typename VectorType>
254 bool
256 {
257  return false;
258 }
259 
260 
261 
262 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
263 void
265 compute_shapes_virtual (const std::vector<Point<dim> > &unit_points,
267 {
268  const unsigned int n_points=unit_points.size();
269 
270  for (unsigned int point=0; point<n_points; ++point)
271  {
272  if (data.shape_values.size()!=0)
273  for (unsigned int i=0; i<data.n_shape_functions; ++i)
274  data.shape(point, i) = fe->shape_value(i, unit_points[point]);
275 
276  if (data.shape_derivatives.size()!=0)
277  for (unsigned int i=0; i<data.n_shape_functions; ++i)
278  data.derivative(point, i) = fe->shape_grad(i, unit_points[point]);
279 
280  if (data.shape_second_derivatives.size()!=0)
281  for (unsigned int i=0; i<data.n_shape_functions; ++i)
282  data.second_derivative(point, i) = fe->shape_grad_grad(i, unit_points[point]);
283 
284  if (data.shape_third_derivatives.size()!=0)
285  for (unsigned int i=0; i<data.n_shape_functions; ++i)
286  data.third_derivative(point, i) = fe->shape_3rd_derivative(i, unit_points[point]);
287 
288  if (data.shape_fourth_derivatives.size()!=0)
289  for (unsigned int i=0; i<data.n_shape_functions; ++i)
290  data.fourth_derivative(point, i) = fe->shape_4th_derivative(i, unit_points[point]);
291  }
292 }
293 
294 
295 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
298 {
299  // add flags if the respective quantities are necessary to compute
300  // what we need. note that some flags appear in both conditions and
301  // in subsequent set operations. this leads to some circular
302  // logic. the only way to treat this is to iterate. since there are
303  // 5 if-clauses in the loop, it will take at most 4 iterations to
304  // converge. do them:
305  UpdateFlags out = in;
306  for (unsigned int i=0; i<5; ++i)
307  {
308  // The following is a little incorrect:
309  // If not applied on a face,
310  // update_boundary_forms does not
311  // make sense. On the other hand,
312  // it is necessary on a
313  // face. Currently,
314  // update_boundary_forms is simply
315  // ignored for the interior of a
316  // cell.
317  if (out & (update_JxW_values
319  out |= update_boundary_forms;
320 
328 
329  if (out & (update_inverse_jacobians
334 
335  // The contravariant transformation
336  // is a Piola transformation, which
337  // requires the determinant of the
338  // Jacobi matrix of the transformation.
339  // Therefore these values have to be
340  // updated for each cell.
342  out |= update_JxW_values;
343 
344  if (out & update_normal_vectors)
345  out |= update_JxW_values;
346  }
347 
348  return out;
349 }
350 
351 
352 
353 
354 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
355 void
357 (const UpdateFlags update_flags,
358  const Quadrature<dim> &q,
359  const unsigned int n_original_q_points,
360  InternalData &data) const
361 {
362  // store the flags in the internal data object so we can access them
363  // in fill_fe_*_values(). use the transitive hull of the required
364  // flags
365  data.update_each = requires_update_flags(update_flags);
366 
367  const unsigned int n_q_points = q.size();
368 
369  // see if we need the (transformation) shape function values
370  // and/or gradients and resize the necessary arrays
371  if (data.update_each & update_quadrature_points)
372  data.shape_values.resize(data.n_shape_functions * n_q_points);
373 
374  if (data.update_each & (update_covariant_transformation
382  data.shape_derivatives.resize(data.n_shape_functions * n_q_points);
383 
384  if (data.update_each & update_covariant_transformation)
385  data.covariant.resize(n_original_q_points);
386 
387  if (data.update_each & update_contravariant_transformation)
388  data.contravariant.resize(n_original_q_points);
389 
390  if (data.update_each & update_volume_elements)
391  data.volume_elements.resize(n_original_q_points);
392 
394  data.shape_second_derivatives.resize(data.n_shape_functions * n_q_points);
395 
397  data.shape_third_derivatives.resize(data.n_shape_functions * n_q_points);
398 
400  data.shape_fourth_derivatives.resize(data.n_shape_functions * n_q_points);
401 
403 }
404 
405 
406 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
407 void
409 (const UpdateFlags update_flags,
410  const Quadrature<dim> &q,
411  const unsigned int n_original_q_points,
412  InternalData &data) const
413 {
414  compute_data (update_flags, q, n_original_q_points, data);
415 
416  if (dim > 1)
417  {
418  if (data.update_each & update_boundary_forms)
419  {
420  data.aux.resize (dim-1, std::vector<Tensor<1,spacedim> > (n_original_q_points));
421 
422  // Compute tangentials to the unit cell.
423  for (unsigned int i=0; i<data.unit_tangentials.size(); ++i)
424  data.unit_tangentials[i].resize (n_original_q_points);
425 
426  if (dim==2)
427  {
428  // ensure a counterclockwise
429  // orientation of tangentials
430  static const int tangential_orientation[4]= {-1,1,1,-1};
431  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
432  {
433  Tensor<1,dim> tang;
434  tang[1-i/2]=tangential_orientation[i];
435  std::fill (data.unit_tangentials[i].begin(),
436  data.unit_tangentials[i].end(),
437  tang);
438  }
439  }
440  else if (dim==3)
441  {
442  for (unsigned int i=0; i<GeometryInfo<dim>::faces_per_cell; ++i)
443  {
444  Tensor<1,dim> tang1, tang2;
445 
446  const unsigned int nd=
448 
449  // first tangential
450  // vector in direction
451  // of the (nd+1)%3 axis
452  // and inverted in case
453  // of unit inward normal
454  tang1[(nd+1)%dim]=GeometryInfo<dim>::unit_normal_orientation[i];
455  // second tangential
456  // vector in direction
457  // of the (nd+2)%3 axis
458  tang2[(nd+2)%dim]=1.;
459 
460  // same unit tangents
461  // for all quadrature
462  // points on this face
463  std::fill (data.unit_tangentials[i].begin(),
464  data.unit_tangentials[i].end(),
465  tang1);
466  std::fill (data.unit_tangentials[GeometryInfo<dim>::faces_per_cell+i].begin(),
467  data.unit_tangentials[GeometryInfo<dim>::faces_per_cell+i].end(),
468  tang2);
469  }
470  }
471  }
472  }
473 }
474 
475 
476 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
477 typename
480  const Quadrature<dim> &quadrature) const
481 {
482  InternalData *data = new InternalData(*fe, fe_mask);
483  this->compute_data (update_flags, quadrature,
484  quadrature.size(), *data);
485  return data;
486 }
487 
488 
489 
490 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
493 (const UpdateFlags update_flags,
494  const Quadrature<dim-1> &quadrature) const
495 {
496  InternalData *data = new InternalData(*fe, fe_mask);
498  this->compute_face_data (update_flags, q,
499  quadrature.size(), *data);
500 
501  return data;
502 }
503 
504 
505 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
508 (const UpdateFlags update_flags,
509  const Quadrature<dim-1> &quadrature) const
510 {
511  InternalData *data = new InternalData(*fe, fe_mask);
513  this->compute_face_data (update_flags, q,
514  quadrature.size(), *data);
515 
516  return data;
517 }
518 
519 
520 
521 namespace internal
522 {
523  namespace
524  {
531  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
532  void
533  maybe_compute_q_points (const typename ::QProjector<dim>::DataSetDescriptor data_set,
534  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
536  const ComponentMask &fe_mask,
537  const std::vector<unsigned int> &fe_to_real,
538  std::vector<Point<spacedim> > &quadrature_points)
539  {
540  const UpdateFlags update_flags = data.update_each;
541 
542  if (update_flags & update_quadrature_points)
543  {
544  for (unsigned int point=0; point<quadrature_points.size(); ++point)
545  {
546  Point<spacedim> result;
547  const double *shape = &data.shape(point+data_set,0);
548 
549  for (unsigned int k=0; k<data.n_shape_functions; ++k)
550  {
551  unsigned int comp_k = fe.system_to_component_index(k).first;
552  if (fe_mask[comp_k])
553  result[fe_to_real[comp_k]] += data.local_dof_values[k] * shape[k];
554  }
555 
556  quadrature_points[point] = result;
557  }
558  }
559  }
560 
568  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
569  void
570  maybe_update_Jacobians (const CellSimilarity::Similarity cell_similarity,
571  const typename ::QProjector<dim>::DataSetDescriptor data_set,
572  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
574  const ComponentMask &fe_mask,
575  const std::vector<unsigned int> &fe_to_real)
576  {
577  const UpdateFlags update_flags = data.update_each;
578 
579  // then Jacobians
580  if (update_flags & update_contravariant_transformation)
581  {
582 
583  // if the current cell is just a translation of the previous one, no
584  // need to recompute jacobians...
585  if (cell_similarity != CellSimilarity::translation)
586  {
587  const unsigned int n_q_points = data.contravariant.size();
588 
589  Assert (data.n_shape_functions > 0, ExcInternalError());
590 
591  for (unsigned int point=0; point<n_q_points; ++point)
592  {
593  const Tensor<1,dim> *data_derv =
594  &data.derivative(point+data_set, 0);
595 
596  Tensor<1, dim> result[spacedim];
597 
598  for (unsigned int k=0; k<data.n_shape_functions; ++k)
599  {
600  unsigned int comp_k = fe.system_to_component_index(k).first;
601  if (fe_mask[comp_k])
602  result[fe_to_real[comp_k]] += data.local_dof_values[k] * data_derv[k];
603  }
604 
605  // write result into contravariant data
606  for (unsigned int i=0; i<spacedim; ++i)
607  {
608  data.contravariant[point][i] = result[i];
609  }
610  }
611  }
612  }
613 
614  if (update_flags & update_covariant_transformation)
615  {
616  AssertDimension(data.covariant.size(), data.contravariant.size());
617  if (cell_similarity != CellSimilarity::translation)
618  for (unsigned int point=0; point<data.contravariant.size(); ++point)
619  data.covariant[point] = (data.contravariant[point]).covariant_form();
620  }
621 
622  if (update_flags & update_volume_elements)
623  {
624  AssertDimension(data.covariant.size(), data.volume_elements.size());
625  if (cell_similarity != CellSimilarity::translation)
626  for (unsigned int point=0; point<data.contravariant.size(); ++point)
627  data.volume_elements[point] = data.contravariant[point].determinant();
628  }
629  }
630 
637  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
638  void
639  maybe_update_jacobian_grads (const CellSimilarity::Similarity cell_similarity,
640  const typename ::QProjector<dim>::DataSetDescriptor data_set,
641  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
643  const ComponentMask &fe_mask,
644  const std::vector<unsigned int> &fe_to_real,
645  std::vector<DerivativeForm<2,dim,spacedim> > &jacobian_grads)
646  {
647  const UpdateFlags update_flags = data.update_each;
648  if (update_flags & update_jacobian_grads)
649  {
650  const unsigned int n_q_points = jacobian_grads.size();
651 
652  if (cell_similarity != CellSimilarity::translation)
653  {
654  for (unsigned int point=0; point<n_q_points; ++point)
655  {
656  const Tensor<2,dim> *second =
657  &data.second_derivative(point+data_set, 0);
658 
660 
661  for (unsigned int k=0; k<data.n_shape_functions; ++k)
662  {
663  unsigned int comp_k = fe.system_to_component_index(k).first;
664  if (fe_mask[comp_k])
665  for (unsigned int j=0; j<dim; ++j)
666  for (unsigned int l=0; l<dim; ++l)
667  result[fe_to_real[comp_k]][j][l] += (second[k][j][l]
668  * data.local_dof_values[k]);
669  }
670 
671  // never touch any data for j=dim in case dim<spacedim, so
672  // it will always be zero as it was initialized
673  for (unsigned int i=0; i<spacedim; ++i)
674  for (unsigned int j=0; j<dim; ++j)
675  for (unsigned int l=0; l<dim; ++l)
676  jacobian_grads[point][i][j][l] = result[i][j][l];
677  }
678  }
679  }
680  }
681 
688  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
689  void
690  maybe_update_jacobian_pushed_forward_grads (
691  const CellSimilarity::Similarity cell_similarity,
692  const typename ::QProjector<dim>::DataSetDescriptor data_set,
693  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
695  const ComponentMask &fe_mask,
696  const std::vector<unsigned int> &fe_to_real,
697  std::vector<Tensor<3,spacedim> > &jacobian_pushed_forward_grads )
698  {
699  const UpdateFlags update_flags = data.update_each;
700  if (update_flags & update_jacobian_pushed_forward_grads)
701  {
702  const unsigned int n_q_points = jacobian_pushed_forward_grads.size();
703 
704  if (cell_similarity != CellSimilarity::translation)
705  {
706  double tmp[spacedim][spacedim][spacedim];
707  for (unsigned int point=0; point<n_q_points; ++point)
708  {
709  const Tensor<2,dim> *second =
710  &data.second_derivative(point+data_set, 0);
711 
713 
714  for (unsigned int k=0; k<data.n_shape_functions; ++k)
715  {
716  unsigned int comp_k = fe.system_to_component_index(k).first;
717  if (fe_mask[comp_k])
718  for (unsigned int j=0; j<dim; ++j)
719  for (unsigned int l=0; l<dim; ++l)
720  result[fe_to_real[comp_k]][j][l] += (second[k][j][l]
721  * data.local_dof_values[k]);
722  }
723 
724  // first push forward the j-components
725  for (unsigned int i=0; i<spacedim; ++i)
726  for (unsigned int j=0; j<spacedim; ++j)
727  for (unsigned int l=0; l<dim; ++l)
728  {
729  tmp[i][j][l] = result[i][0][l] *
730  data.covariant[point][j][0];
731  for (unsigned int jr=1; jr<dim; ++jr)
732  {
733  tmp[i][j][l] += result[i][jr][l] *
734  data.covariant[point][j][jr];
735  }
736  }
737 
738  // now, pushing forward the l-components
739  for (unsigned int i=0; i<spacedim; ++i)
740  for (unsigned int j=0; j<spacedim; ++j)
741  for (unsigned int l=0; l<spacedim; ++l)
742  {
743  jacobian_pushed_forward_grads[point][i][j][l] = tmp[i][j][0] *
744  data.covariant[point][l][0];
745  for (unsigned int lr=1; lr<dim; ++lr)
746  {
747  jacobian_pushed_forward_grads[point][i][j][l] += tmp[i][j][lr] *
748  data.covariant[point][l][lr];
749  }
750 
751  }
752  }
753  }
754  }
755  }
756 
763  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
764  void
765  maybe_update_jacobian_2nd_derivatives (const CellSimilarity::Similarity cell_similarity,
766  const typename ::QProjector<dim>::DataSetDescriptor data_set,
767  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
769  const ComponentMask &fe_mask,
770  const std::vector<unsigned int> &fe_to_real,
771  std::vector<DerivativeForm<3,dim,spacedim> > &jacobian_2nd_derivatives)
772  {
773  const UpdateFlags update_flags = data.update_each;
774  if (update_flags & update_jacobian_2nd_derivatives)
775  {
776  const unsigned int n_q_points = jacobian_2nd_derivatives.size();
777 
778  if (cell_similarity != CellSimilarity::translation)
779  {
780  for (unsigned int point=0; point<n_q_points; ++point)
781  {
782  const Tensor<3,dim> *third =
783  &data.third_derivative(point+data_set, 0);
784 
786 
787  for (unsigned int k=0; k<data.n_shape_functions; ++k)
788  {
789  unsigned int comp_k = fe.system_to_component_index(k).first;
790  if (fe_mask[comp_k])
791  for (unsigned int j=0; j<dim; ++j)
792  for (unsigned int l=0; l<dim; ++l)
793  for (unsigned int m=0; m<dim; ++m)
794  result[fe_to_real[comp_k]][j][l][m] += (third[k][j][l][m]
795  * data.local_dof_values[k]);
796  }
797 
798  // never touch any data for j=dim in case dim<spacedim, so
799  // it will always be zero as it was initialized
800  for (unsigned int i=0; i<spacedim; ++i)
801  for (unsigned int j=0; j<dim; ++j)
802  for (unsigned int l=0; l<dim; ++l)
803  for (unsigned int m=0; m<dim; ++m)
804  jacobian_2nd_derivatives[point][i][j][l][m] = result[i][j][l][m];
805  }
806  }
807  }
808  }
809 
816  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
817  void
818  maybe_update_jacobian_pushed_forward_2nd_derivatives (
819  const CellSimilarity::Similarity cell_similarity,
820  const typename ::QProjector<dim>::DataSetDescriptor data_set,
821  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
823  const ComponentMask &fe_mask,
824  const std::vector<unsigned int> &fe_to_real,
825  std::vector<Tensor<4,spacedim> > &jacobian_pushed_forward_2nd_derivatives )
826  {
827  const UpdateFlags update_flags = data.update_each;
829  {
830  const unsigned int n_q_points = jacobian_pushed_forward_2nd_derivatives.size();
831 
832  if (cell_similarity != CellSimilarity::translation)
833  {
834  double tmp[spacedim][spacedim][spacedim][spacedim];
835  for (unsigned int point=0; point<n_q_points; ++point)
836  {
837  const Tensor<3,dim> *third =
838  &data.third_derivative(point+data_set, 0);
839 
841 
842  for (unsigned int k=0; k<data.n_shape_functions; ++k)
843  {
844  unsigned int comp_k = fe.system_to_component_index(k).first;
845  if (fe_mask[comp_k])
846  for (unsigned int j=0; j<dim; ++j)
847  for (unsigned int l=0; l<dim; ++l)
848  for (unsigned int m=0; m<dim; ++m)
849  result[fe_to_real[comp_k]][j][l][m] += (third[k][j][l][m]
850  * data.local_dof_values[k]);
851  }
852 
853  // push forward the j-coordinate
854  for (unsigned int i=0; i<spacedim; ++i)
855  for (unsigned int j=0; j<spacedim; ++j)
856  for (unsigned int l=0; l<dim; ++l)
857  for (unsigned int m=0; m<dim; ++m)
858  {
859  jacobian_pushed_forward_2nd_derivatives[point][i][j][l][m]
860  = result[i][0][l][m]*
861  data.covariant[point][j][0];
862  for (unsigned int jr=1; jr<dim; ++jr)
863  jacobian_pushed_forward_2nd_derivatives[point][i][j][l][m]
864  += result[i][jr][l][m]*
865  data.covariant[point][j][jr];
866  }
867 
868  // push forward the l-coordinate
869  for (unsigned int i=0; i<spacedim; ++i)
870  for (unsigned int j=0; j<spacedim; ++j)
871  for (unsigned int l=0; l<spacedim; ++l)
872  for (unsigned int m=0; m<dim; ++m)
873  {
874  tmp[i][j][l][m]
875  = jacobian_pushed_forward_2nd_derivatives[point][i][j][0][m]*
876  data.covariant[point][l][0];
877  for (unsigned int lr=1; lr<dim; ++lr)
878  tmp[i][j][l][m]
879  += jacobian_pushed_forward_2nd_derivatives[point][i][j][lr][m]*
880  data.covariant[point][l][lr];
881  }
882 
883  // push forward the m-coordinate
884  for (unsigned int i=0; i<spacedim; ++i)
885  for (unsigned int j=0; j<spacedim; ++j)
886  for (unsigned int l=0; l<spacedim; ++l)
887  for (unsigned int m=0; m<spacedim; ++m)
888  {
889  jacobian_pushed_forward_2nd_derivatives[point][i][j][l][m]
890  = tmp[i][j][l][0]*
891  data.covariant[point][m][0];
892  for (unsigned int mr=1; mr<dim; ++mr)
893  jacobian_pushed_forward_2nd_derivatives[point][i][j][l][m]
894  += tmp[i][j][l][mr]*
895  data.covariant[point][m][mr];
896  }
897  }
898  }
899  }
900  }
901  }
902 
909  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
910  void
912  const typename ::QProjector<dim>::DataSetDescriptor data_set,
913  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
915  const ComponentMask &fe_mask,
916  const std::vector<unsigned int> &fe_to_real,
917  std::vector<DerivativeForm<4,dim,spacedim> > &jacobian_3rd_derivatives)
918  {
919  const UpdateFlags update_flags = data.update_each;
920  if (update_flags & update_jacobian_3rd_derivatives)
921  {
922  const unsigned int n_q_points = jacobian_3rd_derivatives.size();
923 
924  if (cell_similarity != CellSimilarity::translation)
925  {
926  for (unsigned int point=0; point<n_q_points; ++point)
927  {
928  const Tensor<4,dim> *fourth =
929  &data.fourth_derivative(point+data_set, 0);
930 
932 
933  for (unsigned int k=0; k<data.n_shape_functions; ++k)
934  {
935  unsigned int comp_k = fe.system_to_component_index(k).first;
936  if (fe_mask[comp_k])
937  for (unsigned int j=0; j<dim; ++j)
938  for (unsigned int l=0; l<dim; ++l)
939  for (unsigned int m=0; m<dim; ++m)
940  for (unsigned int n=0; n<dim; ++n)
941  result[fe_to_real[comp_k]][j][l][m][n] += (fourth[k][j][l][m][n]
942  * data.local_dof_values[k]);
943  }
944 
945  // never touch any data for j,l,m,n=dim in case dim<spacedim, so
946  // it will always be zero as it was initialized
947  for (unsigned int i=0; i<spacedim; ++i)
948  for (unsigned int j=0; j<dim; ++j)
949  for (unsigned int l=0; l<dim; ++l)
950  for (unsigned int m=0; m<dim; ++m)
951  for (unsigned int n=0; n<dim; ++n)
952  jacobian_3rd_derivatives[point][i][j][l][m][n] = result[i][j][l][m][n];
953  }
954  }
955  }
956  }
957 
965  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
966  void
968  const CellSimilarity::Similarity cell_similarity,
969  const typename ::QProjector<dim>::DataSetDescriptor data_set,
970  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
972  const ComponentMask &fe_mask,
973  const std::vector<unsigned int> &fe_to_real,
974  std::vector<Tensor<5,spacedim> > &jacobian_pushed_forward_3rd_derivatives )
975  {
976  const UpdateFlags update_flags = data.update_each;
978  {
979  const unsigned int n_q_points = jacobian_pushed_forward_3rd_derivatives.size();
980 
981  if (cell_similarity != CellSimilarity::translation)
982  {
983  double tmp[spacedim][spacedim][spacedim][spacedim][spacedim];
984  for (unsigned int point=0; point<n_q_points; ++point)
985  {
986  const Tensor<4,dim> *fourth =
987  &data.fourth_derivative(point+data_set, 0);
988 
990 
991  for (unsigned int k=0; k<data.n_shape_functions; ++k)
992  {
993  unsigned int comp_k = fe.system_to_component_index(k).first;
994  if (fe_mask[comp_k])
995  for (unsigned int j=0; j<dim; ++j)
996  for (unsigned int l=0; l<dim; ++l)
997  for (unsigned int m=0; m<dim; ++m)
998  for (unsigned int n=0; n<dim; ++n)
999  result[fe_to_real[comp_k]][j][l][m][n]
1000  += (fourth[k][j][l][m][n]
1001  * data.local_dof_values[k]);
1002  }
1003 
1004  // push-forward the j-coordinate
1005  for (unsigned int i=0; i<spacedim; ++i)
1006  for (unsigned int j=0; j<spacedim; ++j)
1007  for (unsigned int l=0; l<dim; ++l)
1008  for (unsigned int m=0; m<dim; ++m)
1009  for (unsigned int n=0; n<dim; ++n)
1010  {
1011  tmp[i][j][l][m][n] = result[i][0][l][m][n] *
1012  data.covariant[point][j][0];
1013  for (unsigned int jr=1; jr<dim; ++jr)
1014  tmp[i][j][l][m][n] += result[i][jr][l][m][n] *
1015  data.covariant[point][j][jr];
1016  }
1017 
1018  // push-forward the l-coordinate
1019  for (unsigned int i=0; i<spacedim; ++i)
1020  for (unsigned int j=0; j<spacedim; ++j)
1021  for (unsigned int l=0; l<spacedim; ++l)
1022  for (unsigned int m=0; m<dim; ++m)
1023  for (unsigned int n=0; n<dim; ++n)
1024  {
1025  jacobian_pushed_forward_3rd_derivatives[point][i][j][l][m][n]
1026  = tmp[i][j][0][m][n] *
1027  data.covariant[point][l][0];
1028  for (unsigned int lr=1; lr<dim; ++lr)
1029  jacobian_pushed_forward_3rd_derivatives[point][i][j][l][m][n]
1030  += tmp[i][j][lr][m][n] *
1031  data.covariant[point][l][lr];
1032  }
1033 
1034  // push-forward the m-coordinate
1035  for (unsigned int i=0; i<spacedim; ++i)
1036  for (unsigned int j=0; j<spacedim; ++j)
1037  for (unsigned int l=0; l<spacedim; ++l)
1038  for (unsigned int m=0; m<spacedim; ++m)
1039  for (unsigned int n=0; n<dim; ++n)
1040  {
1041  tmp[i][j][l][m][n]
1042  = jacobian_pushed_forward_3rd_derivatives[point][i][j][l][0][n] *
1043  data.covariant[point][m][0];
1044  for (unsigned int mr=1; mr<dim; ++mr)
1045  tmp[i][j][l][m][n]
1046  += jacobian_pushed_forward_3rd_derivatives[point][i][j][l][mr][n] *
1047  data.covariant[point][m][mr];
1048  }
1049 
1050  // push-forward the n-coordinate
1051  for (unsigned int i=0; i<spacedim; ++i)
1052  for (unsigned int j=0; j<spacedim; ++j)
1053  for (unsigned int l=0; l<spacedim; ++l)
1054  for (unsigned int m=0; m<spacedim; ++m)
1055  for (unsigned int n=0; n<spacedim; ++n)
1056  {
1057  jacobian_pushed_forward_3rd_derivatives[point][i][j][l][m][n]
1058  = tmp[i][j][l][m][0] *
1059  data.covariant[point][n][0];
1060  for (unsigned int nr=1; nr<dim; ++nr)
1061  jacobian_pushed_forward_3rd_derivatives[point][i][j][l][m][n]
1062  += tmp[i][j][l][m][nr] *
1063  data.covariant[point][n][nr];
1064  }
1065  }
1066  }
1067  }
1068  }
1069 
1070 
1080  template <int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1081  void
1082  maybe_compute_face_data (const ::Mapping<dim,spacedim> &mapping,
1083  const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
1084  const unsigned int face_no,
1085  const unsigned int subface_no,
1086  const std::vector<double> &weights,
1087  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
1089  {
1090  const UpdateFlags update_flags = data.update_each;
1091 
1092  if (update_flags & update_boundary_forms)
1093  {
1094  const unsigned int n_q_points = output_data.boundary_forms.size();
1095  if (update_flags & update_normal_vectors)
1096  AssertDimension (output_data.normal_vectors.size(), n_q_points);
1097  if (update_flags & update_JxW_values)
1098  AssertDimension (output_data.JxW_values.size(), n_q_points);
1099 
1100  // map the unit tangentials to the real cell. checking for d!=dim-1
1101  // eliminates compiler warnings regarding unsigned int expressions <
1102  // 0.
1103  for (unsigned int d=0; d!=dim-1; ++d)
1104  {
1106  data.unit_tangentials.size(),
1107  ExcInternalError());
1108  Assert (data.aux[d].size() <=
1109  data.unit_tangentials[face_no+GeometryInfo<dim>::faces_per_cell*d].size(),
1110  ExcInternalError());
1111 
1112  mapping.transform (make_array_view(data.unit_tangentials[face_no+GeometryInfo<dim>::faces_per_cell*d]),
1114  data,
1115  make_array_view(data.aux[d]));
1116  }
1117 
1118  // if dim==spacedim, we can use the unit tangentials to compute the
1119  // boundary form by simply taking the cross product
1120  if (dim == spacedim)
1121  {
1122  for (unsigned int i=0; i<n_q_points; ++i)
1123  switch (dim)
1124  {
1125  case 1:
1126  // in 1d, we don't have access to any of the data.aux
1127  // fields (because it has only dim-1 components), but we
1128  // can still compute the boundary form by simply looking
1129  // at the number of the face
1130  output_data.boundary_forms[i][0] = (face_no == 0 ?
1131  -1 : +1);
1132  break;
1133  case 2:
1134  output_data.boundary_forms[i] = cross_product_2d(data.aux[0][i]);
1135  break;
1136  case 3:
1137  output_data.boundary_forms[i] =
1138  cross_product_3d(data.aux[0][i], data.aux[1][i]);
1139  break;
1140  default:
1141  Assert(false, ExcNotImplemented());
1142  }
1143  }
1144  else //(dim < spacedim)
1145  {
1146  // in the codim-one case, the boundary form results from the
1147  // cross product of all the face tangential vectors and the cell
1148  // normal vector
1149  //
1150  // to compute the cell normal, use the same method used in
1151  // fill_fe_values for cells above
1152  AssertDimension (data.contravariant.size(), n_q_points);
1153 
1154  for (unsigned int point=0; point<n_q_points; ++point)
1155  {
1156  if (dim==1)
1157  {
1158  // J is a tangent vector
1159  output_data.boundary_forms[point] = data.contravariant[point].transpose()[0];
1160  output_data.boundary_forms[point] /=
1161  (face_no == 0 ? -1. : +1.) * output_data.boundary_forms[point].norm();
1162 
1163  }
1164 
1165  if (dim==2)
1166  {
1167  const DerivativeForm<1,spacedim,dim> DX_t =
1168  data.contravariant[point].transpose();
1169 
1170  Tensor<1, spacedim> cell_normal =
1171  cross_product_3d(DX_t[0], DX_t[1]);
1172  cell_normal /= cell_normal.norm();
1173 
1174  // then compute the face normal from the face tangent
1175  // and the cell normal:
1176  output_data.boundary_forms[point] =
1177  cross_product_3d(data.aux[0][point], cell_normal);
1178  }
1179 
1180  }
1181  }
1182 
1183  if (update_flags & (update_normal_vectors | update_JxW_values))
1184  for (unsigned int i=0; i<output_data.boundary_forms.size(); ++i)
1185  {
1186  if (update_flags & update_JxW_values)
1187  {
1188  output_data.JxW_values[i] = output_data.boundary_forms[i].norm() * weights[i];
1189 
1190  if (subface_no != numbers::invalid_unsigned_int)
1191  {
1192  const double area_ratio=GeometryInfo<dim>::subface_ratio(
1193  cell->subface_case(face_no), subface_no);
1194  output_data.JxW_values[i] *= area_ratio;
1195  }
1196  }
1197 
1198  if (update_flags & update_normal_vectors)
1199  output_data.normal_vectors[i] = Point<spacedim>(output_data.boundary_forms[i] / output_data.boundary_forms[i].norm());
1200  }
1201 
1202  if (update_flags & update_jacobians)
1203  for (unsigned int point=0; point<n_q_points; ++point)
1204  output_data.jacobians[point] = data.contravariant[point];
1205 
1206  if (update_flags & update_inverse_jacobians)
1207  for (unsigned int point=0; point<n_q_points; ++point)
1208  output_data.inverse_jacobians[point] = data.covariant[point].transpose();
1209  }
1210  }
1211 
1218  template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1219  void
1220  do_fill_fe_face_values (const ::Mapping<dim,spacedim> &mapping,
1221  const typename ::Triangulation<dim,spacedim>::cell_iterator &cell,
1222  const unsigned int face_no,
1223  const unsigned int subface_no,
1224  const typename ::QProjector<dim>::DataSetDescriptor data_set,
1225  const Quadrature<dim-1> &quadrature,
1226  const typename ::MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &data,
1228  const ComponentMask &fe_mask,
1229  const std::vector<unsigned int> &fe_to_real,
1231  {
1232  maybe_compute_q_points<dim,spacedim,VectorType,DoFHandlerType>
1233  (data_set,
1234  data,
1235  fe, fe_mask, fe_to_real,
1236  output_data.quadrature_points);
1237 
1238  maybe_update_Jacobians<dim,spacedim,VectorType,DoFHandlerType>
1240  data_set,
1241  data,
1242  fe, fe_mask, fe_to_real);
1243 
1244  maybe_update_jacobian_grads<dim,spacedim,VectorType,DoFHandlerType>
1246  data_set,
1247  data,
1248  fe, fe_mask, fe_to_real,
1249  output_data.jacobian_grads);
1250 
1251  maybe_update_jacobian_pushed_forward_grads<dim,spacedim,VectorType,DoFHandlerType>
1253  data_set,
1254  data,
1255  fe, fe_mask, fe_to_real,
1256  output_data.jacobian_pushed_forward_grads);
1257 
1258  maybe_update_jacobian_2nd_derivatives<dim,spacedim,VectorType,DoFHandlerType>
1260  data_set,
1261  data,
1262  fe, fe_mask, fe_to_real,
1263  output_data.jacobian_2nd_derivatives);
1264 
1265  maybe_update_jacobian_pushed_forward_2nd_derivatives<dim,spacedim,VectorType,DoFHandlerType>
1267  data_set,
1268  data,
1269  fe, fe_mask, fe_to_real,
1271 
1272  maybe_update_jacobian_3rd_derivatives<dim,spacedim,VectorType,DoFHandlerType>
1274  data_set,
1275  data,
1276  fe, fe_mask, fe_to_real,
1277  output_data.jacobian_3rd_derivatives);
1278 
1279  maybe_update_jacobian_pushed_forward_3rd_derivatives<dim,spacedim,VectorType,DoFHandlerType>
1281  data_set,
1282  data,
1283  fe, fe_mask, fe_to_real,
1285 
1286  maybe_compute_face_data<dim,spacedim,VectorType,DoFHandlerType>
1287  (mapping,
1288  cell, face_no, subface_no,
1289  quadrature.get_weights(), data,
1290  output_data);
1291  }
1292 }
1293 
1294 
1295 // Note that the CellSimilarity flag is modifiable, since MappingFEField can need to
1296 // recalculate data even when cells are similar.
1297 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1301  const CellSimilarity::Similarity cell_similarity,
1302  const Quadrature<dim> &quadrature,
1303  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
1305 {
1306  // convert data object to internal data for this class. fails with an
1307  // exception if that is not possible
1308  Assert (dynamic_cast<const InternalData *> (&internal_data) != 0, ExcInternalError());
1309  const InternalData &data = static_cast<const InternalData &> (internal_data);
1310 
1311  const unsigned int n_q_points=quadrature.size();
1312  const CellSimilarity::Similarity updated_cell_similarity
1313  = (get_degree() == 1
1314  ?
1315  cell_similarity
1316  :
1318 
1319  update_internal_dofs(cell, data);
1320 
1321  internal::maybe_compute_q_points<dim,spacedim,VectorType,DoFHandlerType>
1323  data, *fe, fe_mask, fe_to_real,
1324  output_data.quadrature_points);
1325 
1326  internal::maybe_update_Jacobians<dim,spacedim,VectorType,DoFHandlerType>
1327  (cell_similarity,
1329  data, *fe, fe_mask, fe_to_real);
1330 
1331  const UpdateFlags update_flags = data.update_each;
1332  const std::vector<double> &weights=quadrature.get_weights();
1333 
1334  // Multiply quadrature weights by absolute value of Jacobian determinants or
1335  // the area element g=sqrt(DX^t DX) in case of codim > 0
1336 
1337  if (update_flags & (update_normal_vectors | update_JxW_values))
1338  {
1339  AssertDimension (output_data.JxW_values.size(), n_q_points);
1340 
1341  Assert( !(update_flags & update_normal_vectors ) ||
1342  (output_data.normal_vectors.size() == n_q_points),
1343  ExcDimensionMismatch(output_data.normal_vectors.size(), n_q_points));
1344 
1345 
1346  if (cell_similarity != CellSimilarity::translation)
1347  for (unsigned int point=0; point<n_q_points; ++point)
1348  {
1349  if (dim == spacedim)
1350  {
1351  const double det = data.contravariant[point].determinant();
1352 
1353  // check for distorted cells.
1354 
1355  // TODO: this allows for anisotropies of up to 1e6 in 3D and
1356  // 1e12 in 2D. might want to find a finer
1357  // (dimension-independent) criterion
1358  Assert (det > 1e-12*Utilities::fixed_power<dim>(cell->diameter()/
1359  std::sqrt(double(dim))),
1360  (typename Mapping<dim,spacedim>::ExcDistortedMappedCell(cell->center(), det, point)));
1361  output_data.JxW_values[point] = weights[point] * det;
1362  }
1363  // if dim==spacedim, then there is no cell normal to
1364  // compute. since this is for FEValues (and not FEFaceValues),
1365  // there are also no face normals to compute
1366  else //codim>0 case
1367  {
1368  Tensor<1, spacedim> DX_t [dim];
1369  for (unsigned int i=0; i<spacedim; ++i)
1370  for (unsigned int j=0; j<dim; ++j)
1371  DX_t[j][i] = data.contravariant[point][i][j];
1372 
1373  Tensor<2, dim> G; //First fundamental form
1374  for (unsigned int i=0; i<dim; ++i)
1375  for (unsigned int j=0; j<dim; ++j)
1376  G[i][j] = DX_t[i] * DX_t[j];
1377 
1378  output_data.JxW_values[point] = sqrt(determinant(G)) * weights[point];
1379 
1380  if (cell_similarity == CellSimilarity::inverted_translation)
1381  {
1382  // we only need to flip the normal
1383  if (update_flags & update_normal_vectors)
1384  output_data.normal_vectors[point] *= -1.;
1385  }
1386  else
1387  {
1388  if (update_flags & update_normal_vectors)
1389  {
1390  Assert (spacedim - dim == 1,
1391  ExcMessage("There is no cell normal in codim 2."));
1392 
1393  if (dim==1)
1394  output_data.normal_vectors[point] =
1395  cross_product_2d(-DX_t[0]);
1396  else //dim == 2
1397  output_data.normal_vectors[point] =
1398  cross_product_3d(DX_t[0], DX_t[1]);
1399 
1400  output_data.normal_vectors[point] /= output_data.normal_vectors[point].norm();
1401 
1402  if (cell->direction_flag() == false)
1403  output_data.normal_vectors[point] *= -1.;
1404  }
1405 
1406  }
1407  } //codim>0 case
1408  }
1409  }
1410 
1411  // copy values from InternalData to vector given by reference
1412  if (update_flags & update_jacobians)
1413  {
1414  AssertDimension (output_data.jacobians.size(), n_q_points);
1415  if (cell_similarity != CellSimilarity::translation)
1416  for (unsigned int point=0; point<n_q_points; ++point)
1417  output_data.jacobians[point] = data.contravariant[point];
1418  }
1419 
1420  // copy values from InternalData to vector given by reference
1421  if (update_flags & update_inverse_jacobians)
1422  {
1423  AssertDimension (output_data.inverse_jacobians.size(), n_q_points);
1424  if (cell_similarity != CellSimilarity::translation)
1425  for (unsigned int point=0; point<n_q_points; ++point)
1426  output_data.inverse_jacobians[point] = data.covariant[point].transpose();
1427  }
1428 
1429  // calculate derivatives of the Jacobians
1430  internal::maybe_update_jacobian_grads<dim,spacedim,VectorType,DoFHandlerType>
1431  (cell_similarity,
1433  data, *fe, fe_mask, fe_to_real,
1434  output_data.jacobian_grads);
1435 
1436  // calculate derivatives of the Jacobians pushed forward to real cell coordinates
1437  internal::maybe_update_jacobian_pushed_forward_grads<dim,spacedim,VectorType,DoFHandlerType>
1438  (cell_similarity,
1440  data, *fe, fe_mask, fe_to_real,
1441  output_data.jacobian_pushed_forward_grads);
1442 
1443  // calculate hessians of the Jacobians
1444  internal::maybe_update_jacobian_2nd_derivatives<dim,spacedim,VectorType,DoFHandlerType>
1445  (cell_similarity,
1447  data, *fe, fe_mask, fe_to_real,
1448  output_data.jacobian_2nd_derivatives);
1449 
1450  // calculate hessians of the Jacobians pushed forward to real cell coordinates
1451  internal::maybe_update_jacobian_pushed_forward_2nd_derivatives<dim,spacedim,VectorType,DoFHandlerType>
1452  (cell_similarity,
1454  data, *fe, fe_mask, fe_to_real,
1456 
1457  // calculate gradients of the hessians of the Jacobians
1458  internal::maybe_update_jacobian_3rd_derivatives<dim,spacedim,VectorType,DoFHandlerType>
1459  (cell_similarity,
1461  data, *fe, fe_mask, fe_to_real,
1462  output_data.jacobian_3rd_derivatives);
1463 
1464  // calculate gradients of the hessians of the Jacobians pushed forward to real
1465  // cell coordinates
1466  internal::maybe_update_jacobian_pushed_forward_3rd_derivatives<dim,spacedim,VectorType,DoFHandlerType>
1467  (cell_similarity,
1469  data, *fe, fe_mask, fe_to_real,
1471 
1472  return updated_cell_similarity;
1473 }
1474 
1475 
1476 
1477 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1478 void
1481  const unsigned int face_no,
1482  const Quadrature<dim-1> &quadrature,
1483  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
1485 {
1486  // convert data object to internal data for this class. fails with an
1487  // exception if that is not possible
1488  Assert (dynamic_cast<const InternalData *> (&internal_data) != 0,
1489  ExcInternalError());
1490  const InternalData &data = static_cast<const InternalData &> (internal_data);
1491 
1492  update_internal_dofs(cell, data);
1493 
1494  internal::do_fill_fe_face_values<dim,spacedim,VectorType,DoFHandlerType>
1495  (*this,
1496  cell, face_no, numbers::invalid_unsigned_int,
1498  face (face_no,
1499  cell->face_orientation(face_no),
1500  cell->face_flip(face_no),
1501  cell->face_rotation(face_no),
1502  quadrature.size()),
1503  quadrature,
1504  data,
1505  *fe, fe_mask, fe_to_real,
1506  output_data);
1507 }
1508 
1509 
1510 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1511 void
1514  const unsigned int face_no,
1515  const unsigned int subface_no,
1516  const Quadrature<dim-1> &quadrature,
1517  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
1519 {
1520  // convert data object to internal data for this class. fails with an
1521  // exception if that is not possible
1522  Assert (dynamic_cast<const InternalData *> (&internal_data) != 0,
1523  ExcInternalError());
1524  const InternalData &data = static_cast<const InternalData &> (internal_data);
1525 
1526  update_internal_dofs(cell, data);
1527 
1528  internal::do_fill_fe_face_values<dim,spacedim,VectorType,DoFHandlerType>
1529  (*this,
1530  cell, face_no, numbers::invalid_unsigned_int,
1532  subface (face_no, subface_no,
1533  cell->face_orientation(face_no),
1534  cell->face_flip(face_no),
1535  cell->face_rotation(face_no),
1536  quadrature.size(),
1537  cell->subface_case(face_no)),
1538  quadrature,
1539  data,
1540  *fe, fe_mask, fe_to_real,
1541  output_data);
1542 }
1543 
1544 
1545 namespace
1546 {
1547  template<int dim, int spacedim, int rank, typename VectorType, typename DoFHandlerType>
1548  void
1549  transform_fields(const ArrayView<const Tensor<rank,dim> > &input,
1550  const MappingType mapping_type,
1551  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
1552  const ArrayView<Tensor<rank,spacedim> > &output)
1553  {
1554  AssertDimension (input.size(), output.size());
1555  Assert ((dynamic_cast<const typename MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData *>(&mapping_data) != 0),
1556  ExcInternalError());
1558  &data = static_cast<const typename MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &>(mapping_data);
1559 
1560  switch (mapping_type)
1561  {
1562  case mapping_contravariant:
1563  {
1565  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
1566 
1567  for (unsigned int i=0; i<output.size(); ++i)
1568  output[i] = apply_transformation(data.contravariant[i], input[i]);
1569 
1570  return;
1571  }
1572 
1573  case mapping_piola:
1574  {
1576  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
1578  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
1579  Assert (rank==1, ExcMessage("Only for rank 1"));
1580  for (unsigned int i=0; i<output.size(); ++i)
1581  {
1582  output[i] = apply_transformation(data.contravariant[i], input[i]);
1583  output[i] /= data.volume_elements[i];
1584  }
1585  return;
1586  }
1587 
1588 
1589  //We still allow this operation as in the
1590  //reference cell Derivatives are Tensor
1591  //rather than DerivativeForm
1592  case mapping_covariant:
1593  {
1595  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
1596 
1597  for (unsigned int i=0; i<output.size(); ++i)
1598  output[i] = apply_transformation(data.covariant[i], input[i]);
1599 
1600  return;
1601  }
1602 
1603  default:
1604  Assert(false, ExcNotImplemented());
1605  }
1606  }
1607 
1608 
1609  template<int dim, int spacedim, int rank, typename VectorType, typename DoFHandlerType>
1610  void
1611  transform_differential_forms
1612  (const ArrayView<const DerivativeForm<rank, dim,spacedim> > &input,
1613  const MappingType mapping_type,
1614  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
1615  const ArrayView<Tensor<rank+1, spacedim> > &output)
1616  {
1617 
1618  AssertDimension (input.size(), output.size());
1619  Assert ((dynamic_cast<const typename MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData *>(&mapping_data) != 0),
1620  ExcInternalError());
1622  &data = static_cast<const typename MappingFEField<dim,spacedim,VectorType,DoFHandlerType>::InternalData &>(mapping_data);
1623 
1624  switch (mapping_type)
1625  {
1626  case mapping_covariant:
1627  {
1629  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
1630 
1631  for (unsigned int i=0; i<output.size(); ++i)
1632  output[i] = apply_transformation(data.covariant[i], input[i]);
1633 
1634  return;
1635  }
1636  default:
1637  Assert(false, ExcNotImplemented());
1638  }
1639 
1640  }
1641 }
1642 
1643 
1644 
1645 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1646 void
1648 transform (const ArrayView<const Tensor<1,dim> > &input,
1649  const MappingType mapping_type,
1650  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
1651  const ArrayView<Tensor<1,spacedim> > &output) const
1652 {
1653  AssertDimension (input.size(), output.size());
1654 
1655  transform_fields<dim,spacedim,1,VectorType,DoFHandlerType>(input, mapping_type, mapping_data, output);
1656 }
1657 
1658 
1659 
1660 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1661 void
1664  const MappingType mapping_type,
1665  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
1666  const ArrayView<Tensor<2,spacedim> > &output) const
1667 {
1668  AssertDimension (input.size(), output.size());
1669 
1670  transform_differential_forms<dim,spacedim,1,VectorType,DoFHandlerType>(input, mapping_type, mapping_data, output);
1671 }
1672 
1673 
1674 
1675 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1676 void
1678 transform (const ArrayView<const Tensor<2, dim> > &input,
1679  const MappingType ,
1680  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
1681  const ArrayView<Tensor<2,spacedim> > &output) const
1682 {
1683  (void)input;
1684  (void)output;
1685  (void)mapping_data;
1686  AssertDimension (input.size(), output.size());
1687 
1688  AssertThrow(false, ExcNotImplemented());
1689 }
1690 
1691 
1692 
1693 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1694 void
1697  const MappingType mapping_type,
1698  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
1699  const ArrayView<Tensor<3,spacedim> > &output) const
1700 {
1701  AssertDimension (input.size(), output.size());
1702  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
1703  ExcInternalError());
1704  const InternalData &data = static_cast<const InternalData &>(mapping_data);
1705 
1706  switch (mapping_type)
1707  {
1709  {
1711  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
1712 
1713  for (unsigned int q=0; q<output.size(); ++q)
1714  for (unsigned int i=0; i<spacedim; ++i)
1715  for (unsigned int j=0; j<spacedim; ++j)
1716  for (unsigned int k=0; k<spacedim; ++k)
1717  {
1718  output[q][i][j][k] = data.covariant[q][j][0]
1719  * data.covariant[q][k][0]
1720  * input[q][i][0][0];
1721  for (unsigned int J=0; J<dim; ++J)
1722  {
1723  const unsigned int K0 = (0==J)? 1 : 0;
1724  for (unsigned int K=K0; K<dim; ++K)
1725  output[q][i][j][k] += data.covariant[q][j][J]
1726  * data.covariant[q][k][K]
1727  * input[q][i][J][K];
1728  }
1729 
1730  }
1731  return;
1732  }
1733 
1734  default:
1735  Assert(false, ExcNotImplemented());
1736  }
1737 
1738 }
1739 
1740 
1741 
1742 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1743 void
1745 transform (const ArrayView<const Tensor<3,dim> > &input,
1746  const MappingType /*mapping_type*/,
1747  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
1748  const ArrayView<Tensor<3,spacedim> > &output) const
1749 {
1750 
1751  (void)input;
1752  (void)output;
1753  (void)mapping_data;
1754  AssertDimension (input.size(), output.size());
1755 
1756  AssertThrow(false, ExcNotImplemented());
1757 
1758 }
1759 
1760 
1761 
1762 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1766  const Point<dim> &p) const
1767 {
1768 // Use the get_data function to create an InternalData with data vectors of
1769 // the right size and transformation shape values already computed at point
1770 // p.
1771  const Quadrature<dim> point_quadrature(p);
1772  std_cxx11::unique_ptr<InternalData> mdata (get_data(update_quadrature_points | update_jacobians,
1773  point_quadrature));
1774 
1775  update_internal_dofs(cell, *mdata);
1776 
1777  return do_transform_unit_to_real_cell(*mdata);
1778 }
1779 
1780 
1781 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1785 {
1786  Point<spacedim> p_real;
1787 
1788  for (unsigned int i=0; i<data.n_shape_functions; ++i)
1789  {
1790  unsigned int comp_i = fe->system_to_component_index(i).first;
1791  if (fe_mask[comp_i])
1792  p_real[fe_to_real[comp_i]] += data.local_dof_values[i] * data.shape(0,i);
1793  }
1794 
1795  return p_real;
1796 }
1797 
1798 
1799 
1800 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1801 Point<dim>
1804  const Point<spacedim> &p) const
1805 {
1806  // first a Newton iteration based on the real mapping. It uses the center
1807  // point of the cell as a starting point
1808  Point<dim> initial_p_unit;
1809  try
1810  {
1811  initial_p_unit
1812  = StaticMappingQ1<dim,spacedim>::mapping.transform_real_to_unit_cell(cell, p);
1813  }
1814  catch (const typename Mapping<dim,spacedim>::ExcTransformationFailed &)
1815  {
1816  // mirror the conditions of the code below to determine if we need to
1817  // use an arbitrary starting point or if we just need to rethrow the
1818  // exception
1819  for (unsigned int d=0; d<dim; ++d)
1820  initial_p_unit[d] = 0.5;
1821  }
1822 
1823  initial_p_unit = GeometryInfo<dim>::project_to_unit_cell(initial_p_unit);
1824 
1825  // for (unsigned int d=0; d<dim; ++d)
1826  // initial_p_unit[d] = 0.;
1827 
1828  const Quadrature<dim> point_quadrature(initial_p_unit);
1829 
1831  if (spacedim>dim)
1832  update_flags |= update_jacobian_grads;
1833  std_cxx11::unique_ptr<InternalData>
1834  mdata (get_data(update_flags,point_quadrature));
1835 
1836  update_internal_dofs(cell, *mdata);
1837 
1838  return do_transform_real_to_unit_cell(cell, p, initial_p_unit, *mdata);
1839 
1840 }
1841 
1842 
1843 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1844 Point<dim>
1848  const Point<spacedim> &p,
1849  const Point<dim> &initial_p_unit,
1850  InternalData &mdata) const
1851 {
1852  const unsigned int n_shapes=mdata.shape_values.size();
1853  (void)n_shapes;
1854  Assert(n_shapes!=0, ExcInternalError());
1855  AssertDimension (mdata.shape_derivatives.size(), n_shapes);
1856 
1857 
1858  // Newton iteration to solve
1859  // f(x)=p(x)-p=0
1860  // x_{n+1}=x_n-[f'(x)]^{-1}f(x)
1861  // The start value was set to be the
1862  // linear approximation to the cell
1863  // The shape values and derivatives
1864  // of the mapping at this point are
1865  // previously computed.
1866  // f(x)
1867  Point<dim> p_unit = initial_p_unit;
1868  Point<dim> f;
1869  compute_shapes_virtual(std::vector<Point<dim> > (1, p_unit), mdata);
1871  Tensor<1,spacedim> p_minus_F = p - p_real;
1872  const double eps = 1.e-12*cell->diameter();
1873  const unsigned int newton_iteration_limit = 20;
1874  unsigned int newton_iteration=0;
1875  while (p_minus_F.norm_square() > eps*eps)
1876  {
1877  // f'(x)
1878  Point<spacedim> DF[dim];
1879  Tensor<2,dim> df;
1880  for (unsigned int k=0; k<mdata.n_shape_functions; ++k)
1881  {
1882  const Tensor<1,dim> &grad_k = mdata.derivative(0,k);
1883  unsigned int comp_k = fe->system_to_component_index(k).first;
1884  if (fe_mask[comp_k])
1885  for (unsigned int j=0; j<dim; ++j)
1886  DF[j][fe_to_real[comp_k]] += mdata.local_dof_values[k] * grad_k[j];
1887  }
1888  for (unsigned int j=0; j<dim; ++j)
1889  {
1890  f[j] = DF[j] * p_minus_F;
1891  for (unsigned int l=0; l<dim; ++l)
1892  df[j][l] = -DF[j] * DF[l];
1893  }
1894  // Solve [f'(x)]d=f(x)
1895  const Tensor<1, dim> delta =
1896  invert(df) * static_cast<const Tensor<1, dim> &>(f);
1897  // do a line search
1898  double step_length = 1;
1899  do
1900  {
1901  // update of p_unit. The
1902  // spacedimth component of
1903  // transformed point is simply
1904  // ignored in codimension one
1905  // case. When this component is
1906  // not zero, then we are
1907  // projecting the point to the
1908  // surface or curve identified
1909  // by the cell.
1910  Point<dim> p_unit_trial = p_unit;
1911  for (unsigned int i=0; i<dim; ++i)
1912  p_unit_trial[i] -= step_length * delta[i];
1913  // shape values and derivatives
1914  // at new p_unit point
1915  compute_shapes_virtual(std::vector<Point<dim> > (1, p_unit_trial), mdata);
1916  // f(x)
1917  Point<spacedim> p_real_trial = do_transform_unit_to_real_cell(mdata);
1918  const Tensor<1,spacedim> f_trial = p - p_real_trial;
1919  // see if we are making progress with the current step length
1920  // and if not, reduce it by a factor of two and try again
1921  if (f_trial.norm() < p_minus_F.norm())
1922  {
1923  p_real = p_real_trial;
1924  p_unit = p_unit_trial;
1925  p_minus_F = f_trial;
1926  break;
1927  }
1928  else if (step_length > 0.05)
1929  step_length /= 2;
1930  else
1931  goto failure;
1932  }
1933  while (true);
1934  ++newton_iteration;
1935  if (newton_iteration > newton_iteration_limit)
1936  goto failure;
1937  }
1938  return p_unit;
1939  // if we get to the following label, then we have either run out
1940  // of Newton iterations, or the line search has not converged.
1941  // in either case, we need to give up, so throw an exception that
1942  // can then be caught
1943 failure:
1945  // ...the compiler wants us to return something, though we can
1946  // of course never get here...
1947  return Point<dim>();
1948 }
1949 
1950 
1951 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1952 unsigned int
1954 {
1955  return fe->degree;
1956 }
1957 
1958 
1959 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1962 {
1963  return this->fe_mask;
1964 }
1965 
1966 
1967 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1970 {
1972 }
1973 
1974 
1975 template<int dim, int spacedim, typename VectorType, typename DoFHandlerType>
1976 void
1980 {
1981  Assert(euler_dof_handler != 0, ExcMessage("euler_dof_handler is empty"));
1982 
1983  typename DoFHandlerType::cell_iterator dof_cell(*cell, euler_dof_handler);
1984  Assert (dof_cell->active() == true, ExcInactiveCell());
1985 
1986  dof_cell->get_dof_indices(data.local_dof_indices);
1987 
1988  for (unsigned int i=0; i<data.local_dof_values.size(); ++i)
1989  data.local_dof_values[i] = (*euler_vector)(data.local_dof_indices[i]);
1990 }
1991 
1992 // explicit instantiations
1993 #define SPLIT_INSTANTIATIONS_COUNT 2
1994 #ifndef SPLIT_INSTANTIATIONS_INDEX
1995 #define SPLIT_INSTANTIATIONS_INDEX 0
1996 #endif
1997 #include "mapping_fe_field.inst"
1998 
1999 
2000 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
static const unsigned int invalid_unsigned_int
Definition: types.h:170
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
virtual Mapping< dim, spacedim >::InternalDataBase * get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const
Contravariant transformation.
const std::vector< Point< dim > > & get_points() const
const std::vector< double > & get_weights() const
virtual std::size_t memory_consumption() const
MappingType
Definition: mapping.h:50
Volume element.
Outer normal vector, not normalized.
std::vector< unsigned int > fe_to_real
Determinant of the Jacobian.
Transformed quadrature points.
SmartPointer< const FiniteElement< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > fe
#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 ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static DataSetDescriptor cell()
Definition: qprojector.h:348
Point< spacedim > do_transform_unit_to_real_cell(const InternalData &mdata) const
std::vector< Tensor< 1, spacedim > > boundary_forms
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
unsigned int get_degree() const
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
virtual Mapping< dim, spacedim > * clone() const
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const
InternalData(const FiniteElement< dim, spacedim > &fe, const ComponentMask mask)
static ::ExceptionBase & ExcMessage(std::string arg1)
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
std::vector< Tensor< 3, dim > > shape_third_derivatives
#define Assert(cond, exc)
Definition: exceptions.h:313
UpdateFlags
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
Abstract base class for mapping classes.
Definition: dof_tools.h:46
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
std::vector< Tensor< 1, spacedim > > normal_vectors
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const
std::vector< double > volume_elements
numbers::NumberTraits< Number >::real_type norm_square() const
Definition: tensor.h:1000
unsigned int size() const
const double & shape(const unsigned int qpoint, const unsigned int shape_nr) const
virtual InternalData * get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const
DerivativeForm< 1, spacedim, dim, Number > transpose() const
virtual Mapping< dim, spacedim >::InternalDataBase * get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) 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)
SmartPointer< const DoFHandlerType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_dof_handler
virtual bool preserves_vertex_locations() const
Gradient of volume element.
std::vector< double > local_dof_values
const Tensor< 3, dim > & third_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
unsigned int size() const
const unsigned int dofs_per_cell
Definition: fe_base.h:295
std::vector< double > shape_values
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
friend class MappingFEField
std::vector< types::global_dof_index > local_dof_indices
virtual void compute_shapes_virtual(const std::vector< Point< dim > > &unit_points, typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
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)
const Tensor< 2, dim > & second_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
SmartPointer< const VectorType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > euler_vector
std::vector< Point< spacedim > > quadrature_points
Normal vectors.
const Tensor< 4, dim > & fourth_derivative(const unsigned int qpoint, const unsigned int shape_nr) const
static Point< dim > project_to_unit_cell(const Point< dim > &p)
static DataSetDescriptor subface(const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
const Tensor< 1, dim > & derivative(const unsigned int qpoint, const unsigned int shape_nr) const
static ::ExceptionBase & ExcNotImplemented()
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const
ComponentMask get_component_mask() const
std::vector< Tensor< 2, dim > > shape_second_derivatives
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 update_internal_dofs(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const
static ::ExceptionBase & ExcInactiveCell()
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
Point< 3 > point(const gp_Pnt &p)
Definition: utilities.cc:176
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
std::vector< Tensor< 1, dim > > shape_derivatives
Point< dim > do_transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit, InternalData &mdata) const
static double subface_ratio(const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
Covariant transformation.
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()
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Definition: quadrature.cc:1081