ASPECT
compat.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2015 - 2024 by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 #ifndef _aspect_compat_h
22 #define _aspect_compat_h
23 
24 #include <deal.II/base/config.h>
25 #include <deal.II/base/mpi.h>
26 
27 // C++11 related includes.
28 #include <functional>
29 #include <memory>
30 
31 #if !DEAL_II_VERSION_GTE(9,6,0)
32 # include <deal.II/multigrid/mg_transfer_matrix_free.h>
33 # include <deal.II/grid/manifold.h>
34 # include <deal.II/grid/manifold_lib.h>
35 #else
36 // We have a 'using' declaration for new enough deal.II versions. When
37 // we remove the backward compatibility, we should also remove the
38 // 'using' declaration and then can also remove the following
39 // #include:
40 # include <deal.II/grid/manifold_lib.h>
41 #endif
42 
43 #if !DEAL_II_VERSION_GTE(9,7,0)
44 # include <deal.II/grid/grid_generator.h>
45 #endif
46 
47 
48 namespace aspect
49 {
50  namespace big_mpi
51  {
52 
53  using ::Utilities::MPI::broadcast;
54 
55  }
56 
57 
58 // deal.II 9.6 introduces the new MGTransferMF class as a replacement
59 // for MGTransferMatrixFree. Instead of putting an ifdef in every place,
60 // do this in one central location:
61 #if !DEAL_II_VERSION_GTE(9,6,0)
62  template <int dim, class NumberType>
63  using MGTransferMF = ::MGTransferMatrixFree<dim,NumberType>;
64 #endif
65 
66 
67 #if !DEAL_II_VERSION_GTE(9,7,0)
68 
74  template <typename T, typename P = void>
75  using ObserverPointer = ::SmartPointer<T, P>;
76 
81 #endif
82 
83 
84 // deal.II versions up to 9.5 had a poorly designed interface of the
85 // SphericalManifold class that made it impossible for us to use.
86 // This file thus contains a copy of it.
87 #if !DEAL_II_VERSION_GTE(9,6,0)
88  using namespace dealii;
89 
98  template <int dim, int spacedim = dim>
99  class SphericalManifold : public Manifold<dim, spacedim>
100  {
101  public:
108  SphericalManifold(const Point<spacedim> center = Point<spacedim>());
109 
113  virtual std::unique_ptr<Manifold<dim, spacedim>>
114  clone() const override;
115 
123  virtual Point<spacedim>
124  get_intermediate_point(const Point<spacedim> &p1,
125  const Point<spacedim> &p2,
126  const double w) const override;
127 
132  virtual Tensor<1, spacedim>
133  get_tangent_vector(const Point<spacedim> &x1,
134  const Point<spacedim> &x2) const override;
135 
139  virtual Tensor<1, spacedim>
140  normal_vector(
141  const typename Triangulation<dim, spacedim>::face_iterator &face,
142  const Point<spacedim> &p) const override;
143 
147  virtual void
148  get_normals_at_vertices(
149  const typename Triangulation<dim, spacedim>::face_iterator &face,
150  typename Manifold<dim, spacedim>::FaceVertexNormals &face_vertex_normals)
151  const override;
152 
167  virtual void
168  get_new_points(const ArrayView<const Point<spacedim>> &surrounding_points,
169  const Table<2, double> &weights,
170  ArrayView<Point<spacedim>> new_points) const override;
171 
176  virtual Point<spacedim>
177  get_new_point(const ArrayView<const Point<spacedim>> &vertices,
178  const ArrayView<const double> &weights) const override;
179 
183  const Point<spacedim> center;
184 
185  private:
192  std::pair<double, Tensor<1, spacedim>>
193  guess_new_point(const ArrayView<const Tensor<1, spacedim>> &directions,
194  const ArrayView<const double> &distances,
195  const ArrayView<const double> &weights) const;
196 
214  void
215  do_get_new_points(const ArrayView<const Point<spacedim>> &surrounding_points,
216  const ArrayView<const double> &weights,
217  ArrayView<Point<spacedim>> new_points) const;
218 
222  const PolarManifold<spacedim> polar_manifold;
223  };
224 
225 #else
226 
227 // For sufficiently new deal.II versions, we can use the deal.II class, but to
228 // avoid name clashes, we have to import the class into namespace aspect. Once
229 // we rely on these sufficiently new versions of deal.II, we can not only remove
230 // the code above, but also the following lines, and in all places where we
231 // reference 'aspect::SphericalManifold' simply use 'SphericalManifold' instead
232 // (which then refers to the deal.II class).
233 
234  using ::SphericalManifold;
235 
236 #endif
237 
238 
239 // deal.II version 9.7 introduces a new class VectorFunctionFromTensorFunctionObject
240 // that we would like to use also for earlier versions
241 #if !DEAL_II_VERSION_GTE(9,7,0)
242 
243  using namespace dealii;
244 
279  template <int dim, typename RangeNumberType = double>
281  : public Function<dim, RangeNumberType>
282  {
283  public:
301  const std::function<Tensor<1, dim, RangeNumberType>(const Point<dim> &)>
302  &tensor_function_object,
303  const unsigned int selected_component = 0,
304  const unsigned int n_components = dim);
305 
310  virtual ~VectorFunctionFromTensorFunctionObject() override = default;
311 
315  virtual RangeNumberType
316  value(const Point<dim> &p, const unsigned int component = 0) const override;
317 
323  virtual void
324  vector_value(const Point<dim> &p,
325  Vector<RangeNumberType> &values) const override;
326 
334  virtual void
335  vector_value_list(
336  const std::vector<Point<dim>> &points,
337  std::vector<Vector<RangeNumberType>> &value_list) const override;
338 
339  private:
344  const std::function<Tensor<1, dim, RangeNumberType>(const Point<dim> &)>
345  tensor_function_object;
346 
353  const unsigned int selected_component;
354  };
355 
356 
357  template <int dim, typename RangeNumberType>
360  const std::function<Tensor<1, dim, RangeNumberType>(const Point<dim> &)>
361  &tensor_function_object,
362  const unsigned int selected_component,
363  const unsigned int n_components)
364  : Function<dim, RangeNumberType>(n_components)
365  , tensor_function_object(tensor_function_object)
366  , selected_component(selected_component)
367  {
368  // Verify that the Tensor<1,dim,RangeNumberType> will fit in the given length
369  // selected_components and not hang over the end of the vector.
370  AssertIndexRange(selected_component + dim - 1, this->n_components);
371  }
372 
373 
374 
375  template <int dim, typename RangeNumberType>
376  inline RangeNumberType
378  const Point<dim> &p,
379  const unsigned int component) const
380  {
381  AssertIndexRange(component, this->n_components);
382 
383  // if the requested component is out of the range selected, then we can
384  // return early
385  if ((component < selected_component) ||
386  (component >= selected_component + dim))
387  return 0;
388 
389  // otherwise retrieve the values from the <tt>tensor_function</tt> to be
390  // placed at the <tt>selected_component</tt> to
391  // <tt>selected_component + dim - 1</tt> elements of the <tt>Vector</tt>
392  // values and pick the correct one
393  const Tensor<1, dim, RangeNumberType> tensor_value =
394  tensor_function_object(p);
395 
396  return tensor_value[component - selected_component];
397  }
398 
399 
400  template <int dim, typename RangeNumberType>
401  inline void
403  const Point<dim> &p,
404  Vector<RangeNumberType> &values) const
405  {
406  Assert(values.size() == this->n_components,
407  ExcDimensionMismatch(values.size(), this->n_components));
408 
409  // Retrieve the values from the <tt>tensor_function</tt> to be placed at
410  // the <tt>selected_component</tt> to
411  // <tt>selected_component + dim - 1</tt> elements of the <tt>Vector</tt>
412  // values.
413  const Tensor<1, dim, RangeNumberType> tensor_value =
414  tensor_function_object(p);
415 
416  // First we make all elements of values = 0
417  values = 0;
418 
419  // Second we adjust the desired components to take on the values in
420  // <tt>tensor_value</tt>.
421  for (unsigned int i = 0; i < dim; ++i)
422  values(i + selected_component) = tensor_value[i];
423  }
424 
425 
433  template <int dim, typename RangeNumberType>
434  void
436  const std::vector<Point<dim>> &points,
437  std::vector<Vector<RangeNumberType>> &value_list) const
438  {
439  Assert(value_list.size() == points.size(),
440  ExcDimensionMismatch(value_list.size(), points.size()));
441 
442  const unsigned int n_points = points.size();
443 
444  for (unsigned int p = 0; p < n_points; ++p)
446  points[p], value_list[p]);
447  }
448 
449 #endif
450 
451 // deal.II versions up to 9.6 had a bug for very thin shell geometries.
452 // This function contains a fixed version.
453 #if !DEAL_II_VERSION_GTE(9,7,0)
454 
455  template <int dim>
456  void
457  colorize_quarter_hyper_shell(Triangulation<dim> &tria,
458  const Point<dim> &center,
459  const double inner_radius,
460  const double outer_radius);
461 #endif
462 
463  // deal.II 9.8 made ReferenceCell a template class, whereas older versions
464  // had it as a non-template class. This is a problem.
465  // Rather than litter our own code base with #ifdefs, we can just define the
466  // templated class variant here for older deal.II versions, and then we can
467  // use the same code in all versions.
468 #if !DEAL_II_VERSION_GTE(9,8,0)
469  template <int dim> using ReferenceCell = ::ReferenceCell;
470 #endif
471 
472 }
473 
474 #endif
::SmartPointer< T, P > ObserverPointer
Definition: compat.h:75
VectorFunctionFromTensorFunctionObject(const std::function< Tensor< 1, dim, RangeNumberType >(const Point< dim > &)> &tensor_function_object, const unsigned int selected_component=0, const unsigned int n_components=dim)
Definition: compat.h:359
::MGTransferMatrixFree< dim, NumberType > MGTransferMF
Definition: compat.h:63
const PolarManifold< spacedim > polar_manifold
Definition: compat.h:222
void colorize_quarter_hyper_shell(Triangulation< dim > &tria, const Point< dim > &center, const double inner_radius, const double outer_radius)
const Point< spacedim > center
Definition: compat.h:183
::ReferenceCell ReferenceCell
Definition: compat.h:469