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