Reference documentation for deal.II version Git 24041777c5 2021-12-02 20:43:25 -0700
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mapping_q_eulerian.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2021 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
17 #include <deal.II/base/utilities.h>
18 
21 
22 #include <deal.II/fe/fe.h>
23 #include <deal.II/fe/fe_tools.h>
26 
28 
32 #include <deal.II/lac/la_vector.h>
37 #include <deal.II/lac/vector.h>
38 
39 #include <memory>
40 
42 
43 
44 
45 // .... MAPPING Q EULERIAN CONSTRUCTOR
46 
47 
48 template <int dim, class VectorType, int spacedim>
50  const unsigned int degree,
51  const DoFHandler<dim, spacedim> &euler_dof_handler,
52  const VectorType & euler_vector,
53  const unsigned int level)
54  : MappingQ<dim, spacedim>(degree)
55  , euler_vector(&euler_vector)
56  , euler_dof_handler(&euler_dof_handler)
57  , level(level)
58  , support_quadrature(degree)
59  , fe_values(euler_dof_handler.get_fe(),
60  support_quadrature,
62 {}
63 
64 
65 
66 template <int dim, class VectorType, int spacedim>
67 std::unique_ptr<Mapping<dim, spacedim>>
69 {
70  return std::make_unique<MappingQEulerian<dim, VectorType, spacedim>>(
71  this->get_degree(), *euler_dof_handler, *euler_vector, this->level);
72 }
73 
74 
75 
76 // .... SUPPORT QUADRATURE CONSTRUCTOR
77 
78 template <int dim, class VectorType, int spacedim>
80  SupportQuadrature(const unsigned int map_degree)
81  : Quadrature<dim>(Utilities::fixed_power<dim>(map_degree + 1))
82 {
83  // first we determine the support points on the unit cell in lexicographic
84  // order, which are (in accordance with MappingQ) the support points of
85  // QGaussLobatto.
86  const QGaussLobatto<dim> q_iterated(map_degree + 1);
87  const unsigned int n_q_points = q_iterated.size();
88 
89  // we then need to define a renumbering vector that allows us to go from a
90  // lexicographic numbering scheme to a hierarchic one.
91  const std::vector<unsigned int> renumber =
92  FETools::lexicographic_to_hierarchic_numbering<dim>(map_degree);
93 
94  // finally we assign the quadrature points in the required order.
95  for (unsigned int q = 0; q < n_q_points; ++q)
96  this->quadrature_points[renumber[q]] = q_iterated.point(q);
97 }
98 
99 
100 
101 // .... COMPUTE MAPPING SUPPORT POINTS
102 
103 template <int dim, class VectorType, int spacedim>
104 boost::container::small_vector<Point<spacedim>,
107  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
108 {
109  // get the vertices as the first 2^dim mapping support points
110  const std::vector<Point<spacedim>> a = compute_mapping_support_points(cell);
111 
112  boost::container::small_vector<Point<spacedim>,
114  vertex_locations(a.begin(),
116 
117  return vertex_locations;
118 }
119 
120 
121 
122 template <int dim, class VectorType, int spacedim>
123 std::vector<Point<spacedim>>
125  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
126 {
127  const bool mg_vector = level != numbers::invalid_unsigned_int;
128 
129  const types::global_dof_index n_dofs =
130  mg_vector ? euler_dof_handler->n_dofs(level) : euler_dof_handler->n_dofs();
131  const types::global_dof_index vector_size = euler_vector->size();
132 
133  (void)n_dofs;
134  (void)vector_size;
135 
136  AssertDimension(vector_size, n_dofs);
137 
138  // we then transform our tria iterator into a dof iterator so we can access
139  // data not associated with triangulations
140  typename DoFHandler<dim, spacedim>::cell_iterator dof_cell(*cell,
142 
143  Assert(mg_vector || dof_cell->is_active() == true, ExcInactiveCell());
144 
145  // our quadrature rule is chosen so that each quadrature point corresponds
146  // to a support point in the undeformed configuration. We can then query
147  // the given displacement field at these points to determine the shift
148  // vector that maps the support points to the deformed configuration.
149 
150  // We assume that the given field contains dim displacement components, but
151  // that there may be other solution components as well (e.g. pressures).
152  // this class therefore assumes that the first dim components represent the
153  // actual shift vector we need, and simply ignores any components after
154  // that. This implies that the user should order components appropriately,
155  // or create a separate dof handler for the displacements.
156  const unsigned int n_support_pts = support_quadrature.size();
157  const unsigned int n_components = euler_dof_handler->get_fe(0).n_components();
158 
159  Assert(n_components >= spacedim,
160  ExcDimensionMismatch(n_components, spacedim));
161 
162  std::vector<Vector<typename VectorType::value_type>> shift_vector(
163  n_support_pts, Vector<typename VectorType::value_type>(n_components));
164 
165  std::vector<types::global_dof_index> dof_indices(
166  euler_dof_handler->get_fe(0).n_dofs_per_cell());
167  // fill shift vector for each support point using an fe_values object. make
168  // sure that the fe_values variable isn't used simultaneously from different
169  // threads
170  std::lock_guard<std::mutex> lock(fe_values_mutex);
171  fe_values.reinit(dof_cell);
172  if (mg_vector)
173  {
174  dof_cell->get_mg_dof_indices(dof_indices);
175  fe_values.get_function_values(*euler_vector, dof_indices, shift_vector);
176  }
177  else
178  fe_values.get_function_values(*euler_vector, shift_vector);
179 
180  // and finally compute the positions of the support points in the deformed
181  // configuration.
182  std::vector<Point<spacedim>> a(n_support_pts);
183  for (unsigned int q = 0; q < n_support_pts; ++q)
184  {
185  a[q] = fe_values.quadrature_point(q);
186  for (unsigned int d = 0; d < spacedim; ++d)
187  a[q](d) += shift_vector[q](d);
188  }
189 
190  return a;
191 }
192 
193 
194 
195 template <int dim, class VectorType, int spacedim>
198  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
200  const Quadrature<dim> & quadrature,
201  const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
203  &output_data) const
204 {
205  // call the function of the base class, but ignoring
206  // any potentially detected cell similarity between
207  // the current and the previous cell
210  quadrature,
211  internal_data,
212  output_data);
213  // also return the updated flag since any detected
214  // similarity wasn't based on the mapped field, but
215  // the original vertices which are meaningless
217 }
218 
219 
220 // explicit instantiations
221 #include "mapping_q_eulerian.inst"
222 
223 
unsigned int get_degree() const
Definition: mapping_q.cc:444
Shape function values.
static const unsigned int invalid_unsigned_int
Definition: types.h:196
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1655
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim > > euler_dof_handler
const unsigned int level
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
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
T fixed_power(const T t)
Definition: utilities.h:1082
FEValues< dim, spacedim > fe_values
#define Assert(cond, exc)
Definition: exceptions.h:1461
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:404
unsigned int level
Definition: grid_out.cc:4614
const SupportQuadrature support_quadrature
MappingQEulerian(const unsigned int degree, const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector, const unsigned int level=numbers::invalid_unsigned_int)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< Point< dim > > quadrature_points
Definition: quadrature.h:289
unsigned int size() const
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:403
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:466
SupportQuadrature(const unsigned int map_degree)
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
Definition: mapping_q.cc:945
SmartPointer< const VectorType, MappingQEulerian< dim, VectorType, spacedim > > euler_vector
Threads::Mutex fe_values_mutex