Reference documentation for deal.II version 9.0.0
mapping_q_eulerian.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/utilities.h>
17 #include <deal.II/base/quadrature_lib.h>
18 #include <deal.II/base/std_cxx14/memory.h>
19 #include <deal.II/lac/vector.h>
20 #include <deal.II/lac/block_vector.h>
21 #include <deal.II/lac/la_vector.h>
22 #include <deal.II/lac/parallel_vector.h>
23 #include <deal.II/lac/parallel_block_vector.h>
24 #include <deal.II/lac/petsc_parallel_vector.h>
25 #include <deal.II/lac/petsc_parallel_block_vector.h>
26 #include <deal.II/lac/trilinos_vector.h>
27 #include <deal.II/lac/trilinos_parallel_block_vector.h>
28 #include <deal.II/lac/trilinos_parallel_block_vector.h>
29 #include <deal.II/grid/tria_iterator.h>
30 #include <deal.II/dofs/dof_handler.h>
31 #include <deal.II/dofs/dof_accessor.h>
32 #include <deal.II/fe/fe.h>
33 #include <deal.II/fe/fe_tools.h>
34 #include <deal.II/fe/mapping_q_eulerian.h>
35 #include <deal.II/fe/mapping_q1_eulerian.h>
36 
37 DEAL_II_NAMESPACE_OPEN
38 
39 
40 
41 // .... MAPPING Q EULERIAN CONSTRUCTOR
42 
43 template <int dim, class VectorType, int spacedim>
45 MappingQEulerianGeneric (const unsigned int degree,
46  const MappingQEulerian<dim, VectorType, spacedim> &mapping_q_eulerian)
47  :
48  MappingQGeneric<dim,spacedim>(degree),
49  mapping_q_eulerian (mapping_q_eulerian),
50  support_quadrature(degree),
51  fe_values(mapping_q_eulerian.euler_dof_handler->get_fe(),
52  support_quadrature,
54 {}
55 
56 
57 
58 template <int dim, class VectorType, int spacedim>
60 MappingQEulerian (const unsigned int degree,
62  const VectorType &euler_vector,
63  const unsigned int level)
64  :
65  MappingQ<dim,spacedim>(degree, true),
68  level(level)
69 {
70  // reset the q1 mapping we use for interior cells (and previously
71  // set by the MappingQ constructor) to a MappingQ1Eulerian with the
72  // current vector
73  this->q1_mapping = std::make_shared<MappingQ1Eulerian<dim,VectorType,spacedim>>
75 
76  // also reset the qp mapping pointer with our own class
77  this->qp_mapping = std::make_shared<MappingQEulerianGeneric>(degree,*this);
78 }
79 
80 
81 
82 template <int dim, class VectorType, int spacedim>
83 std::unique_ptr<Mapping<dim,spacedim> >
85 {
86  return std_cxx14::make_unique<MappingQEulerian<dim,VectorType,spacedim>>
87  (this->get_degree(), *euler_dof_handler, *euler_vector);
88 }
89 
90 
91 
92 // .... SUPPORT QUADRATURE CONSTRUCTOR
93 
94 template <int dim, class VectorType, int spacedim>
97 SupportQuadrature (const unsigned int map_degree)
98  :
99  Quadrature<dim>(Utilities::fixed_power<dim>(map_degree+1))
100 {
101  // first we determine the support points on the unit cell in lexicographic
102  // order, which are (in accordance with MappingQ) the support points of
103  // QGaussLobatto.
104  const QGaussLobatto<dim> q_iterated(map_degree+1);
105  const unsigned int n_q_points = q_iterated.size();
106 
107  // we then need to define a renumbering vector that allows us to go from a
108  // lexicographic numbering scheme to a hierarchic one. this fragment is
109  // taking almost verbatim from the MappingQ class.
110  std::vector<unsigned int> renumber(n_q_points);
111  std::vector<unsigned int> dpo(dim+1, 1U);
112  for (unsigned int i=1; i<dpo.size(); ++i)
113  dpo[i]=dpo[i-1]*(map_degree-1);
114 
116  renumber);
117 
118  // finally we assign the quadrature points in the required order.
119  for (unsigned int q=0; q<n_q_points; ++q)
120  this->quadrature_points[renumber[q]] = q_iterated.point(q);
121 }
122 
123 
124 
125 // .... COMPUTE MAPPING SUPPORT POINTS
126 
127 template <int dim, class VectorType, int spacedim>
128 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
131 (const typename Triangulation<dim,spacedim>::cell_iterator &cell) const
132 {
133  // get the vertices as the first 2^dim mapping support points
134  const std::vector<Point<spacedim> > a
135  = dynamic_cast<const MappingQEulerianGeneric &>(*this->qp_mapping).compute_mapping_support_points(cell);
136 
137  std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell> vertex_locations;
138  std::copy (a.begin(),
140  vertex_locations.begin());
141 
142  return vertex_locations;
143 }
144 
145 
146 
147 template <int dim, class VectorType, int spacedim>
148 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
151 (const typename Triangulation<dim,spacedim>::cell_iterator &cell) const
152 {
153  return mapping_q_eulerian.get_vertices (cell);
154 }
155 
156 
157 
158 template <int dim, class VectorType, int spacedim>
159 bool
162 {
163  return false;
164 }
165 
166 
167 
168 template <int dim, class VectorType, int spacedim>
169 std::vector<Point<spacedim> >
172 {
173  const bool mg_vector = mapping_q_eulerian.level != numbers::invalid_unsigned_int;
174 
175  const types::global_dof_index n_dofs = mg_vector ?
176  mapping_q_eulerian.euler_dof_handler->n_dofs(mapping_q_eulerian.level) :
177  mapping_q_eulerian.euler_dof_handler->n_dofs();
178  const types::global_dof_index vector_size = mapping_q_eulerian.euler_vector->size();
179 
180  (void)n_dofs;
181  (void)vector_size;
182 
183  AssertDimension(vector_size,n_dofs);
184 
185  // we then transform our tria iterator into a dof iterator so we can access
186  // data not associated with triangulations
187  typename DoFHandler<dim,spacedim>::cell_iterator dof_cell(*cell,
188  mapping_q_eulerian.euler_dof_handler);
189 
190  Assert (mg_vector || dof_cell->active() == true, ExcInactiveCell());
191 
192  // our quadrature rule is chosen so that each quadrature point corresponds
193  // to a support point in the undeformed configuration. We can then query
194  // the given displacement field at these points to determine the shift
195  // vector that maps the support points to the deformed configuration.
196 
197  // We assume that the given field contains dim displacement components, but
198  // that there may be other solution components as well (e.g. pressures).
199  // this class therefore assumes that the first dim components represent the
200  // actual shift vector we need, and simply ignores any components after
201  // that. This implies that the user should order components appropriately,
202  // or create a separate dof handler for the displacements.
203  const unsigned int n_support_pts = support_quadrature.size();
204  const unsigned int n_components = mapping_q_eulerian.euler_dof_handler->get_fe(0).n_components();
205 
206  Assert (n_components >= spacedim, ExcDimensionMismatch(n_components, spacedim) );
207 
208  std::vector<Vector<typename VectorType::value_type> >
209  shift_vector(n_support_pts,
211 
212  std::vector<types::global_dof_index> dof_indices(mapping_q_eulerian.euler_dof_handler->get_fe(0).dofs_per_cell);
213  // fill shift vector for each support point using an fe_values object. make
214  // sure that the fe_values variable isn't used simultaneously from different
215  // threads
217  fe_values.reinit(dof_cell);
218  if (mg_vector)
219  {
220  dof_cell->get_mg_dof_indices(dof_indices);
221  fe_values.get_function_values(*mapping_q_eulerian.euler_vector,
222  dof_indices,
223  shift_vector);
224  }
225  else
226  fe_values.get_function_values(*mapping_q_eulerian.euler_vector, shift_vector);
227 
228  // and finally compute the positions of the support points in the deformed
229  // configuration.
230  std::vector<Point<spacedim> > a(n_support_pts);
231  for (unsigned int q=0; q<n_support_pts; ++q)
232  {
233  a[q] = fe_values.quadrature_point(q);
234  for (unsigned int d=0; d<spacedim; ++d)
235  a[q](d) += shift_vector[q](d);
236  }
237 
238  return a;
239 }
240 
241 
242 
243 template <int dim, class VectorType, int spacedim>
248  const Quadrature<dim> &quadrature,
249  const typename Mapping<dim,spacedim>::InternalDataBase &internal_data,
251 {
252  // call the function of the base class, but ignoring
253  // any potentially detected cell similarity between
254  // the current and the previous cell
257  quadrature,
258  internal_data,
259  output_data);
260  // also return the updated flag since any detected
261  // similarity wasn't based on the mapped field, but
262  // the original vertices which are meaningless
264 }
265 
266 
267 
268 // explicit instantiations
269 #include "mapping_q_eulerian.inst"
270 
271 
272 DEAL_II_NAMESPACE_CLOSE
Shape function values.
static const unsigned int invalid_unsigned_int
Definition: types.h:173
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
const unsigned int level
static ::ExceptionBase & ExcInactiveCell()
const Point< dim > & point(const unsigned int i) const
Transformed quadrature points.
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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
std::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
Definition: mapping_q.h:359
ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:257
MappingQEulerianGeneric(const unsigned int degree, const MappingQEulerian< dim, VectorType, spacedim > &mapping_q_eulerian)
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
MappingQEulerian(const unsigned int degree, const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector, const unsigned int level=numbers::invalid_unsigned_int)
SmartPointer< const VectorType, MappingQEulerian< dim, VectorType, spacedim > > euler_vector
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
std::vector< Point< dim > > quadrature_points
Definition: quadrature.h:252
unsigned int size() const
void lexicographic_to_hierarchic_numbering(const FiniteElementData< dim > &fe_data, std::vector< unsigned int > &l2h)
Definition: cuda.h:31
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:375
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim > > euler_dof_handler
virtual bool preserves_vertex_locations() const override
const MappingQEulerian< dim, VectorType, spacedim > & mapping_q_eulerian