Reference documentation for deal.II version 8.5.1
mapping_cartesian.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2016 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/signaling_nan.h>
21 #include <deal.II/base/memory_consumption.h>
22 #include <deal.II/lac/full_matrix.h>
23 #include <deal.II/grid/tria.h>
24 #include <deal.II/grid/tria_iterator.h>
25 #include <deal.II/dofs/dof_accessor.h>
26 #include <deal.II/fe/mapping_cartesian.h>
27 #include <deal.II/fe/fe_values.h>
28 
29 #include <cmath>
30 #include <algorithm>
31 
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 
36 template<int dim, int spacedim>
38 
39 
40 
41 template<int dim, int spacedim>
43  :
44  cell_extents (numbers::signaling_nan<Tensor<1,dim> >()),
45  volume_element (numbers::signaling_nan<double>()),
46  quadrature_points (q.get_points ())
47 {}
48 
49 
50 
51 template<int dim, int spacedim>
52 std::size_t
54 {
57  MemoryConsumption::memory_consumption (volume_element) +
58  MemoryConsumption::memory_consumption (quadrature_points));
59 }
60 
61 
62 
63 template <int dim, int spacedim>
64 bool
66 {
67  return true;
68 }
69 
70 
71 
72 template<int dim, int spacedim>
75 {
76  // this mapping is pretty simple in that it can basically compute
77  // every piece of information wanted by FEValues without requiring
78  // computing any other quantities. boundary forms are one exception
79  // since they can be computed from the normal vectors without much
80  // further ado
81  UpdateFlags out = in;
82  if (out & update_boundary_forms)
83  out |= update_normal_vectors;
84 
85  return out;
86 }
87 
88 
89 
90 template<int dim, int spacedim>
93  const Quadrature<dim> &q) const
94 {
95  InternalData *data = new InternalData (q);
96 
97  // store the flags in the internal data object so we can access them
98  // in fill_fe_*_values(). use the transitive hull of the required
99  // flags
100  data->update_each = requires_update_flags(update_flags);
101 
102  return data;
103 }
104 
105 
106 
107 template<int dim, int spacedim>
110  const Quadrature<dim-1>& quadrature) const
111 {
112  InternalData *data
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 data;
125 }
126 
127 
128 
129 template<int dim, int spacedim>
132  const Quadrature<dim-1> &quadrature) const
133 {
134  InternalData *data
136 
137  // verify that we have computed the transitive hull of the required
138  // flags and that FEValues has faithfully passed them on to us
139  Assert (update_flags == requires_update_flags (update_flags),
140  ExcInternalError());
141 
142  // store the flags in the internal data object so we can access them
143  // in fill_fe_*_values()
144  data->update_each = update_flags;
145 
146  return data;
147 }
148 
149 
150 
151 
152 template<int dim, int spacedim>
153 void
155  const unsigned int face_no,
156  const unsigned int sub_no,
157  const CellSimilarity::Similarity cell_similarity,
158  const InternalData &data,
159  std::vector<Point<dim> > &quadrature_points,
160  std::vector<Tensor<1,dim> > &normal_vectors) const
161 {
162  const UpdateFlags update_flags = data.update_each;
163 
164  // some more sanity checks
165  if (face_no != invalid_face_number)
166  {
167  // Add 1 on both sides of
168  // assertion to avoid compiler
169  // warning about testing
170  // unsigned int < 0 in 1d.
173 
174  // We would like to check for
175  // sub_no < cell->face(face_no)->n_children(),
176  // but unfortunately the current
177  // function is also called for
178  // faces without children (see
179  // tests/fe/mapping.cc). Therefore,
180  // we must use following workaround
181  // of two separate assertions
182  Assert ((sub_no == invalid_face_number) ||
183  cell->face(face_no)->has_children() ||
185  ExcIndexRange (sub_no, 0,
187  Assert ((sub_no == invalid_face_number) ||
188  !cell->face(face_no)->has_children() ||
189  (sub_no < cell->face(face_no)->n_children()),
190  ExcIndexRange (sub_no, 0, cell->face(face_no)->n_children()));
191  }
192  else
193  // invalid face number, so
194  // subface should be invalid as
195  // well
197 
198  // let @p{start} be the origin of a
199  // local coordinate system. it is
200  // chosen as the (lower) left
201  // vertex
202  const Point<dim> start = cell->vertex(0);
203 
204  // Compute start point and sizes
205  // along axes. Strange vertex
206  // numbering makes this complicated
207  // again.
208  if (cell_similarity != CellSimilarity::translation)
209  {
210  switch (dim)
211  {
212  case 1:
213  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
214  break;
215  case 2:
216  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
217  data.cell_extents[1] = cell->vertex(2)(1) - start(1);
218  break;
219  case 3:
220  data.cell_extents[0] = cell->vertex(1)(0) - start(0);
221  data.cell_extents[1] = cell->vertex(2)(1) - start(1);
222  data.cell_extents[2] = cell->vertex(4)(2) - start(2);
223  break;
224  default:
225  Assert(false, ExcNotImplemented());
226  }
227  }
228 
229 
230  // transform quadrature point. this
231  // is obtained simply by scaling
232  // unit coordinates with lengths in
233  // each direction
234  if (update_flags & update_quadrature_points)
235  {
236  const typename QProjector<dim>::DataSetDescriptor offset
237  = (face_no == invalid_face_number
238  ?
240  :
241  (sub_no == invalid_face_number
242  ?
243  // called from FEFaceValues
245  cell->face_orientation(face_no),
246  cell->face_flip(face_no),
247  cell->face_rotation(face_no),
248  quadrature_points.size())
249  :
250  // called from FESubfaceValues
252  cell->face_orientation(face_no),
253  cell->face_flip(face_no),
254  cell->face_rotation(face_no),
255  quadrature_points.size(),
256  cell->subface_case(face_no))
257  ));
258 
259  for (unsigned int i=0; i<quadrature_points.size(); ++i)
260  {
261  quadrature_points[i] = start;
262  for (unsigned int d=0; d<dim; ++d)
263  quadrature_points[i](d) += data.cell_extents[d] *
264  data.quadrature_points[i+offset](d);
265  }
266  }
267 
268 
269  // compute normal vectors. since
270  // cells are aligned to coordinate
271  // axes, they are simply vectors
272  // with exactly one entry equal to
273  // 1 or -1. Furthermore, all
274  // normals on a face have the same
275  // value
276  if (update_flags & update_normal_vectors)
277  {
279  ExcInternalError());
280 
281  switch (dim)
282  {
283  case 1:
284  {
285  static const Point<dim>
287  = { Point<dim>(-1.),
288  Point<dim>( 1.)
289  };
290  std::fill (normal_vectors.begin(),
291  normal_vectors.end(),
292  normals[face_no]);
293  break;
294  }
295 
296  case 2:
297  {
298  static const Point<dim>
300  = { Point<dim>(-1, 0),
301  Point<dim>( 1, 0),
302  Point<dim>( 0,-1),
303  Point<dim>( 0, 1)
304  };
305  std::fill (normal_vectors.begin(),
306  normal_vectors.end(),
307  normals[face_no]);
308  break;
309  }
310 
311  case 3:
312  {
313  static const Point<dim>
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  };
322  std::fill (normal_vectors.begin(),
323  normal_vectors.end(),
324  normals[face_no]);
325  break;
326  }
327 
328  default:
329  Assert (false, ExcNotImplemented());
330  }
331  }
332 }
333 
334 
335 
336 template<int dim, int spacedim>
340  const CellSimilarity::Similarity cell_similarity,
341  const Quadrature<dim> &quadrature,
342  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
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) != 0, ExcInternalError());
348  const InternalData &data = static_cast<const InternalData &> (internal_data);
349 
350  std::vector<Tensor<1,dim> > dummy;
351 
352  compute_fill (cell, invalid_face_number, invalid_face_number, cell_similarity,
353  data,
354  output_data.quadrature_points,
355  dummy);
356 
357  // compute Jacobian determinant. all values are equal and are the
358  // product of the local lengths in each coordinate direction
359  if (data.update_each & (update_JxW_values | update_volume_elements))
360  if (cell_similarity != CellSimilarity::translation)
361  {
362  double J = data.cell_extents[0];
363  for (unsigned int d=1; d<dim; ++d)
364  J *= data.cell_extents[d];
365  data.volume_element = J;
366  if (data.update_each & update_JxW_values)
367  for (unsigned int i=0; i<output_data.JxW_values.size(); ++i)
368  output_data.JxW_values[i] = J * quadrature.weight(i);
369  }
370  // "compute" Jacobian at the quadrature points, which are all the
371  // same
372  if (data.update_each & update_jacobians)
373  if (cell_similarity != CellSimilarity::translation)
374  for (unsigned int i=0; i<output_data.jacobians.size(); ++i)
375  {
376  output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>();
377  for (unsigned int j=0; j<dim; ++j)
378  output_data.jacobians[i][j][j] = data.cell_extents[j];
379  }
380  // "compute" the derivative of the Jacobian at the quadrature
381  // points, which are all zero of course
382  if (data.update_each & update_jacobian_grads)
383  if (cell_similarity != CellSimilarity::translation)
384  for (unsigned int i=0; i<output_data.jacobian_grads.size(); ++i)
386 
387  if (data.update_each & update_jacobian_pushed_forward_grads)
388  if (cell_similarity != CellSimilarity::translation)
389  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_grads.size(); ++i)
391 
392  // "compute" the hessian of the Jacobian at the quadrature points,
393  // which are all also zero of course
394  if (data.update_each & update_jacobian_2nd_derivatives)
395  if (cell_similarity != CellSimilarity::translation)
396  for (unsigned int i=0; i<output_data.jacobian_2nd_derivatives.size(); ++i)
398 
399  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
400  if (cell_similarity != CellSimilarity::translation)
401  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_2nd_derivatives.size(); ++i)
403 
404  if (data.update_each & update_jacobian_3rd_derivatives)
405  if (cell_similarity != CellSimilarity::translation)
406  for (unsigned int i=0; i<output_data.jacobian_3rd_derivatives.size(); ++i)
408 
409  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
410  if (cell_similarity != CellSimilarity::translation)
411  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_3rd_derivatives.size(); ++i)
413 
414  // "compute" inverse Jacobian at the quadrature points, which are
415  // all the same
416  if (data.update_each & update_inverse_jacobians)
417  if (cell_similarity != CellSimilarity::translation)
418  for (unsigned int i=0; i<output_data.inverse_jacobians.size(); ++i)
419  {
420  output_data.inverse_jacobians[i] = Tensor<2,dim>();
421  for (unsigned int j=0; j<dim; ++j)
422  output_data.inverse_jacobians[j][j]=1./data.cell_extents[j];
423  }
424 
425  return cell_similarity;
426 }
427 
428 
429 
430 template<int dim, int spacedim>
431 void
434  const unsigned int face_no,
435  const Quadrature<dim-1> &quadrature,
436  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
438 {
439  // convert data object to internal
440  // data for this class. fails with
441  // an exception if that is not
442  // possible
443  Assert (dynamic_cast<const InternalData *> (&internal_data) != 0,
444  ExcInternalError());
445  const InternalData &data = static_cast<const InternalData &> (internal_data);
446 
447  compute_fill (cell, face_no, invalid_face_number,
449  data,
450  output_data.quadrature_points,
451  output_data.normal_vectors);
452 
453  // first compute Jacobian determinant, which is simply the product
454  // of the local lengths since the jacobian is diagonal
455  double J = 1.;
456  for (unsigned int d=0; d<dim; ++d)
458  J *= data.cell_extents[d];
459 
460  if (data.update_each & update_JxW_values)
461  for (unsigned int i=0; i<output_data.JxW_values.size(); ++i)
462  output_data.JxW_values[i] = J * quadrature.weight(i);
463 
464  if (data.update_each & update_boundary_forms)
465  for (unsigned int i=0; i<output_data.boundary_forms.size(); ++i)
466  output_data.boundary_forms[i] = J * output_data.normal_vectors[i];
467 
468  if (data.update_each & update_volume_elements)
469  {
470  J = data.cell_extents[0];
471  for (unsigned int d=1; d<dim; ++d)
472  J *= data.cell_extents[d];
473  data.volume_element = J;
474  }
475 
476  if (data.update_each & update_jacobians)
477  for (unsigned int i=0; i<output_data.jacobians.size(); ++i)
478  {
479  output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>();
480  for (unsigned int d=0; d<dim; ++d)
481  output_data.jacobians[i][d][d] = data.cell_extents[d];
482  }
483 
484  if (data.update_each & update_jacobian_grads)
485  for (unsigned int i=0; i<output_data.jacobian_grads.size(); ++i)
487 
488  if (data.update_each & update_jacobian_pushed_forward_grads)
489  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_grads.size(); ++i)
491 
492  if (data.update_each & update_jacobian_2nd_derivatives)
493  for (unsigned int i=0; i<output_data.jacobian_2nd_derivatives.size(); ++i)
495 
496  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
497  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_2nd_derivatives.size(); ++i)
499 
500  if (data.update_each & update_jacobian_3rd_derivatives)
501  for (unsigned int i=0; i<output_data.jacobian_3rd_derivatives.size(); ++i)
503 
504  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
505  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_3rd_derivatives.size(); ++i)
507 
508  if (data.update_each & update_inverse_jacobians)
509  for (unsigned int i=0; i<output_data.inverse_jacobians.size(); ++i)
510  {
512  for (unsigned int d=0; d<dim; ++d)
513  output_data.inverse_jacobians[i][d][d] = 1./data.cell_extents[d];
514  }
515 }
516 
517 
518 
519 template<int dim, int spacedim>
520 void
523  const unsigned int face_no,
524  const unsigned int subface_no,
525  const Quadrature<dim-1> &quadrature,
526  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
528 {
529  // convert data object to internal data for this class. fails with
530  // an exception if that is not possible
531  Assert (dynamic_cast<const InternalData *> (&internal_data) != 0, ExcInternalError());
532  const InternalData &data = static_cast<const InternalData &> (internal_data);
533 
534  compute_fill (cell, face_no, subface_no, CellSimilarity::none,
535  data,
536  output_data.quadrature_points,
537  output_data.normal_vectors);
538 
539  // first compute Jacobian determinant, which is simply the product
540  // of the local lengths since the jacobian is diagonal
541  double J = 1.;
542  for (unsigned int d=0; d<dim; ++d)
544  J *= data.cell_extents[d];
545 
546  if (data.update_each & update_JxW_values)
547  {
548  // Here, cell->face(face_no)->n_children() would be the right
549  // choice, but unfortunately the current function is also called
550  // for faces without children (see tests/fe/mapping.cc). Add
551  // following switch to avoid diffs in tests/fe/mapping.OK
552  const unsigned int n_subfaces=
553  cell->face(face_no)->has_children() ?
554  cell->face(face_no)->n_children() :
556  for (unsigned int i=0; i<output_data.JxW_values.size(); ++i)
557  output_data.JxW_values[i] = J * quadrature.weight(i) / n_subfaces;
558  }
559 
560  if (data.update_each & update_boundary_forms)
561  for (unsigned int i=0; i<output_data.boundary_forms.size(); ++i)
562  output_data.boundary_forms[i] = J * output_data.normal_vectors[i];
563 
564  if (data.update_each & update_volume_elements)
565  {
566  J = data.cell_extents[0];
567  for (unsigned int d=1; d<dim; ++d)
568  J *= data.cell_extents[d];
569  data.volume_element = J;
570  }
571 
572  if (data.update_each & update_jacobians)
573  for (unsigned int i=0; i<output_data.jacobians.size(); ++i)
574  {
575  output_data.jacobians[i] = DerivativeForm<1,dim,spacedim>();
576  for (unsigned int d=0; d<dim; ++d)
577  output_data.jacobians[i][d][d] = data.cell_extents[d];
578  }
579 
580  if (data.update_each & update_jacobian_grads)
581  for (unsigned int i=0; i<output_data.jacobian_grads.size(); ++i)
583 
584  if (data.update_each & update_jacobian_pushed_forward_grads)
585  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_grads.size(); ++i)
587 
588  if (data.update_each & update_jacobian_2nd_derivatives)
589  for (unsigned int i=0; i<output_data.jacobian_2nd_derivatives.size(); ++i)
591 
592  if (data.update_each & update_jacobian_pushed_forward_2nd_derivatives)
593  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_2nd_derivatives.size(); ++i)
595 
596  if (data.update_each & update_jacobian_3rd_derivatives)
597  for (unsigned int i=0; i<output_data.jacobian_3rd_derivatives.size(); ++i)
599 
600  if (data.update_each & update_jacobian_pushed_forward_3rd_derivatives)
601  for (unsigned int i=0; i<output_data.jacobian_pushed_forward_3rd_derivatives.size(); ++i)
603 
604  if (data.update_each & update_inverse_jacobians)
605  for (unsigned int i=0; i<output_data.inverse_jacobians.size(); ++i)
606  {
608  for (unsigned int d=0; d<dim; ++d)
609  output_data.inverse_jacobians[i][d][d] = 1./data.cell_extents[d];
610  }
611 }
612 
613 
614 
615 
616 template<int dim, int spacedim>
617 void
619 transform (const ArrayView<const Tensor<1,dim> > &input,
620  const MappingType mapping_type,
621  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
622  const ArrayView<Tensor<1,spacedim> > &output) const
623 {
624  AssertDimension (input.size(), output.size());
625  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
626  ExcInternalError());
627  const InternalData &data = static_cast<const InternalData &> (mapping_data);
628 
629  switch (mapping_type)
630  {
631  case mapping_covariant:
632  {
634  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
635 
636  for (unsigned int i=0; i<output.size(); ++i)
637  for (unsigned int d=0; d<dim; ++d)
638  output[i][d] = input[i][d]/data.cell_extents[d];
639  return;
640  }
641 
643  {
645  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
646 
647  for (unsigned int i=0; i<output.size(); ++i)
648  for (unsigned int d=0; d<dim; ++d)
649  output[i][d] = input[i][d]*data.cell_extents[d];
650  return;
651  }
652  case mapping_piola:
653  {
655  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
657  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
658 
659  for (unsigned int i=0; i<output.size(); ++i)
660  for (unsigned int d=0; d<dim; ++d)
661  output[i][d] = input[i][d] * data.cell_extents[d] / data.volume_element;
662  return;
663  }
664  default:
665  Assert(false, ExcNotImplemented());
666  }
667 }
668 
669 
670 
671 template<int dim, int spacedim>
672 void
675  const MappingType mapping_type,
676  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
677  const ArrayView<Tensor<2,spacedim> > &output) const
678 {
679  AssertDimension (input.size(), output.size());
680  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
681  ExcInternalError());
682  const InternalData &data = static_cast<const InternalData &>(mapping_data);
683 
684  switch (mapping_type)
685  {
686  case mapping_covariant:
687  {
689  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
690 
691  for (unsigned int i=0; i<output.size(); ++i)
692  for (unsigned int d1=0; d1<dim; ++d1)
693  for (unsigned int d2=0; d2<dim; ++d2)
694  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2];
695  return;
696  }
697 
699  {
701  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
702 
703  for (unsigned int i=0; i<output.size(); ++i)
704  for (unsigned int d1=0; d1<dim; ++d1)
705  for (unsigned int d2=0; d2<dim; ++d2)
706  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
707  return;
708  }
709 
711  {
713  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
714 
715  for (unsigned int i=0; i<output.size(); ++i)
716  for (unsigned int d1=0; d1<dim; ++d1)
717  for (unsigned int d2=0; d2<dim; ++d2)
718  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2] / data.cell_extents[d1];
719  return;
720  }
721 
723  {
725  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
726 
727  for (unsigned int i=0; i<output.size(); ++i)
728  for (unsigned int d1=0; d1<dim; ++d1)
729  for (unsigned int d2=0; d2<dim; ++d2)
730  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] / data.cell_extents[d1];
731  return;
732  }
733 
734  case mapping_piola:
735  {
737  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
739  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
740 
741  for (unsigned int i=0; i<output.size(); ++i)
742  for (unsigned int d1=0; d1<dim; ++d1)
743  for (unsigned int d2=0; d2<dim; ++d2)
744  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2]
745  / data.volume_element;
746  return;
747  }
748 
750  {
752  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
754  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
755 
756  for (unsigned int i=0; i<output.size(); ++i)
757  for (unsigned int d1=0; d1<dim; ++d1)
758  for (unsigned int d2=0; d2<dim; ++d2)
759  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2]
760  / data.cell_extents[d1] / data.volume_element;
761  return;
762  }
763 
764  default:
765  Assert(false, ExcNotImplemented());
766  }
767 }
768 
769 
770 
771 
772 template<int dim, int spacedim>
773 void
775 transform (const ArrayView<const Tensor<2, dim> > &input,
776  const MappingType mapping_type,
777  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
778  const ArrayView<Tensor<2, spacedim> > &output) const
779 {
780 
781  AssertDimension (input.size(), output.size());
782  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
783  ExcInternalError());
784  const InternalData &data = static_cast<const InternalData &>(mapping_data);
785 
786  switch (mapping_type)
787  {
788  case mapping_covariant:
789  {
791  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
792 
793  for (unsigned int i=0; i<output.size(); ++i)
794  for (unsigned int d1=0; d1<dim; ++d1)
795  for (unsigned int d2=0; d2<dim; ++d2)
796  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2];
797  return;
798  }
799 
801  {
803  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
804 
805  for (unsigned int i=0; i<output.size(); ++i)
806  for (unsigned int d1=0; d1<dim; ++d1)
807  for (unsigned int d2=0; d2<dim; ++d2)
808  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2];
809  return;
810  }
811 
813  {
815  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
816 
817  for (unsigned int i=0; i<output.size(); ++i)
818  for (unsigned int d1=0; d1<dim; ++d1)
819  for (unsigned int d2=0; d2<dim; ++d2)
820  output[i][d1][d2] = input[i][d1][d2] / data.cell_extents[d2] / data.cell_extents[d1];
821  return;
822  }
823 
825  {
827  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
828 
829  for (unsigned int i=0; i<output.size(); ++i)
830  for (unsigned int d1=0; d1<dim; ++d1)
831  for (unsigned int d2=0; d2<dim; ++d2)
832  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2] / data.cell_extents[d1];
833  return;
834  }
835 
836  case mapping_piola:
837  {
839  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
841  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
842 
843  for (unsigned int i=0; i<output.size(); ++i)
844  for (unsigned int d1=0; d1<dim; ++d1)
845  for (unsigned int d2=0; d2<dim; ++d2)
846  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2]
847  / data.volume_element;
848  return;
849  }
850 
852  {
854  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
856  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
857 
858  for (unsigned int i=0; i<output.size(); ++i)
859  for (unsigned int d1=0; d1<dim; ++d1)
860  for (unsigned int d2=0; d2<dim; ++d2)
861  output[i][d1][d2] = input[i][d1][d2] * data.cell_extents[d2]
862  / data.cell_extents[d1] / data.volume_element;
863  return;
864  }
865 
866  default:
867  Assert(false, ExcNotImplemented());
868  }
869 
870 }
871 
872 
873 template<int dim, int spacedim>
874 void
877  const MappingType mapping_type,
878  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
879  const ArrayView<Tensor<3,spacedim> > &output) const
880 {
881 
882  AssertDimension (input.size(), output.size());
883  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
884  ExcInternalError());
885  const InternalData &data = static_cast<const InternalData &>(mapping_data);
886 
887  switch (mapping_type)
888  {
890  {
892  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
893 
894  for (unsigned int q=0; q<output.size(); ++q)
895  for (unsigned int i=0; i<spacedim; ++i)
896  for (unsigned int j=0; j<spacedim; ++j)
897  for (unsigned int k=0; k<spacedim; ++k)
898  {
899 
900  output[q][i][j][k] = input[q][i][j][k] / data.cell_extents[j] / data.cell_extents[k];
901 
902  }
903  return;
904  }
905  default:
906  Assert(false, ExcNotImplemented());
907  }
908 
909 }
910 
911 template<int dim, int spacedim>
912 void
914 transform (const ArrayView<const Tensor<3,dim> > &input,
915  const MappingType mapping_type,
916  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_data,
917  const ArrayView<Tensor<3,spacedim> > &output) const
918 {
919 
920  AssertDimension (input.size(), output.size());
921  Assert (dynamic_cast<const InternalData *>(&mapping_data) != 0,
922  ExcInternalError());
923  const InternalData &data = static_cast<const InternalData &>(mapping_data);
924 
925  switch (mapping_type)
926  {
928  {
930  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
932  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
933 
934  for (unsigned int q=0; q<output.size(); ++q)
935  for (unsigned int i=0; i<spacedim; ++i)
936  for (unsigned int j=0; j<spacedim; ++j)
937  for (unsigned int k=0; k<spacedim; ++k)
938  {
939  output[q][i][j][k] = input[q][i][j][k]
940  * data.cell_extents[i]
941  / data.cell_extents[j]
942  / data.cell_extents[k];
943  }
944  return;
945  }
946 
948  {
950  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
951 
952  for (unsigned int q=0; q<output.size(); ++q)
953  for (unsigned int i=0; i<spacedim; ++i)
954  for (unsigned int j=0; j<spacedim; ++j)
955  for (unsigned int k=0; k<spacedim; ++k)
956  {
957  output[q][i][j][k] = input[q][i][j][k]
958  / data.cell_extents[i]
959  / data.cell_extents[j]
960  / data.cell_extents[k];
961  }
962 
963  return;
964  }
965 
967  {
969  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_covariant_transformation"));
971  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_contravariant_transformation"));
973  typename FEValuesBase<dim>::ExcAccessToUninitializedField("update_volume_elements"));
974 
975  for (unsigned int q=0; q<output.size(); ++q)
976  for (unsigned int i=0; i<spacedim; ++i)
977  for (unsigned int j=0; j<spacedim; ++j)
978  for (unsigned int k=0; k<spacedim; ++k)
979  {
980  output[q][i][j][k] = input[q][i][j][k]
981  * data.cell_extents[i]
982  / data.volume_element
983  / data.cell_extents[j]
984  / data.cell_extents[k];
985  }
986 
987  return;
988  }
989 
990  default:
991  Assert(false, ExcNotImplemented());
992  }
993 }
994 
995 
996 template <int dim, int spacedim>
999  const typename Triangulation<dim,spacedim>::cell_iterator &cell,
1000  const Point<dim> &p) const
1001 {
1002  Tensor<1,dim> length;
1003  const Point<dim> start = cell->vertex(0);
1004  switch (dim)
1005  {
1006  case 1:
1007  length[0] = cell->vertex(1)(0) - start(0);
1008  break;
1009  case 2:
1010  length[0] = cell->vertex(1)(0) - start(0);
1011  length[1] = cell->vertex(2)(1) - start(1);
1012  break;
1013  case 3:
1014  length[0] = cell->vertex(1)(0) - start(0);
1015  length[1] = cell->vertex(2)(1) - start(1);
1016  length[2] = cell->vertex(4)(2) - start(2);
1017  break;
1018  default:
1019  Assert(false, ExcNotImplemented());
1020  }
1021 
1022  Point<dim> p_real = cell->vertex(0);
1023  for (unsigned int d=0; d<dim; ++d)
1024  p_real(d) += length[d]*p(d);
1025 
1026  return p_real;
1027 }
1028 
1029 
1030 
1031 template<int dim, int spacedim>
1032 Point<dim>
1034  const typename Triangulation<dim,spacedim>::cell_iterator &cell,
1035  const Point<spacedim> &p) const
1036 {
1037 
1038  if (dim != spacedim)
1039  Assert(false, ExcNotImplemented());
1040  const Point<dim> &start = cell->vertex(0);
1041  Point<dim> real = p;
1042  real -= start;
1043 
1044  switch (dim)
1045  {
1046  case 1:
1047  real(0) /= cell->vertex(1)(0) - start(0);
1048  break;
1049  case 2:
1050  real(0) /= cell->vertex(1)(0) - start(0);
1051  real(1) /= cell->vertex(2)(1) - start(1);
1052  break;
1053  case 3:
1054  real(0) /= cell->vertex(1)(0) - start(0);
1055  real(1) /= cell->vertex(2)(1) - start(1);
1056  real(2) /= cell->vertex(4)(2) - start(2);
1057  break;
1058  default:
1059  Assert(false, ExcNotImplemented());
1060  }
1061  return real;
1062 }
1063 
1064 
1065 template<int dim, int spacedim>
1068 {
1069  return new MappingCartesian<dim, spacedim>(*this);
1070 }
1071 
1072 
1073 //---------------------------------------------------------------------------
1074 // explicit instantiations
1075 #include "mapping_cartesian.inst"
1076 
1077 
1078 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:1146
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 ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
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:313
UpdateFlags
std::vector< DerivativeForm< 1, dim, spacedim > > jacobians
Abstract base class for mapping classes.
Definition: dof_tools.h:46
std::vector< Tensor< 1, spacedim > > normal_vectors
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
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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)
Definition: mpi.h:41
std::vector< DerivativeForm< 4, dim, spacedim > > jacobian_3rd_derivatives
std::vector< Point< spacedim > > quadrature_points
Normal vectors.
static ::ExceptionBase & ExcNotImplemented()
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 ::ExceptionBase & ExcInternalError()
static DataSetDescriptor face(const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
Definition: quadrature.cc:1081