Reference documentation for deal.II version 9.1.0-pre
mapping_cartesian.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2018 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/array_view.h>
17 #include <deal.II/base/memory_consumption.h>
18 #include <deal.II/base/qprojector.h>
19 #include <deal.II/base/quadrature.h>
20 #include <deal.II/base/signaling_nan.h>
21 #include <deal.II/base/std_cxx14/memory.h>
22 #include <deal.II/base/tensor.h>
23 
24 #include <deal.II/dofs/dof_accessor.h>
25 
26 #include <deal.II/fe/fe_values.h>
27 #include <deal.II/fe/mapping_cartesian.h>
28 
29 #include <deal.II/grid/tria.h>
30 #include <deal.II/grid/tria_iterator.h>
31 
32 #include <deal.II/lac/full_matrix.h>
33 
34 #include <algorithm>
35 #include <cmath>
36 
37 
38 DEAL_II_NAMESPACE_OPEN
39 
40 
41 template <int dim, int spacedim>
43 
44 
45 
46 template <int dim, int spacedim>
48  const Quadrature<dim> &q)
49  : cell_extents(numbers::signaling_nan<Tensor<1, dim>>())
50  , volume_element(numbers::signaling_nan<double>())
51  , quadrature_points(q.get_points())
52 {}
53 
54 
55 
56 template <int dim, int spacedim>
57 std::size_t
59 {
64 }
65 
66 
67 
68 template <int dim, int spacedim>
69 bool
71 {
72  return true;
73 }
74 
75 
76 
77 template <int dim, int spacedim>
80  const UpdateFlags in) const
81 {
82  // this mapping is pretty simple in that it can basically compute
83  // every piece of information wanted by FEValues without requiring
84  // computing any other quantities. boundary forms are one exception
85  // since they can be computed from the normal vectors without much
86  // further ado
87  UpdateFlags out = in;
88  if (out & update_boundary_forms)
89  out |= update_normal_vectors;
90 
91  return out;
92 }
93 
94 
95 
96 template <int dim, int spacedim>
97 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
99  const Quadrature<dim> &q) const
100 {
101  auto data = std_cxx14::make_unique<InternalData>(q);
102 
103  // store the flags in the internal data object so we can access them
104  // in fill_fe_*_values(). use the transitive hull of the required
105  // flags
106  data->update_each = requires_update_flags(update_flags);
107 
108  return std::move(data);
109 }
110 
111 
112 
113 template <int dim, int spacedim>
114 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
116  const UpdateFlags update_flags,
117  const Quadrature<dim - 1> &quadrature) const
118 {
119  auto data = std_cxx14::make_unique<InternalData>(
121 
122  // verify that we have computed the transitive hull of the required
123  // flags and that FEValues has faithfully passed them on to us
124  Assert(update_flags == requires_update_flags(update_flags),
125  ExcInternalError());
126 
127  // store the flags in the internal data object so we can access them
128  // in fill_fe_*_values()
129  data->update_each = update_flags;
130 
131  return std::move(data);
132 }
133 
134 
135 
136 template <int dim, int spacedim>
137 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
139  const UpdateFlags update_flags,
140  const Quadrature<dim - 1> &quadrature) const
141 {
142  auto data = std_cxx14::make_unique<InternalData>(
144 
145  // verify that we have computed the transitive hull of the required
146  // flags and that FEValues has faithfully passed them on to us
147  Assert(update_flags == requires_update_flags(update_flags),
148  ExcInternalError());
149 
150  // store the flags in the internal data object so we can access them
151  // in fill_fe_*_values()
152  data->update_each = update_flags;
153 
154  return std::move(data);
155 }
156 
157 
158 
159 template <int dim, int spacedim>
160 void
162  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
163  const unsigned int face_no,
164  const unsigned int sub_no,
165  const CellSimilarity::Similarity cell_similarity,
166  const InternalData & data,
167  std::vector<Point<dim>> & quadrature_points,
168  std::vector<Tensor<1, dim>> &normal_vectors) const
169 {
170  const UpdateFlags update_flags = data.update_each;
171 
172  // some more sanity checks
173  if (face_no != invalid_face_number)
174  {
175  // Add 1 on both sides of
176  // assertion to avoid compiler
177  // warning about testing
178  // unsigned int < 0 in 1d.
179  Assert(face_no + 1 < GeometryInfo<dim>::faces_per_cell + 1,
181 
182  // We would like to check for
183  // sub_no < cell->face(face_no)->n_children(),
184  // but unfortunately the current
185  // function is also called for
186  // faces without children (see
187  // tests/fe/mapping.cc). Therefore,
188  // we must use following workaround
189  // of two separate assertions
190  Assert((sub_no == invalid_face_number) ||
191  cell->face(face_no)->has_children() ||
192  (sub_no + 1 < GeometryInfo<dim>::max_children_per_face + 1),
193  ExcIndexRange(sub_no,
194  0,
196  Assert((sub_no == invalid_face_number) ||
197  !cell->face(face_no)->has_children() ||
198  (sub_no < cell->face(face_no)->n_children()),
199  ExcIndexRange(sub_no, 0, cell->face(face_no)->n_children()));
200  }
201  else
202  // invalid face number, so
203  // subface should be invalid as
204  // well
206 
207  // let @p{start} be the origin of a
208  // local coordinate system. it is
209  // chosen as the (lower) left
210  // vertex
211  const Point<dim> start = cell->vertex(0);
212 
213  // Compute start point and sizes
214  // along axes. Strange vertex
215  // numbering makes this complicated
216  // again.
217  if (cell_similarity != CellSimilarity::translation)
218  {
219  switch (dim)
220  {
221  case 1:
222  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
223  break;
224  case 2:
225  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
226  data.cell_extents[1] = cell->vertex(2)(1) - start(1);
227  break;
228  case 3:
229  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
230  data.cell_extents[1] = cell->vertex(2)(1) - start(1);
231  data.cell_extents[2] = cell->vertex(4)(2) - start(2);
232  break;
233  default:
234  Assert(false, ExcNotImplemented());
235  }
236  }
237 
238 
239  // transform quadrature point. this
240  // is obtained simply by scaling
241  // unit coordinates with lengths in
242  // each direction
243  if (update_flags & update_quadrature_points)
244  {
245  const typename QProjector<dim>::DataSetDescriptor offset =
246  (face_no == invalid_face_number ?
248  (sub_no == invalid_face_number ?
249  // called from FEFaceValues
251  face_no,
252  cell->face_orientation(face_no),
253  cell->face_flip(face_no),
254  cell->face_rotation(face_no),
255  quadrature_points.size()) :
256  // called from FESubfaceValues
258  face_no,
259  sub_no,
260  cell->face_orientation(face_no),
261  cell->face_flip(face_no),
262  cell->face_rotation(face_no),
263  quadrature_points.size(),
264  cell->subface_case(face_no))));
265 
266  for (unsigned int i = 0; i < quadrature_points.size(); ++i)
267  {
268  quadrature_points[i] = start;
269  for (unsigned int d = 0; d < dim; ++d)
270  quadrature_points[i](d) +=
271  data.cell_extents[d] * data.quadrature_points[i + offset](d);
272  }
273  }
274 
275 
276  // compute normal vectors. since
277  // cells are aligned to coordinate
278  // axes, they are simply vectors
279  // with exactly one entry equal to
280  // 1 or -1. Furthermore, all
281  // normals on a face have the same
282  // value
283  if (update_flags & update_normal_vectors)
284  {
286 
287  switch (dim)
288  {
289  case 1:
290  {
291  static const Point<dim> normals[GeometryInfo<1>::faces_per_cell] =
292  {Point<dim>(-1.), Point<dim>(1.)};
293  std::fill(normal_vectors.begin(),
294  normal_vectors.end(),
295  normals[face_no]);
296  break;
297  }
298 
299  case 2:
300  {
301  static const Point<dim> normals[GeometryInfo<2>::faces_per_cell] =
302  {Point<dim>(-1, 0),
303  Point<dim>(1, 0),
304  Point<dim>(0, -1),
305  Point<dim>(0, 1)};
306  std::fill(normal_vectors.begin(),
307  normal_vectors.end(),
308  normals[face_no]);
309  break;
310  }
311 
312  case 3:
313  {
314  static const Point<dim> normals[GeometryInfo<3>::faces_per_cell] =
315  {Point<dim>(-1, 0, 0),
316  Point<dim>(1, 0, 0),
317  Point<dim>(0, -1, 0),
318  Point<dim>(0, 1, 0),
319  Point<dim>(0, 0, -1),
320  Point<dim>(0, 0, 1)};
321  std::fill(normal_vectors.begin(),
322  normal_vectors.end(),
323  normals[face_no]);
324  break;
325  }
326 
327  default:
328  Assert(false, ExcNotImplemented());
329  }
330  }
331 }
332 
333 
334 
335 template <int dim, int spacedim>
338  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
339  const CellSimilarity::Similarity cell_similarity,
340  const Quadrature<dim> & quadrature,
341  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
343  &output_data) const
344 {
345  // convert data object to internal data for this class. fails with
346  // an exception if that is not possible
347  Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
348  ExcInternalError());
349  const InternalData &data = static_cast<const InternalData &>(internal_data);
350 
351  std::vector<Tensor<1, dim>> dummy;
352 
353  compute_fill(cell,
356  cell_similarity,
357  data,
358  output_data.quadrature_points,
359  dummy);
360 
361  // compute Jacobian determinant. all values are equal and are the
362  // product of the local lengths in each coordinate direction
364  if (cell_similarity != CellSimilarity::translation)
365  {
366  double J = data.cell_extents[0];
367  for (unsigned int d = 1; d < dim; ++d)
368  J *= data.cell_extents[d];
369  data.volume_element = J;
370  if (data.update_each & update_JxW_values)
371  for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
372  output_data.JxW_values[i] = J * quadrature.weight(i);
373  }
374  // "compute" Jacobian at the quadrature points, which are all the
375  // same
376  if (data.update_each & update_jacobians)
377  if (cell_similarity != CellSimilarity::translation)
378  for (unsigned int i = 0; i < output_data.jacobians.size(); ++i)
379  {
380  output_data.jacobians[i] = DerivativeForm<1, dim, spacedim>();
381  for (unsigned int j = 0; j < dim; ++j)
382  output_data.jacobians[i][j][j] = data.cell_extents[j];
383  }
384  // "compute" the derivative of the Jacobian at the quadrature
385  // points, which are all zero of course
387  if (cell_similarity != CellSimilarity::translation)
388  for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
390 
392  if (cell_similarity != CellSimilarity::translation)
393  for (unsigned int i = 0;
394  i < output_data.jacobian_pushed_forward_grads.size();
395  ++i)
397 
398  // "compute" the hessian of the Jacobian at the quadrature points,
399  // which are all also zero of course
401  if (cell_similarity != CellSimilarity::translation)
402  for (unsigned int i = 0; i < output_data.jacobian_2nd_derivatives.size();
403  ++i)
404  output_data.jacobian_2nd_derivatives[i] =
406 
408  if (cell_similarity != CellSimilarity::translation)
409  for (unsigned int i = 0;
410  i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
411  ++i)
414 
416  if (cell_similarity != CellSimilarity::translation)
417  for (unsigned int i = 0; i < output_data.jacobian_3rd_derivatives.size();
418  ++i)
419  output_data.jacobian_3rd_derivatives[i] =
421 
423  if (cell_similarity != CellSimilarity::translation)
424  for (unsigned int i = 0;
425  i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
426  ++i)
429 
430  // "compute" inverse Jacobian at the quadrature points, which are
431  // all the same
433  if (cell_similarity != CellSimilarity::translation)
434  for (unsigned int i = 0; i < output_data.inverse_jacobians.size(); ++i)
435  {
436  output_data.inverse_jacobians[i] = Tensor<2, dim>();
437  for (unsigned int j = 0; j < dim; ++j)
438  output_data.inverse_jacobians[i][j][j] = 1. / data.cell_extents[j];
439  }
440 
441  return cell_similarity;
442 }
443 
444 
445 
446 template <int dim, int spacedim>
447 void
449  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
450  const unsigned int face_no,
451  const Quadrature<dim - 1> & quadrature,
452  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
454  &output_data) const
455 {
456  // convert data object to internal
457  // data for this class. fails with
458  // an exception if that is not
459  // possible
460  Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
461  ExcInternalError());
462  const InternalData &data = static_cast<const InternalData &>(internal_data);
463 
464  compute_fill(cell,
465  face_no,
468  data,
469  output_data.quadrature_points,
470  output_data.normal_vectors);
471 
472  // first compute Jacobian determinant, which is simply the product
473  // of the local lengths since the jacobian is diagonal
474  double J = 1.;
475  for (unsigned int d = 0; d < dim; ++d)
477  J *= data.cell_extents[d];
478 
479  if (data.update_each & update_JxW_values)
480  for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
481  output_data.JxW_values[i] = J * quadrature.weight(i);
482 
484  for (unsigned int i = 0; i < output_data.boundary_forms.size(); ++i)
485  output_data.boundary_forms[i] = J * output_data.normal_vectors[i];
486 
488  {
489  J = data.cell_extents[0];
490  for (unsigned int d = 1; d < dim; ++d)
491  J *= data.cell_extents[d];
492  data.volume_element = J;
493  }
494 
495  if (data.update_each & update_jacobians)
496  for (unsigned int i = 0; i < output_data.jacobians.size(); ++i)
497  {
498  output_data.jacobians[i] = DerivativeForm<1, dim, spacedim>();
499  for (unsigned int d = 0; d < dim; ++d)
500  output_data.jacobians[i][d][d] = data.cell_extents[d];
501  }
502 
504  for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
506 
508  for (unsigned int i = 0;
509  i < output_data.jacobian_pushed_forward_grads.size();
510  ++i)
512 
514  for (unsigned int i = 0; i < output_data.jacobian_2nd_derivatives.size();
515  ++i)
516  output_data.jacobian_2nd_derivatives[i] =
518 
520  for (unsigned int i = 0;
521  i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
522  ++i)
525 
527  for (unsigned int i = 0; i < output_data.jacobian_3rd_derivatives.size();
528  ++i)
529  output_data.jacobian_3rd_derivatives[i] =
531 
533  for (unsigned int i = 0;
534  i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
535  ++i)
538 
540  for (unsigned int i = 0; i < output_data.inverse_jacobians.size(); ++i)
541  {
543  for (unsigned int d = 0; d < dim; ++d)
544  output_data.inverse_jacobians[i][d][d] = 1. / data.cell_extents[d];
545  }
546 }
547 
548 
549 
550 template <int dim, int spacedim>
551 void
553  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
554  const unsigned int face_no,
555  const unsigned int subface_no,
556  const Quadrature<dim - 1> & quadrature,
557  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
559  &output_data) const
560 {
561  // convert data object to internal data for this class. fails with
562  // an exception if that is not possible
563  Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
564  ExcInternalError());
565  const InternalData &data = static_cast<const InternalData &>(internal_data);
566 
567  compute_fill(cell,
568  face_no,
569  subface_no,
571  data,
572  output_data.quadrature_points,
573  output_data.normal_vectors);
574 
575  // first compute Jacobian determinant, which is simply the product
576  // of the local lengths since the jacobian is diagonal
577  double J = 1.;
578  for (unsigned int d = 0; d < dim; ++d)
580  J *= data.cell_extents[d];
581 
582  if (data.update_each & update_JxW_values)
583  {
584  // Here, cell->face(face_no)->n_children() would be the right
585  // choice, but unfortunately the current function is also called
586  // for faces without children (see tests/fe/mapping.cc). Add
587  // following switch to avoid diffs in tests/fe/mapping.OK
588  const unsigned int n_subfaces =
589  cell->face(face_no)->has_children() ?
590  cell->face(face_no)->n_children() :
592  for (unsigned int i = 0; i < output_data.JxW_values.size(); ++i)
593  output_data.JxW_values[i] = J * quadrature.weight(i) / n_subfaces;
594  }
595 
597  for (unsigned int i = 0; i < output_data.boundary_forms.size(); ++i)
598  output_data.boundary_forms[i] = J * output_data.normal_vectors[i];
599 
601  {
602  J = data.cell_extents[0];
603  for (unsigned int d = 1; d < dim; ++d)
604  J *= data.cell_extents[d];
605  data.volume_element = J;
606  }
607 
608  if (data.update_each & update_jacobians)
609  for (unsigned int i = 0; i < output_data.jacobians.size(); ++i)
610  {
611  output_data.jacobians[i] = DerivativeForm<1, dim, spacedim>();
612  for (unsigned int d = 0; d < dim; ++d)
613  output_data.jacobians[i][d][d] = data.cell_extents[d];
614  }
615 
617  for (unsigned int i = 0; i < output_data.jacobian_grads.size(); ++i)
619 
621  for (unsigned int i = 0;
622  i < output_data.jacobian_pushed_forward_grads.size();
623  ++i)
625 
627  for (unsigned int i = 0; i < output_data.jacobian_2nd_derivatives.size();
628  ++i)
629  output_data.jacobian_2nd_derivatives[i] =
631 
633  for (unsigned int i = 0;
634  i < output_data.jacobian_pushed_forward_2nd_derivatives.size();
635  ++i)
638 
640  for (unsigned int i = 0; i < output_data.jacobian_3rd_derivatives.size();
641  ++i)
642  output_data.jacobian_3rd_derivatives[i] =
644 
646  for (unsigned int i = 0;
647  i < output_data.jacobian_pushed_forward_3rd_derivatives.size();
648  ++i)
651 
653  for (unsigned int i = 0; i < output_data.inverse_jacobians.size(); ++i)
654  {
656  for (unsigned int d = 0; d < dim; ++d)
657  output_data.inverse_jacobians[i][d][d] = 1. / data.cell_extents[d];
658  }
659 }
660 
661 
662 
663 template <int dim, int spacedim>
664 void
666  const ArrayView<const Tensor<1, dim>> & input,
667  const MappingType mapping_type,
668  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
669  const ArrayView<Tensor<1, spacedim>> & output) const
670 {
671  AssertDimension(input.size(), output.size());
672  Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
673  ExcInternalError());
674  const InternalData &data = static_cast<const InternalData &>(mapping_data);
675 
676  switch (mapping_type)
677  {
678  case mapping_covariant:
679  {
682  "update_covariant_transformation"));
683 
684  for (unsigned int i = 0; i < output.size(); ++i)
685  for (unsigned int d = 0; d < dim; ++d)
686  output[i][d] = input[i][d] / data.cell_extents[d];
687  return;
688  }
689 
691  {
694  "update_contravariant_transformation"));
695 
696  for (unsigned int i = 0; i < output.size(); ++i)
697  for (unsigned int d = 0; d < dim; ++d)
698  output[i][d] = input[i][d] * data.cell_extents[d];
699  return;
700  }
701  case mapping_piola:
702  {
705  "update_contravariant_transformation"));
708  "update_volume_elements"));
709 
710  for (unsigned int i = 0; i < output.size(); ++i)
711  for (unsigned int d = 0; d < dim; ++d)
712  output[i][d] =
713  input[i][d] * data.cell_extents[d] / data.volume_element;
714  return;
715  }
716  default:
717  Assert(false, ExcNotImplemented());
718  }
719 }
720 
721 
722 
723 template <int dim, int spacedim>
724 void
726  const ArrayView<const DerivativeForm<1, dim, spacedim>> &input,
727  const MappingType mapping_type,
728  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
729  const ArrayView<Tensor<2, spacedim>> & output) const
730 {
731  AssertDimension(input.size(), output.size());
732  Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
733  ExcInternalError());
734  const InternalData &data = static_cast<const InternalData &>(mapping_data);
735 
736  switch (mapping_type)
737  {
738  case mapping_covariant:
739  {
742  "update_covariant_transformation"));
743 
744  for (unsigned int i = 0; i < output.size(); ++i)
745  for (unsigned int d1 = 0; d1 < dim; ++d1)
746  for (unsigned int d2 = 0; d2 < dim; ++d2)
747  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2];
748  return;
749  }
750 
752  {
755  "update_contravariant_transformation"));
756 
757  for (unsigned int i = 0; i < output.size(); ++i)
758  for (unsigned int d1 = 0; d1 < dim; ++d1)
759  for (unsigned int d2 = 0; d2 < dim; ++d2)
760  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
761  return;
762  }
763 
765  {
768  "update_covariant_transformation"));
769 
770  for (unsigned int i = 0; i < output.size(); ++i)
771  for (unsigned int d1 = 0; d1 < dim; ++d1)
772  for (unsigned int d2 = 0; d2 < dim; ++d2)
773  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2] /
774  data.cell_extents[d1];
775  return;
776  }
777 
779  {
782  "update_contravariant_transformation"));
783 
784  for (unsigned int i = 0; i < output.size(); ++i)
785  for (unsigned int d1 = 0; d1 < dim; ++d1)
786  for (unsigned int d2 = 0; d2 < dim; ++d2)
787  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
788  data.cell_extents[d1];
789  return;
790  }
791 
792  case mapping_piola:
793  {
796  "update_contravariant_transformation"));
799  "update_volume_elements"));
800 
801  for (unsigned int i = 0; i < output.size(); ++i)
802  for (unsigned int d1 = 0; d1 < dim; ++d1)
803  for (unsigned int d2 = 0; d2 < dim; ++d2)
804  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
805  data.volume_element;
806  return;
807  }
808 
810  {
813  "update_contravariant_transformation"));
816  "update_volume_elements"));
817 
818  for (unsigned int i = 0; i < output.size(); ++i)
819  for (unsigned int d1 = 0; d1 < dim; ++d1)
820  for (unsigned int d2 = 0; d2 < dim; ++d2)
821  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
822  data.cell_extents[d1] / data.volume_element;
823  return;
824  }
825 
826  default:
827  Assert(false, ExcNotImplemented());
828  }
829 }
830 
831 
832 
833 template <int dim, int spacedim>
834 void
836  const ArrayView<const Tensor<2, dim>> & input,
837  const MappingType mapping_type,
838  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
839  const ArrayView<Tensor<2, spacedim>> & output) const
840 {
841  AssertDimension(input.size(), output.size());
842  Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
843  ExcInternalError());
844  const InternalData &data = static_cast<const InternalData &>(mapping_data);
845 
846  switch (mapping_type)
847  {
848  case mapping_covariant:
849  {
852  "update_covariant_transformation"));
853 
854  for (unsigned int i = 0; i < output.size(); ++i)
855  for (unsigned int d1 = 0; d1 < dim; ++d1)
856  for (unsigned int d2 = 0; d2 < dim; ++d2)
857  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2];
858  return;
859  }
860 
862  {
865  "update_contravariant_transformation"));
866 
867  for (unsigned int i = 0; i < output.size(); ++i)
868  for (unsigned int d1 = 0; d1 < dim; ++d1)
869  for (unsigned int d2 = 0; d2 < dim; ++d2)
870  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
871  return;
872  }
873 
875  {
878  "update_covariant_transformation"));
879 
880  for (unsigned int i = 0; i < output.size(); ++i)
881  for (unsigned int d1 = 0; d1 < dim; ++d1)
882  for (unsigned int d2 = 0; d2 < dim; ++d2)
883  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2] /
884  data.cell_extents[d1];
885  return;
886  }
887 
889  {
892  "update_contravariant_transformation"));
893 
894  for (unsigned int i = 0; i < output.size(); ++i)
895  for (unsigned int d1 = 0; d1 < dim; ++d1)
896  for (unsigned int d2 = 0; d2 < dim; ++d2)
897  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
898  data.cell_extents[d1];
899  return;
900  }
901 
902  case mapping_piola:
903  {
906  "update_contravariant_transformation"));
909  "update_volume_elements"));
910 
911  for (unsigned int i = 0; i < output.size(); ++i)
912  for (unsigned int d1 = 0; d1 < dim; ++d1)
913  for (unsigned int d2 = 0; d2 < dim; ++d2)
914  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
915  data.volume_element;
916  return;
917  }
918 
920  {
923  "update_contravariant_transformation"));
926  "update_volume_elements"));
927 
928  for (unsigned int i = 0; i < output.size(); ++i)
929  for (unsigned int d1 = 0; d1 < dim; ++d1)
930  for (unsigned int d2 = 0; d2 < dim; ++d2)
931  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] /
932  data.cell_extents[d1] / data.volume_element;
933  return;
934  }
935 
936  default:
937  Assert(false, ExcNotImplemented());
938  }
939 }
940 
941 
942 template <int dim, int spacedim>
943 void
945  const ArrayView<const DerivativeForm<2, dim, spacedim>> &input,
946  const MappingType mapping_type,
947  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
948  const ArrayView<Tensor<3, spacedim>> & output) const
949 {
950  AssertDimension(input.size(), output.size());
951  Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
952  ExcInternalError());
953  const InternalData &data = static_cast<const InternalData &>(mapping_data);
954 
955  switch (mapping_type)
956  {
958  {
961  "update_covariant_transformation"));
962 
963  for (unsigned int q = 0; q < output.size(); ++q)
964  for (unsigned int i = 0; i < spacedim; ++i)
965  for (unsigned int j = 0; j < spacedim; ++j)
966  for (unsigned int k = 0; k < spacedim; ++k)
967  {
968  output[q][i][j][k] = input[q][i][j][k] /
969  data.cell_extents[j] /
970  data.cell_extents[k];
971  }
972  return;
973  }
974  default:
975  Assert(false, ExcNotImplemented());
976  }
977 }
978 
979 template <int dim, int spacedim>
980 void
982  const ArrayView<const Tensor<3, dim>> & input,
983  const MappingType mapping_type,
984  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
985  const ArrayView<Tensor<3, spacedim>> & output) const
986 {
987  AssertDimension(input.size(), output.size());
988  Assert(dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
989  ExcInternalError());
990  const InternalData &data = static_cast<const InternalData &>(mapping_data);
991 
992  switch (mapping_type)
993  {
995  {
998  "update_covariant_transformation"));
1001  "update_contravariant_transformation"));
1002 
1003  for (unsigned int q = 0; q < output.size(); ++q)
1004  for (unsigned int i = 0; i < spacedim; ++i)
1005  for (unsigned int j = 0; j < spacedim; ++j)
1006  for (unsigned int k = 0; k < spacedim; ++k)
1007  {
1008  output[q][i][j][k] =
1009  input[q][i][j][k] * data.cell_extents[i] /
1010  data.cell_extents[j] / data.cell_extents[k];
1011  }
1012  return;
1013  }
1014 
1016  {
1019  "update_covariant_transformation"));
1020 
1021  for (unsigned int q = 0; q < output.size(); ++q)
1022  for (unsigned int i = 0; i < spacedim; ++i)
1023  for (unsigned int j = 0; j < spacedim; ++j)
1024  for (unsigned int k = 0; k < spacedim; ++k)
1025  {
1026  output[q][i][j][k] =
1027  input[q][i][j][k] / data.cell_extents[i] /
1028  data.cell_extents[j] / data.cell_extents[k];
1029  }
1030 
1031  return;
1032  }
1033 
1034  case mapping_piola_hessian:
1035  {
1038  "update_covariant_transformation"));
1041  "update_contravariant_transformation"));
1044  "update_volume_elements"));
1045 
1046  for (unsigned int q = 0; q < output.size(); ++q)
1047  for (unsigned int i = 0; i < spacedim; ++i)
1048  for (unsigned int j = 0; j < spacedim; ++j)
1049  for (unsigned int k = 0; k < spacedim; ++k)
1050  {
1051  output[q][i][j][k] =
1052  input[q][i][j][k] * data.cell_extents[i] /
1053  data.volume_element / data.cell_extents[j] /
1054  data.cell_extents[k];
1055  }
1056 
1057  return;
1058  }
1059 
1060  default:
1061  Assert(false, ExcNotImplemented());
1062  }
1063 }
1064 
1065 
1066 template <int dim, int spacedim>
1069  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
1070  const Point<dim> & p) const
1071 {
1072  Tensor<1, dim> length;
1073  const Point<dim> start = cell->vertex(0);
1074  switch (dim)
1075  {
1076  case 1:
1077  length[0] = cell->vertex(1)(0) - start(0);
1078  break;
1079  case 2:
1080  length[0] = cell->vertex(1)(0) - start(0);
1081  length[1] = cell->vertex(2)(1) - start(1);
1082  break;
1083  case 3:
1084  length[0] = cell->vertex(1)(0) - start(0);
1085  length[1] = cell->vertex(2)(1) - start(1);
1086  length[2] = cell->vertex(4)(2) - start(2);
1087  break;
1088  default:
1089  Assert(false, ExcNotImplemented());
1090  }
1091 
1092  Point<dim> p_real = cell->vertex(0);
1093  for (unsigned int d = 0; d < dim; ++d)
1094  p_real(d) += length[d] * p(d);
1095 
1096  return p_real;
1097 }
1098 
1099 
1100 
1101 template <int dim, int spacedim>
1102 Point<dim>
1104  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
1105  const Point<spacedim> & p) const
1106 {
1107  if (dim != spacedim)
1108  Assert(false, ExcNotImplemented());
1109  const Point<dim> &start = cell->vertex(0);
1110  Point<dim> real = p;
1111  real -= start;
1112 
1113  switch (dim)
1114  {
1115  case 1:
1116  real(0) /= cell->vertex(1)(0) - start(0);
1117  break;
1118  case 2:
1119  real(0) /= cell->vertex(1)(0) - start(0);
1120  real(1) /= cell->vertex(2)(1) - start(1);
1121  break;
1122  case 3:
1123  real(0) /= cell->vertex(1)(0) - start(0);
1124  real(1) /= cell->vertex(2)(1) - start(1);
1125  real(2) /= cell->vertex(4)(2) - start(2);
1126  break;
1127  default:
1128  Assert(false, ExcNotImplemented());
1129  }
1130  return real;
1131 }
1132 
1133 
1134 template <int dim, int spacedim>
1135 std::unique_ptr<Mapping<dim, spacedim>>
1137 {
1138  return std_cxx14::make_unique<MappingCartesian<dim, spacedim>>(*this);
1139 }
1140 
1141 
1142 //---------------------------------------------------------------------------
1143 // explicit instantiations
1144 #include "mapping_cartesian.inst"
1145 
1146 
1147 DEAL_II_NAMESPACE_CLOSE
Transformed quadrature weights.
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
Contravariant transformation.
MappingType
Definition: mapping.h:61
std::vector< Tensor< 1, spacedim > > boundary_forms
Volume element.
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
Outer normal vector, not normalized.
std::vector< Point< dim > > quadrature_points
Determinant of the Jacobian.
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
Transformed quadrature points.
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
static::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static DataSetDescriptor cell()
Definition: qprojector.h:344
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
static Quadrature< dim > project_to_all_subfaces(const SubQuadrature &quadrature)
virtual std::size_t memory_consumption() const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
#define Assert(cond, exc)
Definition: exceptions.h:1227
UpdateFlags
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 override
double weight(const unsigned int i) const
Abstract base class for mapping classes.
Definition: dof_tools.h:57
Gradient of volume element.
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
static const unsigned int invalid_face_number
std::vector< Point< spacedim > > quadrature_points
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
InternalData(const Quadrature< dim > &quadrature)
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Definition: mpi.h:55
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual bool preserves_vertex_locations() const override
Normal vectors.
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
static::ExceptionBase & ExcNotImplemented()
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
std::vector< DerivativeForm< 2, dim, spacedim > > jacobian_grads
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
void compute_fill(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const CellSimilarity::Similarity cell_similarity, const InternalData &data, std::vector< Point< dim >> &quadrature_points, std::vector< Tensor< 1, dim >> &normal_vectors) const
Covariant transformation.
std::vector< Tensor< 1, spacedim > > normal_vectors
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
UpdateFlags update_each
Definition: mapping.h:576
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:1139