Reference documentation for deal.II version 9.1.0-pre
mapping_q.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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/array_view.h>
17 #include <deal.II/base/memory_consumption.h>
18 #include <deal.II/base/polynomial.h>
19 #include <deal.II/base/quadrature.h>
20 #include <deal.II/base/quadrature_lib.h>
21 #include <deal.II/base/std_cxx14/memory.h>
22 #include <deal.II/base/tensor_product_polynomials.h>
23 #include <deal.II/base/utilities.h>
24 
25 #include <deal.II/dofs/dof_accessor.h>
26 
27 #include <deal.II/fe/fe_q.h>
28 #include <deal.II/fe/fe_values.h>
29 #include <deal.II/fe/mapping_q.h>
30 
31 #include <deal.II/grid/tria_iterator.h>
32 
33 #include <deal.II/lac/full_matrix.h>
34 
35 #include <memory>
36 #include <numeric>
37 
38 DEAL_II_NAMESPACE_OPEN
39 
40 
41 template <int dim, int spacedim>
43  : use_mapping_q1_on_current_cell(false)
44 {}
45 
46 
47 
48 template <int dim, int spacedim>
49 std::size_t
51 {
52  return (
57 }
58 
59 
60 
61 template <int dim, int spacedim>
62 MappingQ<dim, spacedim>::MappingQ(const unsigned int degree,
63  const bool use_mapping_q_on_all_cells)
64  : polynomial_degree(degree)
65  ,
66 
67  // see whether we want to use *this* mapping objects on *all* cells,
68  // or defer to an explicit Q1 mapping on interior cells. if
69  // degree==1, then we are already that Q1 mapping, so we don't need
70  // it; if dim!=spacedim, there is also no need for anything because
71  // we're most likely on a curved manifold
72  use_mapping_q_on_all_cells(degree == 1 || use_mapping_q_on_all_cells ||
73  (dim != spacedim))
74  ,
75  // create a Q1 mapping for use on interior cells (if necessary)
76  // or to create a good initial guess in transform_real_to_unit_cell()
77  q1_mapping(std::make_shared<MappingQGeneric<dim, spacedim>>(1))
78  ,
79 
80  // create a Q_p mapping; if p=1, simply share the Q_1 mapping already
81  // created via the shared_ptr objects
82  qp_mapping(this->polynomial_degree > 1 ?
83  std::make_shared<const MappingQGeneric<dim, spacedim>>(degree) :
84  q1_mapping)
85 {}
86 
87 
88 
89 template <int dim, int spacedim>
93 {
94  // Note that we really do have to use clone() here, since mapping.q1_mapping
95  // may be MappingQ1Eulerian and mapping.qp_mapping may be MappingQEulerian.
96  std::shared_ptr<const Mapping<dim, spacedim>> other_q1_map =
97  mapping.q1_mapping->clone();
98  q1_mapping = std::dynamic_pointer_cast<const MappingQGeneric<dim, spacedim>>(
99  other_q1_map);
100  Assert(q1_mapping != nullptr, ExcInternalError());
101  Assert(q1_mapping->get_degree() == 1, ExcInternalError());
102 
103  // Same as the other constructor: if possible reuse the Q1 mapping
104  if (this->polynomial_degree == 1)
105  {
107  }
108  else
109  {
110  std::shared_ptr<const Mapping<dim, spacedim>> other_qp_map =
111  mapping.qp_mapping->clone();
112  qp_mapping =
113  std::dynamic_pointer_cast<const MappingQGeneric<dim, spacedim>>(
114  other_qp_map);
115  Assert(qp_mapping != nullptr, ExcInternalError());
116  }
117 }
118 
119 
120 
121 template <int dim, int spacedim>
122 unsigned int
124 {
125  return polynomial_degree;
126 }
127 
128 
129 
130 template <int dim, int spacedim>
131 inline bool
133 {
134  return true;
135 }
136 
137 
138 
139 template <int dim, int spacedim>
142 {
143  return (q1_mapping->requires_update_flags(in) |
144  qp_mapping->requires_update_flags(in));
145 }
146 
147 
148 
149 template <int dim, int spacedim>
150 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
152  const Quadrature<dim> &quadrature) const
153 {
154  auto data = std_cxx14::make_unique<InternalData>();
155 
156  // build the Q1 and Qp internal data objects in parallel
158  std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
160  *qp_mapping,
161  update_flags,
162  quadrature);
163 
165  data->mapping_q1_data = Utilities::dynamic_unique_cast<
167  std::move(q1_mapping->get_data(update_flags, quadrature)));
168 
169  // wait for the task above to finish and use returned value
170  data->mapping_qp_data = Utilities::dynamic_unique_cast<
171  typename MappingQGeneric<dim, spacedim>::InternalData>(
172  std::move(do_get_data.return_value()));
173  return std::move(data);
174 }
175 
176 
177 
178 template <int dim, int spacedim>
179 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
181  const UpdateFlags update_flags,
182  const Quadrature<dim - 1> &quadrature) const
183 {
184  auto data = std_cxx14::make_unique<InternalData>();
185 
186  // build the Q1 and Qp internal data objects in parallel
188  std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
189  do_get_data =
191  *qp_mapping,
192  update_flags,
193  quadrature);
194 
196  data->mapping_q1_data = Utilities::dynamic_unique_cast<
198  std::move(q1_mapping->get_face_data(update_flags, quadrature)));
199 
200  // wait for the task above to finish and use returned value
201  data->mapping_qp_data = Utilities::dynamic_unique_cast<
202  typename MappingQGeneric<dim, spacedim>::InternalData>(
203  std::move(do_get_data.return_value()));
204  return std::move(data);
205 }
206 
207 
208 
209 template <int dim, int spacedim>
210 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
212  const UpdateFlags update_flags,
213  const Quadrature<dim - 1> &quadrature) const
214 {
215  auto data = std_cxx14::make_unique<InternalData>();
216 
217  // build the Q1 and Qp internal data objects in parallel
219  std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>>
220  do_get_data =
222  *qp_mapping,
223  update_flags,
224  quadrature);
225 
227  data->mapping_q1_data = Utilities::dynamic_unique_cast<
229  std::move(q1_mapping->get_subface_data(update_flags, quadrature)));
230 
231  // wait for the task above to finish and use returned value
232  data->mapping_qp_data = Utilities::dynamic_unique_cast<
233  typename MappingQGeneric<dim, spacedim>::InternalData>(
234  std::move(do_get_data.return_value()));
235  return std::move(data);
236 }
237 
238 
239 // Note that the CellSimilarity flag is modifiable, since MappingQ can need to
240 // recalculate data even when cells are similar.
241 template <int dim, int spacedim>
244  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
245  const CellSimilarity::Similarity cell_similarity,
246  const Quadrature<dim> & quadrature,
247  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
249  &output_data) const
250 {
251  // convert data object to internal data for this class. fails with an
252  // exception if that is not possible
253  Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
254  ExcInternalError());
255  const InternalData &data = static_cast<const InternalData &>(internal_data);
256 
257  // check whether this cell needs the full mapping or can be treated by a
258  // reduced Q1 mapping, e.g. if the cell is in the interior of the domain
260  !(use_mapping_q_on_all_cells || cell->has_boundary_lines());
261 
262 
263  // call the base class. we need to ensure that the flag indicating whether
264  // we can use some similarity has to be modified - for a general MappingQ,
265  // the data needs to be recomputed anyway since then the mapping changes the
266  // data. this needs to be known also for later operations, so modify the
267  // variable here. this also affects the calculation of the next cell -- if
268  // we use Q1 data on the next cell, the data will still be invalid.
269  const CellSimilarity::Similarity updated_cell_similarity =
270  ((data.use_mapping_q1_on_current_cell == false) &&
271  (this->polynomial_degree > 1) ?
273  cell_similarity);
274 
275  // depending on the results above, decide whether the Q1 mapping or
276  // the Qp mapping needs to handle this cell
278  q1_mapping->fill_fe_values(cell,
279  updated_cell_similarity,
280  quadrature,
281  *data.mapping_q1_data,
282  output_data);
283  else
284  qp_mapping->fill_fe_values(cell,
285  updated_cell_similarity,
286  quadrature,
287  *data.mapping_qp_data,
288  output_data);
289 
290  return updated_cell_similarity;
291 }
292 
293 
294 
295 template <int dim, int spacedim>
296 void
298  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
299  const unsigned int face_no,
300  const Quadrature<dim - 1> & quadrature,
301  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
303  &output_data) const
304 {
305  // convert data object to internal data for this class. fails with an
306  // exception if that is not possible
307  Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
308  ExcInternalError());
309  const InternalData &data = static_cast<const InternalData &>(internal_data);
310 
311  // check whether this cell needs the full mapping or can be treated by a
312  // reduced Q1 mapping, e.g. if the cell is entirely in the interior of the
313  // domain. note that it is not sufficient to ask whether the present _face_
314  // is in the interior, as the mapping on the face depends on the mapping of
315  // the cell, which in turn depends on the fact whether _any_ of the faces of
316  // this cell is at the boundary, not only the present face
318  !(use_mapping_q_on_all_cells || cell->has_boundary_lines());
319 
320  // depending on the results above, decide whether the Q1 mapping or
321  // the Qp mapping needs to handle this cell
323  q1_mapping->fill_fe_face_values(
324  cell, face_no, quadrature, *data.mapping_q1_data, output_data);
325  else
326  qp_mapping->fill_fe_face_values(
327  cell, face_no, quadrature, *data.mapping_qp_data, output_data);
328 }
329 
330 
331 template <int dim, int spacedim>
332 void
334  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
335  const unsigned int face_no,
336  const unsigned int subface_no,
337  const Quadrature<dim - 1> & quadrature,
338  const typename Mapping<dim, spacedim>::InternalDataBase & internal_data,
340  &output_data) const
341 {
342  // convert data object to internal data for this class. fails with an
343  // exception if that is not possible
344  Assert(dynamic_cast<const InternalData *>(&internal_data) != nullptr,
345  ExcInternalError());
346  const InternalData &data = static_cast<const InternalData &>(internal_data);
347 
348  // check whether this cell needs the full mapping or can be treated by a
349  // reduced Q1 mapping, e.g. if the cell is entirely in the interior of the
350  // domain. note that it is not sufficient to ask whether the present _face_
351  // is in the interior, as the mapping on the face depends on the mapping of
352  // the cell, which in turn depends on the fact whether _any_ of the faces of
353  // this cell is at the boundary, not only the present face
355  !(use_mapping_q_on_all_cells || cell->has_boundary_lines());
356 
357  // depending on the results above, decide whether the Q1 mapping or
358  // the Qp mapping needs to handle this cell
360  q1_mapping->fill_fe_subface_values(cell,
361  face_no,
362  subface_no,
363  quadrature,
364  *data.mapping_q1_data,
365  output_data);
366  else
367  qp_mapping->fill_fe_subface_values(cell,
368  face_no,
369  subface_no,
370  quadrature,
371  *data.mapping_qp_data,
372  output_data);
373 }
374 
375 
376 
377 template <int dim, int spacedim>
378 void
380  const ArrayView<const Tensor<1, dim>> & input,
381  const MappingType mapping_type,
382  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
383  const ArrayView<Tensor<1, spacedim>> & output) const
384 {
385  AssertDimension(input.size(), output.size());
386 
387  const InternalData *data = dynamic_cast<const InternalData *>(&mapping_data);
388  Assert(data != nullptr, ExcInternalError());
389 
390  // check whether we should in fact work on the Q1 portion of it
391  if (data->use_mapping_q1_on_current_cell)
392  q1_mapping->transform(input, mapping_type, *data->mapping_q1_data, output);
393  else
394  // otherwise use the full mapping
395  qp_mapping->transform(input, mapping_type, *data->mapping_qp_data, output);
396 }
397 
398 
399 
400 template <int dim, int spacedim>
401 void
403  const ArrayView<const DerivativeForm<1, dim, spacedim>> &input,
404  const MappingType mapping_type,
405  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
406  const ArrayView<Tensor<2, spacedim>> & output) const
407 {
408  AssertDimension(input.size(), output.size());
409  Assert((dynamic_cast<const typename MappingQ<dim, spacedim>::InternalData *>(
410  &mapping_data) != nullptr),
411  ExcInternalError());
412  const InternalData *data = dynamic_cast<const InternalData *>(&mapping_data);
413 
414  // check whether we should in fact work on the Q1 portion of it
416  q1_mapping->transform(input, mapping_type, *data->mapping_q1_data, output);
417  else
418  // otherwise use the full mapping
419  qp_mapping->transform(input, mapping_type, *data->mapping_qp_data, output);
420 }
421 
422 
423 
424 template <int dim, int spacedim>
425 void
427  const ArrayView<const Tensor<2, dim>> & input,
428  const MappingType mapping_type,
429  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
430  const ArrayView<Tensor<2, spacedim>> & output) const
431 {
432  AssertDimension(input.size(), output.size());
433  Assert((dynamic_cast<const typename MappingQ<dim, spacedim>::InternalData *>(
434  &mapping_data) != nullptr),
435  ExcInternalError());
436  const InternalData *data = dynamic_cast<const InternalData *>(&mapping_data);
437 
438  // check whether we should in fact work on the Q1 portion of it
440  q1_mapping->transform(input, mapping_type, *data->mapping_q1_data, output);
441  else
442  // otherwise use the full mapping
443  qp_mapping->transform(input, mapping_type, *data->mapping_qp_data, output);
444 }
445 
446 
447 
448 template <int dim, int spacedim>
449 void
451  const ArrayView<const DerivativeForm<2, dim, spacedim>> &input,
452  const MappingType mapping_type,
453  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
454  const ArrayView<Tensor<3, spacedim>> & output) const
455 {
456  AssertDimension(input.size(), output.size());
457  Assert((dynamic_cast<const typename MappingQ<dim, spacedim>::InternalData *>(
458  &mapping_data) != nullptr),
459  ExcInternalError());
460  const InternalData *data = dynamic_cast<const InternalData *>(&mapping_data);
461 
462  // check whether we should in fact work on the Q1 portion of it
464  q1_mapping->transform(input, mapping_type, *data->mapping_q1_data, output);
465  else
466  // otherwise use the full mapping
467  qp_mapping->transform(input, mapping_type, *data->mapping_qp_data, output);
468 }
469 
470 
471 
472 template <int dim, int spacedim>
473 void
475  const ArrayView<const Tensor<3, dim>> & input,
476  const MappingType mapping_type,
477  const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
478  const ArrayView<Tensor<3, spacedim>> & output) const
479 {
480  AssertDimension(input.size(), output.size());
481  Assert((dynamic_cast<const typename MappingQ<dim, spacedim>::InternalData *>(
482  &mapping_data) != nullptr),
483  ExcInternalError());
484  const InternalData *data = dynamic_cast<const InternalData *>(&mapping_data);
485 
486  // check whether we should in fact work on the Q1 portion of it
488  q1_mapping->transform(input, mapping_type, *data->mapping_q1_data, output);
489  else
490  // otherwise use the full mapping
491  qp_mapping->transform(input, mapping_type, *data->mapping_qp_data, output);
492 }
493 
494 
495 
496 template <int dim, int spacedim>
499  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
500  const Point<dim> & p) const
501 {
502  // first see, whether we want to use a linear or a higher order
503  // mapping, then either use our own facilities or that of the Q1
504  // mapping we store
505  if (use_mapping_q_on_all_cells || cell->has_boundary_lines())
506  return qp_mapping->transform_unit_to_real_cell(cell, p);
507  else
508  return q1_mapping->transform_unit_to_real_cell(cell, p);
509 }
510 
511 
512 
513 template <int dim, int spacedim>
516  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
517  const Point<spacedim> & p) const
518 {
519  if (cell->has_boundary_lines() || use_mapping_q_on_all_cells ||
520  (dim != spacedim))
521  return qp_mapping->transform_real_to_unit_cell(cell, p);
522  else
523  return q1_mapping->transform_real_to_unit_cell(cell, p);
524 }
525 
526 
527 
528 template <int dim, int spacedim>
529 std::unique_ptr<Mapping<dim, spacedim>>
531 {
532  return std_cxx14::make_unique<MappingQ<dim, spacedim>>(
534 }
535 
536 
537 
538 // explicit instantiations
539 #include "mapping_q.inst"
540 
541 
542 DEAL_II_NAMESPACE_CLOSE
virtual std::size_t memory_consumption() const override
Definition: mapping_q.cc:50
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1366
std::shared_ptr< const MappingQGeneric< dim, spacedim > > qp_mapping
Definition: mapping_q.h:368
unsigned int get_degree() const
Definition: mapping_q.cc:123
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
Definition: mapping_q.cc:151
Task< RT > new_task(const std::function< RT()> &function)
MappingType
Definition: mapping.h:61
bool use_mapping_q1_on_current_cell
Definition: mapping_q.h:251
STL namespace.
const bool use_mapping_q_on_all_cells
Definition: mapping_q.h:333
std::unique_ptr< To > dynamic_unique_cast(std::unique_ptr< From > &&p)
Definition: utilities.h:630
#define Assert(cond, exc)
Definition: exceptions.h:1227
UpdateFlags
Abstract base class for mapping classes.
Definition: dof_tools.h:57
const unsigned int polynomial_degree
Definition: mapping_q.h:327
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_qp_data
Definition: mapping_q.h:265
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
Definition: mapping_q.cc:498
virtual bool preserves_vertex_locations() const override
Definition: mapping_q.cc:132
virtual void transform(const ArrayView< const Tensor< 1, dim >> &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim >> &output) const override
Definition: mapping_q.cc:379
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
Definition: mapping_q.cc:530
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
Definition: mapping_q.cc:180
std::unique_ptr< typename MappingQGeneric< dim, spacedim >::InternalData > mapping_q1_data
Definition: mapping_q.h:258
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
Definition: mapping_q.cc:211
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
Definition: mapping_q.cc:141
MappingQ(const unsigned int polynomial_degree, const bool use_mapping_q_on_all_cells=false)
Definition: mapping_q.cc:62
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
Definition: mapping_q.cc:515
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static::ExceptionBase & ExcInternalError()
std::shared_ptr< const MappingQGeneric< dim, spacedim > > q1_mapping
Definition: mapping_q.h:352