Reference documentation for deal.II version 9.2.0
\(\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 - 2020 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 
18 #include <deal.II/base/utilities.h>
19 
22 
23 #include <deal.II/fe/fe.h>
24 #include <deal.II/fe/fe_tools.h>
27 
29 
33 #include <deal.II/lac/la_vector.h>
38 #include <deal.II/lac/vector.h>
39 
41 
42 
43 
44 // .... MAPPING Q EULERIAN CONSTRUCTOR
45 
46 template <int dim, class VectorType, int spacedim>
49  const unsigned int degree,
50  const MappingQEulerian<dim, VectorType, spacedim> &mapping_q_eulerian)
51  : MappingQGeneric<dim, spacedim>(degree)
52  , mapping_q_eulerian(mapping_q_eulerian)
53  , support_quadrature(degree)
54  , fe_values(mapping_q_eulerian.euler_dof_handler->get_fe(),
55  support_quadrature,
57 {}
58 
59 
60 
61 template <int dim, class VectorType, int spacedim>
63  const unsigned int degree,
65  const VectorType & euler_vector,
66  const unsigned int level)
67  : MappingQ<dim, spacedim>(degree, true)
70  , level(level)
71 {
72  // reset the q1 mapping we use for interior cells (and previously
73  // set by the MappingQ constructor) to a MappingQ1Eulerian with the
74  // current vector
75  this->q1_mapping =
76  std::make_shared<MappingQ1Eulerian<dim, VectorType, spacedim>>(
78 
79  // also reset the qp mapping pointer with our own class
80  this->qp_mapping = std::make_shared<MappingQEulerianGeneric>(degree, *this);
81 }
82 
83 
84 
85 template <int dim, class VectorType, int spacedim>
86 std::unique_ptr<Mapping<dim, spacedim>>
88 {
89  return std_cxx14::make_unique<MappingQEulerian<dim, VectorType, spacedim>>(
90  this->get_degree(), *euler_dof_handler, *euler_vector);
91 }
92 
93 
94 
95 // .... SUPPORT QUADRATURE CONSTRUCTOR
96 
97 template <int dim, class VectorType, int spacedim>
99  SupportQuadrature::SupportQuadrature(const unsigned int map_degree)
100  : Quadrature<dim>(Utilities::fixed_power<dim>(map_degree + 1))
101 {
102  // first we determine the support points on the unit cell in lexicographic
103  // order, which are (in accordance with MappingQ) the support points of
104  // QGaussLobatto.
105  const QGaussLobatto<dim> q_iterated(map_degree + 1);
106  const unsigned int n_q_points = q_iterated.size();
107 
108  // we then need to define a renumbering vector that allows us to go from a
109  // lexicographic numbering scheme to a hierarchic one.
110  const std::vector<unsigned int> renumber =
111  FETools::lexicographic_to_hierarchic_numbering<dim>(map_degree);
112 
113  // finally we assign the quadrature points in the required order.
114  for (unsigned int q = 0; q < n_q_points; ++q)
115  this->quadrature_points[renumber[q]] = q_iterated.point(q);
116 }
117 
118 
119 
120 // .... COMPUTE MAPPING SUPPORT POINTS
121 
122 template <int dim, class VectorType, int spacedim>
123 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
125  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
126 {
127  // get the vertices as the first 2^dim mapping support points
128  const std::vector<Point<spacedim>> a =
129  dynamic_cast<const MappingQEulerianGeneric &>(*this->qp_mapping)
131 
132  std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
133  vertex_locations;
134  std::copy(a.begin(),
136  vertex_locations.begin());
137 
138  return vertex_locations;
139 }
140 
141 
142 
143 template <int dim, class VectorType, int spacedim>
144 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
147  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
148 {
149  return mapping_q_eulerian.get_vertices(cell);
150 }
151 
152 
153 
154 template <int dim, class VectorType, int spacedim>
155 bool
158 {
159  return false;
160 }
161 
162 
163 
164 template <int dim, class VectorType, int spacedim>
165 std::vector<Point<spacedim>>
168  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
169 {
170  const bool mg_vector =
172 
173  const types::global_dof_index n_dofs =
174  mg_vector ?
175  mapping_q_eulerian.euler_dof_handler->n_dofs(mapping_q_eulerian.level) :
176  mapping_q_eulerian.euler_dof_handler->n_dofs();
177  const types::global_dof_index vector_size =
178  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
188  *cell, mapping_q_eulerian.euler_dof_handler);
189 
190  Assert(mg_vector || dof_cell->is_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 =
205  mapping_q_eulerian.euler_dof_handler->get_fe(0).n_components();
206 
207  Assert(n_components >= spacedim,
209 
210  std::vector<Vector<typename VectorType::value_type>> shift_vector(
212 
213  std::vector<types::global_dof_index> dof_indices(
214  mapping_q_eulerian.euler_dof_handler->get_fe(0).dofs_per_cell);
215  // fill shift vector for each support point using an fe_values object. make
216  // sure that the fe_values variable isn't used simultaneously from different
217  // threads
218  std::lock_guard<std::mutex> lock(fe_values_mutex);
219  fe_values.reinit(dof_cell);
220  if (mg_vector)
221  {
222  dof_cell->get_mg_dof_indices(dof_indices);
223  fe_values.get_function_values(*mapping_q_eulerian.euler_vector,
224  dof_indices,
225  shift_vector);
226  }
227  else
228  fe_values.get_function_values(*mapping_q_eulerian.euler_vector,
229  shift_vector);
230 
231  // and finally compute the positions of the support points in the deformed
232  // configuration.
233  std::vector<Point<spacedim>> a(n_support_pts);
234  for (unsigned int q = 0; q < n_support_pts; ++q)
235  {
236  a[q] = fe_values.quadrature_point(q);
237  for (unsigned int d = 0; d < spacedim; ++d)
238  a[q](d) += shift_vector[q](d);
239  }
240 
241  return a;
242 }
243 
244 
245 
246 template <int dim, class VectorType, int spacedim>
249  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
251  const Quadrature<dim> & quadrature,
252  const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
254  &output_data) const
255 {
256  // call the function of the base class, but ignoring
257  // any potentially detected cell similarity between
258  // the current and the previous cell
261  quadrature,
262  internal_data,
263  output_data);
264  // also return the updated flag since any detected
265  // similarity wasn't based on the mapped field, but
266  // the original vertices which are meaningless
268 }
269 
270 
271 
272 // explicit instantiations
273 #include "mapping_q_eulerian.inst"
274 
275 
mapping_q_eulerian.h
la_vector.h
MappingQ
Definition: mapping_manifold.h:33
DoFHandler::cell_iterator
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:357
update_quadrature_points
@ update_quadrature_points
Transformed quadrature points.
Definition: fe_update_flags.h:122
MappingQ::fill_fe_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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
Definition: mapping_q.cc:249
internal::FEValuesImplementation::MappingRelatedData
Definition: fe_update_flags.h:418
MappingQEulerian::ExcInactiveCell
static ::ExceptionBase & ExcInactiveCell()
MappingQEulerian::clone
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_q_eulerian.cc:87
MappingQGeneric
Definition: mapping_q_generic.h:135
MappingQEulerian::MappingQEulerianGeneric::preserves_vertex_locations
virtual bool preserves_vertex_locations() const override
Definition: mapping_q_eulerian.cc:157
Quadrature::quadrature_points
std::vector< Point< dim > > quadrature_points
Definition: quadrature.h:284
utilities.h
tria_iterator.h
MappingQEulerian::euler_dof_handler
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim > > euler_dof_handler
Definition: mapping_q_eulerian.h:179
GeometryInfo
Definition: geometry_info.h:1224
VectorType
CellSimilarity::Similarity
Similarity
Definition: fe_update_flags.h:378
update_values
@ update_values
Shape function values.
Definition: fe_update_flags.h:78
MappingQEulerian::MappingQEulerianGeneric::mapping_q_eulerian
const MappingQEulerian< dim, VectorType, spacedim > & mapping_q_eulerian
Definition: mapping_q_eulerian.h:235
DoFTools::n_components
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
DoFHandler
Definition: dof_handler.h:205
la_parallel_vector.h
quadrature_lib.h
Quadrature::point
const Point< dim > & point(const unsigned int i) const
MappingQEulerian::MappingQEulerian
MappingQEulerian(const unsigned int degree, const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector, const unsigned int level=numbers::invalid_unsigned_int)
Definition: mapping_q_eulerian.cc:62
MappingQ< dim, dim >::qp_mapping
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:368
MappingQEulerian::level
const unsigned int level
Definition: mapping_q_eulerian.h:186
block_vector.h
MappingQEulerian::MappingQEulerianGeneric::fe_values_mutex
Threads::Mutex fe_values_mutex
Definition: mapping_q_eulerian.h:269
la_parallel_block_vector.h
MappingQEulerian::euler_vector
SmartPointer< const VectorType, MappingQEulerian< dim, VectorType, spacedim > > euler_vector
Definition: mapping_q_eulerian.h:172
Utilities::fixed_power
T fixed_power(const T t)
Definition: utilities.h:1072
MappingQEulerian::MappingQEulerianGeneric::SupportQuadrature::SupportQuadrature
SupportQuadrature(const unsigned int map_degree)
Definition: mapping_q_eulerian.cc:99
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
fe.h
petsc_block_vector.h
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
QGaussLobatto
Definition: quadrature_lib.h:76
mapping_q1_eulerian.h
MappingQEulerian::MappingQEulerianGeneric::support_quadrature
const SupportQuadrature support_quadrature
Definition: mapping_q_eulerian.h:254
fe_tools.h
unsigned int
MappingQEulerian::fill_fe_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::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
Definition: mapping_q_eulerian.cc:248
CellSimilarity::invalid_next_cell
@ invalid_next_cell
Definition: fe_update_flags.h:396
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
Mapping::InternalDataBase
Definition: mapping.h:597
vector.h
dof_handler.h
dof_accessor.h
MappingQEulerian::MappingQEulerianGeneric
Definition: mapping_q_eulerian.h:193
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
MappingQEulerian::get_vertices
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
Definition: mapping_q_eulerian.cc:124
StandardExceptions::ExcDimensionMismatch
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
petsc_vector.h
MappingQEulerian
Definition: mapping_q_eulerian.h:95
Quadrature::size
unsigned int size() const
MappingQEulerian::MappingQEulerianGeneric::fe_values
FEValues< dim, spacedim > fe_values
Definition: mapping_q_eulerian.h:264
memory.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
MappingQ< dim, dim >::q1_mapping
std::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
Definition: mapping_q.h:352
Quadrature
Definition: quadrature.h:85
trilinos_vector.h
Vector< typename VectorType::value_type >
TriaIterator
Definition: tria_iterator.h:578
MappingQEulerian::MappingQEulerianGeneric::get_vertices
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
Definition: mapping_q_eulerian.cc:146
MappingQEulerian::MappingQEulerianGeneric::MappingQEulerianGeneric
MappingQEulerianGeneric(const unsigned int degree, const MappingQEulerian< dim, VectorType, spacedim > &mapping_q_eulerian)
Definition: mapping_q_eulerian.cc:48
Utilities
Definition: cuda.h:31
internal::VectorOperations::copy
void copy(const T *begin, const T *end, U *dest)
Definition: vector_operations_internal.h:67
trilinos_parallel_block_vector.h
MappingQEulerian::MappingQEulerianGeneric::compute_mapping_support_points
virtual std::vector< Point< spacedim > > compute_mapping_support_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
Definition: mapping_q_eulerian.cc:167