Reference documentation for deal.II version 8.4.2
mapping_q1_eulerian.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 #include <deal.II/base/std_cxx11/array.h>
17 #include <deal.II/fe/mapping_q1_eulerian.h>
18 #include <deal.II/lac/vector.h>
19 #include <deal.II/lac/petsc_vector.h>
20 #include <deal.II/lac/trilinos_vector.h>
21 #include <deal.II/lac/trilinos_block_vector.h>
22 #include <deal.II/lac/trilinos_parallel_block_vector.h>
23 #include <deal.II/grid/tria_iterator.h>
24 #include <deal.II/dofs/dof_handler.h>
25 #include <deal.II/dofs/dof_accessor.h>
26 #include <deal.II/fe/fe.h>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 
31 template <int dim, class EulerVectorType, int spacedim>
33 MappingQ1Eulerian (const EulerVectorType &euler_transform_vectors,
34  const DoFHandler<dim,spacedim> &shiftmap_dof_handler)
35  :
36  MappingQGeneric<dim,spacedim>(1),
37  euler_transform_vectors(&euler_transform_vectors),
38  shiftmap_dof_handler(&shiftmap_dof_handler)
39 {}
40 
41 
42 
43 template <int dim, class EulerVectorType, int spacedim>
44 std_cxx11::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
47 (const typename Triangulation<dim,spacedim>::cell_iterator &cell) const
48 {
49  std_cxx11::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell> vertices;
50  // The assertions can not be in the constructor, since this would
51  // require to call dof_handler.distribute_dofs(fe) *before* the mapping
52  // object is constructed, which is not necessarily what we want.
53 
54  //TODO: Only one of these two assertions should be relevant
55  AssertDimension (spacedim, shiftmap_dof_handler->get_fe().n_dofs_per_vertex());
56  AssertDimension (shiftmap_dof_handler->get_fe().n_components(), spacedim);
57 
59 
60  // cast the Triangulation<dim>::cell_iterator into a
61  // DoFHandler<dim>::cell_iterator which is necessary for access to
62  // DoFCellAccessor::get_dof_values()
64 
65  // We require the cell to be active since we can only then get nodal
66  // values for the shifts
67  Assert (dof_cell->active() == true, ExcInactiveCell());
68 
69  // now get the values of the shift vectors at the vertices
70  Vector<double> mapping_values (shiftmap_dof_handler->get_fe().dofs_per_cell);
71  dof_cell->get_dof_values (*euler_transform_vectors, mapping_values);
72 
73  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
74  {
75  Point<spacedim> shift_vector;
76 
77  // pick out the value of the shift vector at the present
78  // vertex. since vertex dofs are always numbered first, we can
79  // access them easily
80  for (unsigned int j=0; j<spacedim; ++j)
81  shift_vector[j] = mapping_values(i*spacedim+j);
82 
83  // compute new support point by old (reference) value and added
84  // shift
85  vertices[i] = cell->vertex(i) + shift_vector;
86  }
87  return vertices;
88 }
89 
90 
91 
92 template<int dim, class EulerVectorType, int spacedim>
93 std::vector<Point<spacedim> >
96 {
97  const std_cxx11::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
98  vertices = this->get_vertices(cell);
99 
100  std::vector<Point<spacedim> > a(GeometryInfo<dim>::vertices_per_cell);
101  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
102  a[i] = vertices[i];
103 
104  return a;
105 }
106 
107 
108 
109 
110 
111 template <int dim, class EulerVectorType, int spacedim>
114 {
116 }
117 
118 
119 
120 template<int dim, class EulerVectorType, int spacedim>
121 CellSimilarity::Similarity
124  const CellSimilarity::Similarity ,
125  const Quadrature<dim> &quadrature,
126  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
128 {
129  // call the function of the base class, but ignoring
130  // any potentially detected cell similarity between
131  // the current and the previous cell
133  CellSimilarity::invalid_next_cell,
134  quadrature,
135  internal_data,
136  output_data);
137  // also return the updated flag since any detected
138  // similarity wasn't based on the mapped field, but
139  // the original vertices which are meaningless
140  return CellSimilarity::invalid_next_cell;
141 }
142 
143 
144 
145 // explicit instantiations
146 #include "mapping_q1_eulerian.inst"
147 
148 
149 DEAL_II_NAMESPACE_CLOSE
MappingQ1Eulerian(const VectorType &euler_transform_vectors, const DoFHandler< dim, spacedim > &shiftmap_dof_handler)
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const
virtual MappingQ1Eulerian< dim, VectorType, spacedim > * clone() const
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:249
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
#define Assert(cond, exc)
Definition: exceptions.h:294
SmartPointer< const VectorType, MappingQ1Eulerian< dim, VectorType, spacedim > > euler_transform_vectors
SmartPointer< const DoFHandler< dim, spacedim >, MappingQ1Eulerian< dim, VectorType, spacedim > > shiftmap_dof_handler
virtual std_cxx11::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const