Reference documentation for deal.II version 8.4.2
mapping_q_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/utilities.h>
17 #include <deal.II/base/quadrature_lib.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 #include <deal.II/fe/fe_tools.h>
28 #include <deal.II/fe/mapping_q_eulerian.h>
29 #include <deal.II/fe/mapping_q1_eulerian.h>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 
34 
35 // .... MAPPING Q EULERIAN CONSTRUCTOR
36 
37 template <int dim, class EulerVectorType, int spacedim>
39 MappingQEulerianGeneric (const unsigned int degree,
40  const MappingQEulerian<dim, EulerVectorType, spacedim> &mapping_q_eulerian)
41  :
42  MappingQGeneric<dim,spacedim>(degree),
43  mapping_q_eulerian (mapping_q_eulerian),
44  support_quadrature(degree),
45  fe_values(mapping_q_eulerian.euler_dof_handler->get_fe(),
46  support_quadrature,
48 {}
49 
50 template <int dim, class EulerVectorType, int spacedim>
52 MappingQEulerian (const unsigned int degree,
53  const EulerVectorType &euler_vector,
55  :
56  MappingQ<dim,spacedim>(degree, true),
58  euler_dof_handler(&euler_dof_handler)
59 {
60  // reset the q1 mapping we use for interior cells (and previously
61  // set by the MappingQ constructor) to a MappingQ1Eulerian with the
62  // current vector
64  euler_dof_handler));
65 
66  // also reset the qp mapping pointer with our own class
67  this->qp_mapping.reset (new MappingQEulerianGeneric(degree,*this));
68 }
69 
70 
71 
72 template <int dim, class EulerVectorType, int spacedim>
74 MappingQEulerian (const unsigned int degree,
75  const DoFHandler<dim,spacedim> &euler_dof_handler,
76  const EulerVectorType &euler_vector)
77  :
78  MappingQ<dim,spacedim>(degree, true),
80  euler_dof_handler(&euler_dof_handler)
81 {
82  // reset the q1 mapping we use for interior cells (and previously
83  // set by the MappingQ constructor) to a MappingQ1Eulerian with the
84  // current vector
86  euler_dof_handler));
87 
88  // also reset the qp mapping pointer with our own class
89  this->qp_mapping.reset (new MappingQEulerianGeneric(degree,*this));
90 }
91 
92 
93 
94 template <int dim, class EulerVectorType, int spacedim>
97 {
99  *euler_vector,
101 }
102 
103 
104 
105 // .... SUPPORT QUADRATURE CONSTRUCTOR
106 
107 template <int dim, class EulerVectorType, int spacedim>
110 SupportQuadrature (const unsigned int map_degree)
111  :
112  Quadrature<dim>(Utilities::fixed_power<dim>(map_degree+1))
113 {
114  // first we determine the support points on the unit cell in lexicographic
115  // order, which are (in accordance with MappingQ) the support points of
116  // QGaussLobatto.
117  const QGaussLobatto<dim> q_iterated(map_degree+1);
118  const unsigned int n_q_points = q_iterated.size();
119 
120  // we then need to define a renumbering vector that allows us to go from a
121  // lexicographic numbering scheme to a hierarchic one. this fragment is
122  // taking almost verbatim from the MappingQ class.
123  std::vector<unsigned int> renumber(n_q_points);
124  std::vector<unsigned int> dpo(dim+1, 1U);
125  for (unsigned int i=1; i<dpo.size(); ++i)
126  dpo[i]=dpo[i-1]*(map_degree-1);
127 
129  renumber);
130 
131  // finally we assign the quadrature points in the required order.
132  for (unsigned int q=0; q<n_q_points; ++q)
133  this->quadrature_points[renumber[q]] = q_iterated.point(q);
134 }
135 
136 
137 
138 // .... COMPUTE MAPPING SUPPORT POINTS
139 
140 template <int dim, class EulerVectorType, int spacedim>
141 std_cxx11::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
144 (const typename Triangulation<dim,spacedim>::cell_iterator &cell) const
145 {
146  // get the vertices as the first 2^dim mapping support points
147  const std::vector<Point<spacedim> > a
148  = dynamic_cast<const MappingQEulerianGeneric &>(*this->qp_mapping).compute_mapping_support_points(cell);
149 
150  std_cxx11::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell> vertex_locations;
151  std::copy (a.begin(),
153  vertex_locations.begin());
154 
155  return vertex_locations;
156 }
157 
158 
159 
160 template <int dim, class EulerVectorType, int spacedim>
161 std_cxx11::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
164 (const typename Triangulation<dim,spacedim>::cell_iterator &cell) const
165 {
166  return mapping_q_eulerian.get_vertices (cell);
167 }
168 
169 
170 
171 
172 template <int dim, class EulerVectorType, int spacedim>
173 std::vector<Point<spacedim> >
176 {
177  // first, basic assertion with respect to vector size,
178 
179  const types::global_dof_index n_dofs = mapping_q_eulerian.euler_dof_handler->n_dofs();
180  const types::global_dof_index vector_size = mapping_q_eulerian.euler_vector->size();
181  (void)n_dofs;
182  (void)vector_size;
183 
184  AssertDimension(vector_size,n_dofs);
185 
186  // we then transform our tria iterator into a dof iterator so we can access
187  // data not associated with triangulations
188  typename DoFHandler<dim,spacedim>::cell_iterator dof_cell(*cell,
189  mapping_q_eulerian.euler_dof_handler);
190 
191  Assert (dof_cell->active() == true, ExcInactiveCell());
192 
193  // our quadrature rule is chosen so that each quadrature point corresponds
194  // to a support point in the undeformed configuration. we can then query
195  // the given displacement field at these points to determine the shift
196  // vector that maps the support points to the deformed configuration.
197 
198  // we assume that the given field contains dim displacement components, but
199  // that there may be other solution components as well (e.g. pressures).
200  // this class therefore assumes that the first dim components represent the
201  // actual shift vector we need, and simply ignores any components after
202  // that. this implies that the user should order components appropriately,
203  // or create a separate dof handler for the displacements.
204 
205  const unsigned int n_support_pts = support_quadrature.size();
206  const unsigned int n_components = mapping_q_eulerian.euler_dof_handler->get_fe().n_components();
207 
208  Assert (n_components >= spacedim, ExcDimensionMismatch(n_components, spacedim) );
209 
210  std::vector<Vector<typename EulerVectorType::value_type> >
211  shift_vector(n_support_pts,
213 
214  // fill shift vector for each support point using an fe_values object. make
215  // sure that the fe_values variable isn't used simultaneously from different
216  // threads
218  fe_values.reinit(dof_cell);
219  fe_values.get_function_values(*mapping_q_eulerian.euler_vector, shift_vector);
220 
221  // and finally compute the positions of the support points in the deformed
222  // configuration.
223  std::vector<Point<spacedim> > a(n_support_pts);
224  for (unsigned int q=0; q<n_support_pts; ++q)
225  {
226  a[q] = fe_values.quadrature_point(q);
227  for (unsigned int d=0; d<spacedim; ++d)
228  a[q](d) += shift_vector[q](d);
229  }
230 
231  return a;
232 }
233 
234 
235 
236 template<int dim, class EulerVectorType, int spacedim>
237 CellSimilarity::Similarity
240  const CellSimilarity::Similarity ,
241  const Quadrature<dim> &quadrature,
242  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
244 {
245  // call the function of the base class, but ignoring
246  // any potentially detected cell similarity between
247  // the current and the previous cell
249  CellSimilarity::invalid_next_cell,
250  quadrature,
251  internal_data,
252  output_data);
253  // also return the updated flag since any detected
254  // similarity wasn't based on the mapped field, but
255  // the original vertices which are meaningless
256  return CellSimilarity::invalid_next_cell;
257 }
258 
259 
260 
261 // explicit instantiations
262 #include "mapping_q_eulerian.inst"
263 
264 
265 DEAL_II_NAMESPACE_CLOSE
Shape function values.
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
std_cxx11::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
Definition: mapping_q.h:352
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
virtual std_cxx11::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
const Point< dim > & point(const unsigned int i) const
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:249
virtual Mapping< dim, spacedim > * clone() const
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:294
MappingQEulerianGeneric(const unsigned int degree, const MappingQEulerian< dim, VectorType, spacedim > &mapping_q_eulerian)
Abstract base class for mapping classes.
Definition: dof_tools.h:52
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
SmartPointer< const VectorType, MappingQEulerian< dim, VectorType, spacedim > > euler_vector
MappingQEulerian(const unsigned int degree, const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector)
std::vector< Point< dim > > quadrature_points
Definition: quadrature.h:217
unsigned int size() const
void lexicographic_to_hierarchic_numbering(const FiniteElementData< dim > &fe_data, std::vector< unsigned int > &l2h)
Definition: fe_tools.cc:2086
unsigned int get_degree() const
Definition: mpi.h:55
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim > > euler_dof_handler
std_cxx11::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:368
virtual std_cxx11::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
const MappingQEulerian< dim, VectorType, spacedim > & mapping_q_eulerian