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_q1_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 
17 
20 
21 #include <deal.II/fe/fe.h>
23 
25 
29 #include <deal.II/lac/la_vector.h>
34 #include <deal.II/lac/vector.h>
35 
36 #include <array>
37 
39 
40 
41 template <int dim, class VectorType, int spacedim>
43  const DoFHandler<dim, spacedim> &shiftmap_dof_handler,
44  const VectorType & euler_transform_vectors)
45  : MappingQGeneric<dim, spacedim>(1)
46  , euler_transform_vectors(&euler_transform_vectors)
47  , shiftmap_dof_handler(&shiftmap_dof_handler)
48 {}
49 
50 
51 
52 template <int dim, class VectorType, int spacedim>
53 std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
55  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
56 {
57  std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell> vertices;
58  // The assertions can not be in the constructor, since this would
59  // require to call dof_handler.distribute_dofs(fe) *before* the mapping
60  // object is constructed, which is not necessarily what we want.
61 
62  // TODO: Only one of these two assertions should be relevant
63  AssertDimension(spacedim, shiftmap_dof_handler->get_fe().n_dofs_per_vertex());
64  AssertDimension(shiftmap_dof_handler->get_fe(0).n_components(), spacedim);
65 
66  AssertDimension(shiftmap_dof_handler->n_dofs(),
67  euler_transform_vectors->size());
68 
69  // cast the Triangulation<dim>::cell_iterator into a
70  // DoFHandler<dim>::cell_iterator which is necessary for access to
71  // DoFCellAccessor::get_dof_values()
73  *cell, shiftmap_dof_handler);
74 
75  // We require the cell to be active since we can only then get nodal
76  // values for the shifts
77  Assert(dof_cell->is_active() == true, ExcInactiveCell());
78 
79  // now get the values of the shift vectors at the vertices
81  shiftmap_dof_handler->get_fe().dofs_per_cell);
82  dof_cell->get_dof_values(*euler_transform_vectors, mapping_values);
83 
84  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
85  {
86  Point<spacedim> shift_vector;
87 
88  // pick out the value of the shift vector at the present
89  // vertex. since vertex dofs are always numbered first, we can
90  // access them easily
91  for (unsigned int j = 0; j < spacedim; ++j)
92  shift_vector[j] = mapping_values(i * spacedim + j);
93 
94  // compute new support point by old (reference) value and added
95  // shift
96  vertices[i] = cell->vertex(i) + shift_vector;
97  }
98  return vertices;
99 }
100 
101 
102 
103 template <int dim, class VectorType, int spacedim>
104 std::vector<Point<spacedim>>
106  const typename Triangulation<dim, spacedim>::cell_iterator &cell) const
107 {
108  const std::array<Point<spacedim>, GeometryInfo<dim>::vertices_per_cell>
109  vertices = this->get_vertices(cell);
110 
111  std::vector<Point<spacedim>> a(GeometryInfo<dim>::vertices_per_cell);
112  for (const unsigned int i : GeometryInfo<dim>::vertex_indices())
113  a[i] = vertices[i];
114 
115  return a;
116 }
117 
118 
119 
120 template <int dim, class VectorType, int spacedim>
121 std::unique_ptr<Mapping<dim, spacedim>>
123 {
124  return std_cxx14::make_unique<MappingQ1Eulerian<dim, VectorType, spacedim>>(
125  *this);
126 }
127 
128 
129 
130 template <int dim, class VectorType, int spacedim>
133  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
135  const Quadrature<dim> & quadrature,
136  const typename Mapping<dim, spacedim>::InternalDataBase &internal_data,
138  &output_data) const
139 {
140  // call the function of the base class, but ignoring
141  // any potentially detected cell similarity between
142  // the current and the previous cell
144  cell,
146  quadrature,
147  internal_data,
148  output_data);
149  // also return the updated flag since any detected
150  // similarity wasn't based on the mapped field, but
151  // the original vertices which are meaningless
153 }
154 
155 
156 
157 // explicit instantiations
158 #include "mapping_q1_eulerian.inst"
159 
160 
la_vector.h
DoFHandler::cell_iterator
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:357
internal::FEValuesImplementation::MappingRelatedData
Definition: fe_update_flags.h:418
MappingQGeneric
Definition: mapping_q_generic.h:135
tria_iterator.h
MappingQ1Eulerian::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_q1_eulerian.cc:105
GeometryInfo
Definition: geometry_info.h:1224
VectorType
CellSimilarity::Similarity
Similarity
Definition: fe_update_flags.h:378
DoFHandler
Definition: dof_handler.h:205
la_parallel_vector.h
block_vector.h
la_parallel_block_vector.h
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
MappingQ1Eulerian::clone
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_q1_eulerian.cc:122
mapping_q1_eulerian.h
MappingQ1Eulerian::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_q1_eulerian.cc:54
vertices
Point< 3 > vertices[4]
Definition: data_out_base.cc:174
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
MappingQ1Eulerian::MappingQ1Eulerian
MappingQ1Eulerian(const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector)
Definition: mapping_q1_eulerian.cc:42
Point< spacedim >
petsc_vector.h
MappingQ1Eulerian::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_q1_eulerian.cc:132
memory.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
MappingQGeneric::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_generic.cc:2725
Quadrature
Definition: quadrature.h:85
trilinos_vector.h
Vector< typename VectorType::value_type >
TriaIterator
Definition: tria_iterator.h:578
trilinos_parallel_block_vector.h