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