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 
68 // deal.II versions up to 9.5 had a poorly designed interface of the
69 // SphericalManifold class that made it impossible for us to use.
70 // This file thus contains a copy of it.
71 #if !DEAL_II_VERSION_GTE(9,6,0)
72  using namespace dealii;
73 
82  template <int dim, int spacedim = dim>
83  class SphericalManifold : public Manifold<dim, spacedim>
84  {
85  public:
92  SphericalManifold(const Point<spacedim> center = Point<spacedim>());
93 
97  virtual std::unique_ptr<Manifold<dim, spacedim>>
98  clone() const override;
99 
107  virtual Point<spacedim>
108  get_intermediate_point(const Point<spacedim> &p1,
109  const Point<spacedim> &p2,
110  const double w) const override;
111 
116  virtual Tensor<1, spacedim>
117  get_tangent_vector(const Point<spacedim> &x1,
118  const Point<spacedim> &x2) const override;
119 
123  virtual Tensor<1, spacedim>
124  normal_vector(
125  const typename Triangulation<dim, spacedim>::face_iterator &face,
126  const Point<spacedim> &p) const override;
127 
131  virtual void
132  get_normals_at_vertices(
133  const typename Triangulation<dim, spacedim>::face_iterator &face,
134  typename Manifold<dim, spacedim>::FaceVertexNormals &face_vertex_normals)
135  const override;
136 
151  virtual void
152  get_new_points(const ArrayView<const Point<spacedim>> &surrounding_points,
153  const Table<2, double> &weights,
154  ArrayView<Point<spacedim>> new_points) const override;
155 
160  virtual Point<spacedim>
161  get_new_point(const ArrayView<const Point<spacedim>> &vertices,
162  const ArrayView<const double> &weights) const override;
163 
167  const Point<spacedim> center;
168 
169  private:
176  std::pair<double, Tensor<1, spacedim>>
177  guess_new_point(const ArrayView<const Tensor<1, spacedim>> &directions,
178  const ArrayView<const double> &distances,
179  const ArrayView<const double> &weights) const;
180 
198  void
199  do_get_new_points(const ArrayView<const Point<spacedim>> &surrounding_points,
200  const ArrayView<const double> &weights,
201  ArrayView<Point<spacedim>> new_points) const;
202 
206  const PolarManifold<spacedim> polar_manifold;
207  };
208 
209 #else
210 
211 // For sufficiently new deal.II versions, we can use the deal.II class, but to
212 // avoid name clashes, we have to import the class into namespace aspect. Once
213 // we rely on these sufficiently new versions of deal.II, we can not only remove
214 // the code above, but also the following lines, and in all places where we
215 // reference 'aspect::SphericalManifold' simply use 'SphericalManifold' instead
216 // (which then refers to the deal.II class).
217 
218  using ::SphericalManifold;
219 
220 #endif
221 
222 
223 // deal.II version 9.7 introduces a new class VectorFunctionFromTensorFunctionObject
224 // that we would like to use also for earlier versions
225 #if !DEAL_II_VERSION_GTE(9,7,0)
226 
227  using namespace dealii;
228 
263  template <int dim, typename RangeNumberType = double>
265  : public Function<dim, RangeNumberType>
266  {
267  public:
285  const std::function<Tensor<1, dim, RangeNumberType>(const Point<dim> &)>
286  &tensor_function_object,
287  const unsigned int selected_component = 0,
288  const unsigned int n_components = dim);
289 
294  virtual ~VectorFunctionFromTensorFunctionObject() override = default;
295 
299  virtual RangeNumberType
300  value(const Point<dim> &p, const unsigned int component = 0) const override;
301 
307  virtual void
308  vector_value(const Point<dim> &p,
309  Vector<RangeNumberType> &values) const override;
310 
318  virtual void
319  vector_value_list(
320  const std::vector<Point<dim>> &points,
321  std::vector<Vector<RangeNumberType>> &value_list) const override;
322 
323  private:
328  const std::function<Tensor<1, dim, RangeNumberType>(const Point<dim> &)>
329  tensor_function_object;
330 
337  const unsigned int selected_component;
338  };
339 
340 
341  template <int dim, typename RangeNumberType>
344  const std::function<Tensor<1, dim, RangeNumberType>(const Point<dim> &)>
345  &tensor_function_object,
346  const unsigned int selected_component,
347  const unsigned int n_components)
348  : Function<dim, RangeNumberType>(n_components)
349  , tensor_function_object(tensor_function_object)
350  , selected_component(selected_component)
351  {
352  // Verify that the Tensor<1,dim,RangeNumberType> will fit in the given length
353  // selected_components and not hang over the end of the vector.
354  AssertIndexRange(selected_component + dim - 1, this->n_components);
355  }
356 
357 
358 
359  template <int dim, typename RangeNumberType>
360  inline RangeNumberType
362  const Point<dim> &p,
363  const unsigned int component) const
364  {
365  AssertIndexRange(component, this->n_components);
366 
367  // if the requested component is out of the range selected, then we can
368  // return early
369  if ((component < selected_component) ||
370  (component >= selected_component + dim))
371  return 0;
372 
373  // otherwise retrieve the values from the <tt>tensor_function</tt> to be
374  // placed at the <tt>selected_component</tt> to
375  // <tt>selected_component + dim - 1</tt> elements of the <tt>Vector</tt>
376  // values and pick the correct one
377  const Tensor<1, dim, RangeNumberType> tensor_value =
378  tensor_function_object(p);
379 
380  return tensor_value[component - selected_component];
381  }
382 
383 
384  template <int dim, typename RangeNumberType>
385  inline void
387  const Point<dim> &p,
388  Vector<RangeNumberType> &values) const
389  {
390  Assert(values.size() == this->n_components,
391  ExcDimensionMismatch(values.size(), this->n_components));
392 
393  // Retrieve the values from the <tt>tensor_function</tt> to be placed at
394  // the <tt>selected_component</tt> to
395  // <tt>selected_component + dim - 1</tt> elements of the <tt>Vector</tt>
396  // values.
397  const Tensor<1, dim, RangeNumberType> tensor_value =
398  tensor_function_object(p);
399 
400  // First we make all elements of values = 0
401  values = 0;
402 
403  // Second we adjust the desired components to take on the values in
404  // <tt>tensor_value</tt>.
405  for (unsigned int i = 0; i < dim; ++i)
406  values(i + selected_component) = tensor_value[i];
407  }
408 
409 
417  template <int dim, typename RangeNumberType>
418  void
420  const std::vector<Point<dim>> &points,
421  std::vector<Vector<RangeNumberType>> &value_list) const
422  {
423  Assert(value_list.size() == points.size(),
424  ExcDimensionMismatch(value_list.size(), points.size()));
425 
426  const unsigned int n_points = points.size();
427 
428  for (unsigned int p = 0; p < n_points; ++p)
430  points[p], value_list[p]);
431  }
432 
433 #endif
434 
435 // deal.II versions up to 9.6 had a bug for very thin shell geometries.
436 // This function contains a fixed version.
437 #if !DEAL_II_VERSION_GTE(9,7,0)
438 
439  template <int dim>
440  void
441  colorize_quarter_hyper_shell(Triangulation<dim> &tria,
442  const Point<dim> &center,
443  const double inner_radius,
444  const double outer_radius);
445 #endif
446 
447 }
448 
449 #endif
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:343
::MGTransferMatrixFree< dim, NumberType > MGTransferMF
Definition: compat.h:63
const PolarManifold< spacedim > polar_manifold
Definition: compat.h:206
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:167