Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/std_cxx14/memory.h>
17 #include <deal.II/base/array_view.h>
18 #include <deal.II/base/tensor.h>
19 #include <deal.II/base/quadrature.h>
20 #include <deal.II/base/qprojector.h>
21 #include <deal.II/base/signaling_nan.h>
22 #include <deal.II/base/memory_consumption.h>
23 #include <deal.II/lac/full_matrix.h>
24 #include <deal.II/grid/tria.h>
25 #include <deal.II/grid/tria_iterator.h>
26 #include <deal.II/dofs/dof_accessor.h>
27 #include <deal.II/fe/mapping_cartesian.h>
28 #include <deal.II/fe/fe_values.h>
29 
30 #include <cmath>
31 #include <algorithm>
32 
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 
37 template <int dim, int spacedim>
39 
40 
41 
42 template <int dim, int spacedim>
44  :
45  cell_extents (numbers::signaling_nan<Tensor<1,dim> >()),
46  volume_element (numbers::signaling_nan<double>()),
47  quadrature_points (q.get_points ())
48 {}
49 
50 
51 
52 template <int dim, int spacedim>
53 std::size_t
55 {
58  MemoryConsumption::memory_consumption (volume_element) +
59  MemoryConsumption::memory_consumption (quadrature_points));
60 }
61 
62 
63 
64 template <int dim, int spacedim>
65 bool
67 {
68  return true;
69 }
70 
71 
72 
73 template <int dim, int spacedim>
76 {
77  // this mapping is pretty simple in that it can basically compute
78  // every piece of information wanted by FEValues without requiring
79  // computing any other quantities. boundary forms are one exception
80  // since they can be computed from the normal vectors without much
81  // further ado
82  UpdateFlags out = in;
83  if (out & update_boundary_forms)
84  out |= update_normal_vectors;
85 
86  return out;
87 }
88 
89 
90 
91 template <int dim, int spacedim>
92 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
94  const Quadrature<dim> &q) const
95 {
96  auto data = std_cxx14::make_unique<InternalData> (q);
97 
98  // store the flags in the internal data object so we can access them
99  // in fill_fe_*_values(). use the transitive hull of the required
100  // flags
101  data->update_each = requires_update_flags(update_flags);
102 
103  return std::move(data);
104 }
105 
106 
107 
108 template <int dim, int spacedim>
109 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
111  const Quadrature<dim-1>& quadrature) const
112 {
113  auto data = std_cxx14::make_unique<InternalData> (QProjector<dim>::project_to_all_faces(quadrature));
114 
115  // verify that we have computed the transitive hull of the required
116  // flags and that FEValues has faithfully passed them on to us
117  Assert (update_flags == requires_update_flags (update_flags),
118  ExcInternalError());
119 
120  // store the flags in the internal data object so we can access them
121  // in fill_fe_*_values()
122  data->update_each = update_flags;
123 
124  return std::move(data);
125 }
126 
127 
128 
129 template <int dim, int spacedim>
130 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
132  const Quadrature<dim-1> &quadrature) const
133 {
134  auto data = std_cxx14::make_unique<InternalData> (QProjector<dim>::project_to_all_subfaces(quadrature));
135 
136  // verify that we have computed the transitive hull of the required
137  // flags and that FEValues has faithfully passed them on to us
138  Assert (update_flags == requires_update_flags (update_flags),
139  ExcInternalError());
140 
141  // store the flags in the internal data object so we can access them
142  // in fill_fe_*_values()
143  data->update_each = update_flags;
144 
145  return std::move(data);
146 }
147 
148 
149 
150 
151 template <int dim, int spacedim>
152 void
154  const unsigned int face_no,
155  const unsigned int sub_no,
156  const CellSimilarity::Similarity cell_similarity,
157  const InternalData &data,
158  std::vector<Point<dim> > &quadrature_points,
159  std::vector<Tensor<1,dim> > &normal_vectors) const
160 {
161  const UpdateFlags update_flags = data.update_each;
162 
163  // some more sanity checks
164  if (face_no != invalid_face_number)
165  {
166  // Add 1 on both sides of
167  // assertion to avoid compiler
168  // warning about testing
169  // unsigned int < 0 in 1d.
172 
173  // We would like to check for
174  // sub_no < cell->face(face_no)->n_children(),
175  // but unfortunately the current
176  // function is also called for
177  // faces without children (see
178  // tests/fe/mapping.cc). Therefore,
179  // we must use following workaround
180  // of two separate assertions
181  Assert ((sub_no == invalid_face_number) ||
182  cell->face(face_no)->has_children() ||
184  ExcIndexRange (sub_no, 0,
186  Assert ((sub_no == invalid_face_number) ||
187  !cell->face(face_no)->has_children() ||
188  (sub_no < cell->face(face_no)->n_children()),
189  ExcIndexRange (sub_no, 0, cell->face(face_no)->n_children()));
190  }
191  else
192  // invalid face number, so
193  // subface should be invalid as
194  // well
196 
197  // let @p{start} be the origin of a
198  // local coordinate system. it is
199  // chosen as the (lower) left
200  // vertex
201  const Point<dim> start = cell->vertex(0);
202 
203  // Compute start point and sizes
204  // along axes. Strange vertex
205  // numbering makes this complicated
206  // again.
207  if (cell_similarity != CellSimilarity::translation)
208  {
209  switch (dim)
210  {
211  case 1:
212  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
213  break;
214  case 2:
215  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
216  data.cell_extents[1] = cell->vertex(2)(1) - start(1);
217  break;
218  case 3:
219  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
220  data.cell_extents[1] = cell->vertex(2)(1) - start(1);
221  data.cell_extents[2] = cell->vertex(4)(2) - start(2);
222  break;
223  default:
224  Assert(false, ExcNotImplemented());
225  }
226  }
227 
228 
229  // transform quadrature point. this
230  // is obtained simply by scaling
231  // unit coordinates with lengths in
232  // each direction
233  if (update_flags & update_quadrature_points)
234  {
235  const typename QProjector<dim>::DataSetDescriptor offset
236  = (face_no == invalid_face_number
237  ?
239  :
240  (sub_no == invalid_face_number
241  ?
242  // called from FEFaceValues
244  cell->face_orientation(face_no),
245  cell->face_flip(face_no),
246  cell->face_rotation(face_no),
247  quadrature_points.size())
248  :
249  // called from FESubfaceValues
251  cell->face_orientation(face_no),
252  cell->face_flip(face_no),
253  cell->face_rotation(face_no),
254  quadrature_points.size(),
255  cell->subface_case(face_no))
256  ));
257 
258  for (unsigned int i=0; i<quadrature_points.size(); ++i)
259  {
260  quadrature_points[i] = start;
261  for (unsigned int d=0; d<dim; ++d)
262  quadrature_points[i](d) += data.cell_extents[d] *
263  data.quadrature_points[i+offset](d);
264  }
265  }
266 
267 
268  // compute normal vectors. since
269  // cells are aligned to coordinate
270  // axes, they are simply vectors
271  // with exactly one entry equal to
272  // 1 or -1. Furthermore, all
273  // normals on a face have the same
274  // value
275  if (update_flags & update_normal_vectors)
276  {
278  ExcInternalError());
279 
280  switch (dim)
281  {
282  case 1:
283  {
284  static const Point<dim>
286  = { Point<dim>(-1.),
287  Point<dim>( 1.)
288  };
289  std::fill (normal_vectors.begin(),
290  normal_vectors.end(),
291  normals[face_no]);
292  break;
293  }
294 
295  case 2:
296  {
297  static const Point<dim>
299  = { Point<dim>(-1, 0),
300  Point<dim>( 1, 0),
301  Point<dim>( 0,-1),
302  Point<dim>( 0, 1)
303  };
304  std::fill (normal_vectors.begin(),
305  normal_vectors.end(),
306  normals[face_no]);
307  break;
308  }
309 
310  case 3:
311  {
312  static const Point<dim>
314  = { Point<dim>(-1, 0, 0),
315  Point<dim>( 1, 0, 0),
316  Point<dim>( 0,-1, 0),
317  Point<dim>( 0, 1, 0),
318  Point<dim>( 0, 0,-1),
319  Point<dim>( 0, 0, 1)
320  };
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>
339  const CellSimilarity::Similarity cell_similarity,
340  const Quadrature<dim> &quadrature,
341  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
343 {
344  // convert data object to internal data for this class. fails with
345  // an exception if that is not possible
346  Assert (dynamic_cast<const InternalData *> (&internal_data) != nullptr, ExcInternalError());
347  const InternalData &data = static_cast<const InternalData &> (internal_data);
348 
349  std::vector<Tensor<1,dim> > dummy;
350 
351  compute_fill (cell, invalid_face_number, invalid_face_number, cell_similarity,
352  data,
353  output_data.quadrature_points,
354  dummy);
355 
356  // compute Jacobian determinant. all values are equal and are the
357  // product of the local lengths in each coordinate direction
358  if (data.update_each & (update_JxW_values | update_volume_elements))
359  if (cell_similarity != CellSimilarity::translation)
360  {
361  double J = data.cell_extents[0];
362  for (unsigned int d=1; d<dim; ++d)
363  J *= data.cell_extents[d];
364  data.volume_element = J;
365  if (data.update_each & update_JxW_values)
366  for (unsigned int i=0; i<output_data.JxW_values.size(); ++i)
367  output_data.JxW_values[i] = J * quadrature.weight(i);
368  }
369  // "compute" Jacobian at the quadrature points, which are all the
370  // same
371  if (data.update_each & update_jacobians)
372  if (cell_similarity != CellSimilarity::translation)
373  for (unsigned int i=0; i<output_data.jacobians.size(); ++i)
374  {
375  output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>();
376  for (unsigned int j=0; j<dim; ++j)
377  output_data.jacobians[i][j][j] = data.cell_extents[j];
378  }
379  // "compute" the derivative of the Jacobian at the quadrature
380  // points, which are all zero of course
381  if (data.update_each & update_jacobian_grads)
382  if (cell_similarity != CellSimilarity::translation)
383  for (unsigned int i=0; i<output_data.jacobian_grads.size(); ++i)
385 
386  if (data.update_each & update_jacobian_pushed_forward_grads)
387  if (cell_similarity != CellSimilarity::translation)
388  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_grads.size(); ++i)
390 
391  // "compute" the hessian of the Jacobian at the quadrature points,
392  // which are all also zero of course
393  if (data.update_each & update_jacobian_2nd_derivatives)
394  if (cell_similarity != CellSimilarity::translation)
395  for (unsigned int i=0; i<output_data.jacobian_2nd_derivatives.size(); ++i)
397 
398  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
399  if (cell_similarity != CellSimilarity::translation)
400  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_2nd_derivatives.size(); ++i)
402 
403  if (data.update_each & update_jacobian_3rd_derivatives)
404  if (cell_similarity != CellSimilarity::translation)
405  for (unsigned int i=0; i<output_data.jacobian_3rd_derivatives.size(); ++i)
407 
408  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
409  if (cell_similarity != CellSimilarity::translation)
410  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_3rd_derivatives.size(); ++i)
412 
413  // "compute" inverse Jacobian at the quadrature points, which are
414  // all the same
415  if (data.update_each & update_inverse_jacobians)
416  if (cell_similarity != CellSimilarity::translation)
417  for (unsigned int i=0; i<output_data.inverse_jacobians.size(); ++i)
418  {
419  output_data.inverse_jacobians[i] = Tensor<2,dim>();
420  for (unsigned int j=0; j<dim; ++j)
421  output_data.inverse_jacobians[i][j][j]=1./data.cell_extents[j];
422  }
423 
424  return cell_similarity;
425 }
426 
427 
428 
429 template <int dim, int spacedim>
430 void
433  const unsigned int face_no,
434  const Quadrature<dim-1> &quadrature,
435  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
437 {
438  // convert data object to internal
439  // data for this class. fails with
440  // an exception if that is not
441  // possible
442  Assert (dynamic_cast<const InternalData *> (&internal_data) != nullptr,
443  ExcInternalError());
444  const InternalData &data = static_cast<const InternalData &> (internal_data);
445 
446  compute_fill (cell, face_no, invalid_face_number,
448  data,
449  output_data.quadrature_points,
450  output_data.normal_vectors);
451 
452  // first compute Jacobian determinant, which is simply the product
453  // of the local lengths since the jacobian is diagonal
454  double J = 1.;
455  for (unsigned int d=0; d<dim; ++d)
457  J *= data.cell_extents[d];
458 
459  if (data.update_each & update_JxW_values)
460  for (unsigned int i=0; i<output_data.JxW_values.size(); ++i)
461  output_data.JxW_values[i] = J * quadrature.weight(i);
462 
463  if (data.update_each & update_boundary_forms)
464  for (unsigned int i=0; i<output_data.boundary_forms.size(); ++i)
465  output_data.boundary_forms[i] = J * output_data.normal_vectors[i];
466 
467  if (data.update_each & update_volume_elements)
468  {
469  J = data.cell_extents[0];
470  for (unsigned int d=1; d<dim; ++d)
471  J *= data.cell_extents[d];
472  data.volume_element = J;
473  }
474 
475  if (data.update_each & update_jacobians)
476  for (unsigned int i=0; i<output_data.jacobians.size(); ++i)
477  {
478  output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>();
479  for (unsigned int d=0; d<dim; ++d)
480  output_data.jacobians[i][d][d] = data.cell_extents[d];
481  }
482 
483  if (data.update_each & update_jacobian_grads)
484  for (unsigned int i=0; i<output_data.jacobian_grads.size(); ++i)
486 
487  if (data.update_each & update_jacobian_pushed_forward_grads)
488  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_grads.size(); ++i)
490 
491  if (data.update_each & update_jacobian_2nd_derivatives)
492  for (unsigned int i=0; i<output_data.jacobian_2nd_derivatives.size(); ++i)
494 
495  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
496  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_2nd_derivatives.size(); ++i)
498 
499  if (data.update_each & update_jacobian_3rd_derivatives)
500  for (unsigned int i=0; i<output_data.jacobian_3rd_derivatives.size(); ++i)
502 
503  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
504  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_3rd_derivatives.size(); ++i)
506 
507  if (data.update_each & update_inverse_jacobians)
508  for (unsigned int i=0; i<output_data.inverse_jacobians.size(); ++i)
509  {
511  for (unsigned int d=0; d<dim; ++d)
512  output_data.inverse_jacobians[i][d][d] = 1./data.cell_extents[d];
513  }
514 }
515 
516 
517 
518 template <int dim, int spacedim>
519 void
522  const unsigned int face_no,
523  const unsigned int subface_no,
524  const Quadrature<dim-1> &quadrature,
525  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
527 {
528  // convert data object to internal data for this class. fails with
529  // an exception if that is not possible
530  Assert (dynamic_cast<const InternalData *> (&internal_data) != nullptr, ExcInternalError());
531  const InternalData &data = static_cast<const InternalData &> (internal_data);
532 
533  compute_fill (cell, face_no, subface_no, CellSimilarity::none,
534  data,
535  output_data.quadrature_points,
536  output_data.normal_vectors);
537 
538  // first compute Jacobian determinant, which is simply the product
539  // of the local lengths since the jacobian is diagonal
540  double J = 1.;
541  for (unsigned int d=0; d<dim; ++d)
543  J *= data.cell_extents[d];
544 
545  if (data.update_each & update_JxW_values)
546  {
547  // Here, cell->face(face_no)->n_children() would be the right
548  // choice, but unfortunately the current function is also called
549  // for faces without children (see tests/fe/mapping.cc). Add
550  // following switch to avoid diffs in tests/fe/mapping.OK
551  const unsigned int n_subfaces=
552  cell->face(face_no)->has_children() ?
553  cell->face(face_no)->n_children() :
555  for (unsigned int i=0; i<output_data.JxW_values.size(); ++i)
556  output_data.JxW_values[i] = J * quadrature.weight(i) / n_subfaces;
557  }
558 
559  if (data.update_each & update_boundary_forms)
560  for (unsigned int i=0; i<output_data.boundary_forms.size(); ++i)
561  output_data.boundary_forms[i] = J * output_data.normal_vectors[i];
562 
563  if (data.update_each & update_volume_elements)
564  {
565  J = data.cell_extents[0];
566  for (unsigned int d=1; d<dim; ++d)
567  J *= data.cell_extents[d];
568  data.volume_element = J;
569  }
570 
571  if (data.update_each & update_jacobians)
572  for (unsigned int i=0; i<output_data.jacobians.size(); ++i)
573  {
574  output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>();
575  for (unsigned int d=0; d<dim; ++d)
576  output_data.jacobians[i][d][d] = data.cell_extents[d];
577  }
578 
579  if (data.update_each & update_jacobian_grads)
580  for (unsigned int i=0; i<output_data.jacobian_grads.size(); ++i)
582 
583  if (data.update_each & update_jacobian_pushed_forward_grads)
584  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_grads.size(); ++i)
586 
587  if (data.update_each & update_jacobian_2nd_derivatives)
588  for (unsigned int i=0; i<output_data.jacobian_2nd_derivatives.size(); ++i)
590 
591  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
592  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_2nd_derivatives.size(); ++i)
594 
595  if (data.update_each & update_jacobian_3rd_derivatives)
596  for (unsigned int i=0; i<output_data.jacobian_3rd_derivatives.size(); ++i)
598 
599  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
600  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_3rd_derivatives.size(); ++i)
602 
603  if (data.update_each & update_inverse_jacobians)
604  for (unsigned int i=0; i<output_data.inverse_jacobians.size(); ++i)
605  {
607  for (unsigned int d=0; d<dim; ++d)
608  output_data.inverse_jacobians[i][d][d] = 1./data.cell_extents[d];
609  }
610 }
611 
612 
613 
614 
615 template <int dim, int spacedim>
616 void
618 transform (const ArrayView<const Tensor<1,dim> > &input,
619  const MappingType mapping_type,
620  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
621  const ArrayView<Tensor<1,spacedim> > &output) const
622 {
623  AssertDimension (input.size(), output.size());
624  Assert (dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
625  ExcInternalError());
626  const InternalData &data = static_cast<const InternalData &> (mapping_data);
627 
628  switch (mapping_type)
629  {
630  case mapping_covariant:
631  {
633  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
634 
635  for (unsigned int i=0; i<output.size(); ++i)
636  for (unsigned int d=0; d<dim; ++d)
637  output[i][d] = input[i][d]/data.cell_extents[d];
638  return;
639  }
640 
642  {
644  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
645 
646  for (unsigned int i=0; i<output.size(); ++i)
647  for (unsigned int d=0; d<dim; ++d)
648  output[i][d] = input[i][d]*data.cell_extents[d];
649  return;
650  }
651  case mapping_piola:
652  {
654  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
656  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
657 
658  for (unsigned int i=0; i<output.size(); ++i)
659  for (unsigned int d=0; d<dim; ++d)
660  output[i][d] = input[i][d] * data.cell_extents[d] / data.volume_element;
661  return;
662  }
663  default:
664  Assert(false, ExcNotImplemented());
665  }
666 }
667 
668 
669 
670 template <int dim, int spacedim>
671 void
674  const MappingType mapping_type,
675  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
676  const ArrayView<Tensor<2,spacedim> > &output) const
677 {
678  AssertDimension (input.size(), output.size());
679  Assert (dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
680  ExcInternalError());
681  const InternalData &data = static_cast<const InternalData &>(mapping_data);
682 
683  switch (mapping_type)
684  {
685  case mapping_covariant:
686  {
688  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
689 
690  for (unsigned int i=0; i<output.size(); ++i)
691  for (unsigned int d1=0; d1<dim; ++d1)
692  for (unsigned int d2=0; d2<dim; ++d2)
693  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2];
694  return;
695  }
696 
698  {
700  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
701 
702  for (unsigned int i=0; i<output.size(); ++i)
703  for (unsigned int d1=0; d1<dim; ++d1)
704  for (unsigned int d2=0; d2<dim; ++d2)
705  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
706  return;
707  }
708 
710  {
712  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
713 
714  for (unsigned int i=0; i<output.size(); ++i)
715  for (unsigned int d1=0; d1<dim; ++d1)
716  for (unsigned int d2=0; d2<dim; ++d2)
717  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2] / data.cell_extents[d1];
718  return;
719  }
720 
722  {
724  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
725 
726  for (unsigned int i=0; i<output.size(); ++i)
727  for (unsigned int d1=0; d1<dim; ++d1)
728  for (unsigned int d2=0; d2<dim; ++d2)
729  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] / data.cell_extents[d1];
730  return;
731  }
732 
733  case mapping_piola:
734  {
736  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
738  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
739 
740  for (unsigned int i=0; i<output.size(); ++i)
741  for (unsigned int d1=0; d1<dim; ++d1)
742  for (unsigned int d2=0; d2<dim; ++d2)
743  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2]
744  / data.volume_element;
745  return;
746  }
747 
749  {
751  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
753  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
754 
755  for (unsigned int i=0; i<output.size(); ++i)
756  for (unsigned int d1=0; d1<dim; ++d1)
757  for (unsigned int d2=0; d2<dim; ++d2)
758  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2]
759  / data.cell_extents[d1] / data.volume_element;
760  return;
761  }
762 
763  default:
764  Assert(false, ExcNotImplemented());
765  }
766 }
767 
768 
769 
770 
771 template <int dim, int spacedim>
772 void
774 transform (const ArrayView<const Tensor<2, dim> > &input,
775  const MappingType mapping_type,
776  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
777  const ArrayView<Tensor<2, spacedim> > &output) const
778 {
779 
780  AssertDimension (input.size(), output.size());
781  Assert (dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
782  ExcInternalError());
783  const InternalData &data = static_cast<const InternalData &>(mapping_data);
784 
785  switch (mapping_type)
786  {
787  case mapping_covariant:
788  {
790  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
791 
792  for (unsigned int i=0; i<output.size(); ++i)
793  for (unsigned int d1=0; d1<dim; ++d1)
794  for (unsigned int d2=0; d2<dim; ++d2)
795  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2];
796  return;
797  }
798 
800  {
802  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
803 
804  for (unsigned int i=0; i<output.size(); ++i)
805  for (unsigned int d1=0; d1<dim; ++d1)
806  for (unsigned int d2=0; d2<dim; ++d2)
807  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
808  return;
809  }
810 
812  {
814  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
815 
816  for (unsigned int i=0; i<output.size(); ++i)
817  for (unsigned int d1=0; d1<dim; ++d1)
818  for (unsigned int d2=0; d2<dim; ++d2)
819  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2] / data.cell_extents[d1];
820  return;
821  }
822 
824  {
826  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
827 
828  for (unsigned int i=0; i<output.size(); ++i)
829  for (unsigned int d1=0; d1<dim; ++d1)
830  for (unsigned int d2=0; d2<dim; ++d2)
831  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] / data.cell_extents[d1];
832  return;
833  }
834 
835  case mapping_piola:
836  {
838  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
840  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
841 
842  for (unsigned int i=0; i<output.size(); ++i)
843  for (unsigned int d1=0; d1<dim; ++d1)
844  for (unsigned int d2=0; d2<dim; ++d2)
845  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2]
846  / data.volume_element;
847  return;
848  }
849 
851  {
853  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
855  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
856 
857  for (unsigned int i=0; i<output.size(); ++i)
858  for (unsigned int d1=0; d1<dim; ++d1)
859  for (unsigned int d2=0; d2<dim; ++d2)
860  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2]
861  / data.cell_extents[d1] / data.volume_element;
862  return;
863  }
864 
865  default:
866  Assert(false, ExcNotImplemented());
867  }
868 
869 }
870 
871 
872 template <int dim, int spacedim>
873 void
876  const MappingType mapping_type,
877  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
878  const ArrayView<Tensor<3,spacedim> > &output) const
879 {
880 
881  AssertDimension (input.size(), output.size());
882  Assert (dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
883  ExcInternalError());
884  const InternalData &data = static_cast<const InternalData &>(mapping_data);
885 
886  switch (mapping_type)
887  {
889  {
891  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
892 
893  for (unsigned int q=0; q<output.size(); ++q)
894  for (unsigned int i=0; i<spacedim; ++i)
895  for (unsigned int j=0; j<spacedim; ++j)
896  for (unsigned int k=0; k<spacedim; ++k)
897  {
898 
899  output[q][i][j][k] = input[q][i][j][k] / data.cell_extents[j] / data.cell_extents[k];
900 
901  }
902  return;
903  }
904  default:
905  Assert(false, ExcNotImplemented());
906  }
907 
908 }
909 
910 template <int dim, int spacedim>
911 void
913 transform (const ArrayView<const Tensor<3,dim> > &input,
914  const MappingType mapping_type,
915  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
916  const ArrayView<Tensor<3,spacedim> > &output) const
917 {
918 
919  AssertDimension (input.size(), output.size());
920  Assert (dynamic_cast<const InternalData *>(&mapping_data) != nullptr,
921  ExcInternalError());
922  const InternalData &data = static_cast<const InternalData &>(mapping_data);
923 
924  switch (mapping_type)
925  {
927  {
929  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
931  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
932 
933  for (unsigned int q=0; q<output.size(); ++q)
934  for (unsigned int i=0; i<spacedim; ++i)
935  for (unsigned int j=0; j<spacedim; ++j)
936  for (unsigned int k=0; k<spacedim; ++k)
937  {
938  output[q][i][j][k] = input[q][i][j][k]
939  * data.cell_extents[i]
940  / data.cell_extents[j]
941  / data.cell_extents[k];
942  }
943  return;
944  }
945 
947  {
949  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
950 
951  for (unsigned int q=0; q<output.size(); ++q)
952  for (unsigned int i=0; i<spacedim; ++i)
953  for (unsigned int j=0; j<spacedim; ++j)
954  for (unsigned int k=0; k<spacedim; ++k)
955  {
956  output[q][i][j][k] = input[q][i][j][k]
957  / data.cell_extents[i]
958  / data.cell_extents[j]
959  / data.cell_extents[k];
960  }
961 
962  return;
963  }
964 
966  {
968  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
970  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
972  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
973 
974  for (unsigned int q=0; q<output.size(); ++q)
975  for (unsigned int i=0; i<spacedim; ++i)
976  for (unsigned int j=0; j<spacedim; ++j)
977  for (unsigned int k=0; k<spacedim; ++k)
978  {
979  output[q][i][j][k] = input[q][i][j][k]
980  * data.cell_extents[i]
981  / data.volume_element
982  / data.cell_extents[j]
983  / data.cell_extents[k];
984  }
985 
986  return;
987  }
988 
989  default:
990  Assert(false, ExcNotImplemented());
991  }
992 }
993 
994 
995 template <int dim, int spacedim>
998  const typename Triangulation<dim,spacedim>::cell_iterator &cell,
999  const Point<dim> &p) const
1000 {
1001  Tensor<1,dim> length;
1002  const Point<dim> start = cell->vertex(0);
1003  switch (dim)
1004  {
1005  case 1:
1006  length[0] = cell->vertex(1)(0) - start(0);
1007  break;
1008  case 2:
1009  length[0] = cell->vertex(1)(0) - start(0);
1010  length[1] = cell->vertex(2)(1) - start(1);
1011  break;
1012  case 3:
1013  length[0] = cell->vertex(1)(0) - start(0);
1014  length[1] = cell->vertex(2)(1) - start(1);
1015  length[2] = cell->vertex(4)(2) - start(2);
1016  break;
1017  default:
1018  Assert(false, ExcNotImplemented());
1019  }
1020 
1021  Point<dim> p_real = cell->vertex(0);
1022  for (unsigned int d=0; d<dim; ++d)
1023  p_real(d) += length[d]*p(d);
1024 
1025  return p_real;
1026 }
1027 
1028 
1029 
1030 template <int dim, int spacedim>
1031 Point<dim>
1033  const typename Triangulation<dim,spacedim>::cell_iterator &cell,
1034  const Point<spacedim> &p) const
1035 {
1036 
1037  if (dim != spacedim)
1038  Assert(false, ExcNotImplemented());
1039  const Point<dim> &start = cell->vertex(0);
1040  Point<dim> real = p;
1041  real -= start;
1042 
1043  switch (dim)
1044  {
1045  case 1:
1046  real(0) /= cell->vertex(1)(0) - start(0);
1047  break;
1048  case 2:
1049  real(0) /= cell->vertex(1)(0) - start(0);
1050  real(1) /= cell->vertex(2)(1) - start(1);
1051  break;
1052  case 3:
1053  real(0) /= cell->vertex(1)(0) - start(0);
1054  real(1) /= cell->vertex(2)(1) - start(1);
1055  real(2) /= cell->vertex(4)(2) - start(2);
1056  break;
1057  default:
1058  Assert(false, ExcNotImplemented());
1059  }
1060  return real;
1061 }
1062 
1063 
1064 template <int dim, int spacedim>
1065 std::unique_ptr<Mapping<dim, spacedim> >
1067 {
1068  return std_cxx14::make_unique<MappingCartesian<dim,spacedim>>(*this);
1069 }
1070 
1071 
1072 //---------------------------------------------------------------------------
1073 // explicit instantiations
1074 #include "mapping_cartesian.inst"
1075 
1076 
1077 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
std::vector< Tensor< 1, spacedim > > normal_vectors
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
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
Contravariant transformation.
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
std::vector< Tensor< 4, spacedim > > jacobian_pushed_forward_2nd_derivatives
MappingType
Definition: mapping.h:51
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
Volume element.
Outer normal vector, not normalized.
std::vector< Point< dim > > quadrature_points
Determinant of the Jacobian.
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
std::vector< Tensor< 5, spacedim > > jacobian_pushed_forward_3rd_derivatives
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static DataSetDescriptor cell()
Definition: qprojector.h:348
static Quadrature< dim > project_to_all_subfaces(const SubQuadrature &quadrature)
virtual std::size_t memory_consumption() const override
std::vector< DerivativeForm< 2, dim, spacedim > > jacobian_grads
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
#define Assert(cond, exc)
Definition: exceptions.h:1142
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:46
Gradient of volume element.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< Tensor< 3, spacedim > > jacobian_pushed_forward_grads
std::vector< DerivativeForm< 3, dim, spacedim > > jacobian_2nd_derivatives
static const unsigned int invalid_face_number
std::vector< Point< spacedim > > quadrature_points
static Quadrature< dim > project_to_all_faces(const SubQuadrature &quadrature)
std::vector< DerivativeForm< 1, spacedim, dim > > inverse_jacobians
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:53
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
static ::ExceptionBase & ExcNotImplemented()
std::vector< Tensor< 1, spacedim > > boundary_forms
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
Covariant transformation.
double weight(const unsigned int i) const
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
UpdateFlags update_each
Definition: mapping.h:562
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:1174