]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Describe time by real scalar types in Function and friends
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 26 Oct 2018 09:28:51 +0000 (11:28 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 26 Oct 2018 09:35:29 +0000 (11:35 +0200)
doc/news/changes/incompatibilities/20181026DanielArndt [new file with mode: 0644]
include/deal.II/base/function.h
include/deal.II/base/function.templates.h
include/deal.II/base/tensor_function.h
include/deal.II/base/tensor_function.templates.h
tests/base/functions_13.cc

diff --git a/doc/news/changes/incompatibilities/20181026DanielArndt b/doc/news/changes/incompatibilities/20181026DanielArndt
new file mode 100644 (file)
index 0000000..7a45d4d
--- /dev/null
@@ -0,0 +1,5 @@
+Changed: Time is described by a real-valued scalar in Function, TensorFunction,
+ConstantTensorFunction and ZeroTensorFunction, i.e. all of them derive from
+FunctionTime with a real-valued scalar template parameter.
+<br>
+(Daniel Arndt, 2018/10/26)
index 16fcdbaf8936ff8b38966fdc8ab6a9deec11ced7..952ec34c1e347d5291400602be91b1bd9bc0afcf 100644 (file)
@@ -145,7 +145,9 @@ class TensorFunction;
  * @author Wolfgang Bangerth, 1998, 1999, Luca Heltai 2014
  */
 template <int dim, typename RangeNumberType = double>
-class Function : public FunctionTime<RangeNumberType>, public Subscriptor
+class Function : public FunctionTime<
+                   typename numbers::NumberTraits<RangeNumberType>::real_type>,
+                 public Subscriptor
 {
 public:
   /**
@@ -164,8 +166,9 @@ public:
    * (which defaults to one, i.e. a scalar function), and the time variable,
    * which defaults to zero.
    */
-  Function(const unsigned int    n_components = 1,
-           const RangeNumberType initial_time = 0.0);
+  Function(const unsigned int n_components = 1,
+           const typename numbers::NumberTraits<RangeNumberType>::real_type
+             initial_time = 0.0);
 
   /**
    * Virtual destructor; absolutely necessary in this case.
index 3dbbb7245c3a4a60c4081bc1458710f5502dfabc..a2d57e8c35235a84c26c957a542312a614f008ac 100644 (file)
@@ -33,9 +33,11 @@ const unsigned int Function<dim, RangeNumberType>::dimension;
 
 
 template <int dim, typename RangeNumberType>
-Function<dim, RangeNumberType>::Function(const unsigned int    n_components,
-                                         const RangeNumberType initial_time)
-  : FunctionTime<RangeNumberType>(initial_time)
+Function<dim, RangeNumberType>::Function(
+  const unsigned int                                               n_components,
+  const typename numbers::NumberTraits<RangeNumberType>::real_type initial_time)
+  : FunctionTime<typename numbers::NumberTraits<RangeNumberType>::real_type>(
+      initial_time)
   , n_components(n_components)
 {
   // avoid the construction of function objects that don't return any
index 8b80dae6ca36021fc675076c39ccf96f813e68cc..564b4c59634ae3fa6785f63957f50fc7eb94f0d4 100644 (file)
@@ -54,7 +54,9 @@ DEAL_II_NAMESPACE_OPEN
  * @author Guido Kanschat, 1999
  */
 template <int rank, int dim, typename Number = double>
-class TensorFunction : public FunctionTime<Number>, public Subscriptor
+class TensorFunction
+  : public FunctionTime<typename numbers::NumberTraits<Number>::real_type>,
+    public Subscriptor
 {
 public:
   /**
@@ -68,7 +70,9 @@ public:
    * Constructor. May take an initial value for the time variable, which
    * defaults to zero.
    */
-  TensorFunction(const Number initial_time = Number(0.0));
+  TensorFunction(
+    const typename numbers::NumberTraits<Number>::real_type initial_time =
+      typename numbers::NumberTraits<Number>::real_type(0.0));
 
   /**
    * Virtual destructor; absolutely necessary in this case, as classes are
@@ -128,8 +132,9 @@ public:
    * An initial value for the time variable may be specified, otherwise it
    * defaults to zero.
    */
-  ConstantTensorFunction(const dealii::Tensor<rank, dim, Number> &value,
-                         const Number initial_time = 0.0);
+  ConstantTensorFunction(
+    const dealii::Tensor<rank, dim, Number> &               value,
+    const typename numbers::NumberTraits<Number>::real_type initial_time = 0.0);
 
   virtual ~ConstantTensorFunction() override = default;
 
@@ -175,7 +180,8 @@ public:
    * An initial value for the time variable may be specified, otherwise it
    * defaults to zero.
    */
-  ZeroTensorFunction(const Number initial_time = 0.0);
+  ZeroTensorFunction(
+    const typename numbers::NumberTraits<Number>::real_type initial_time = 0.0);
 };
 
 
index 5e5b0c1992d1ea5b302620abc21a8bd6a64cc4d9..398a7f645272099df113a9dd604a79b7b8d4d14d 100644 (file)
@@ -28,8 +28,10 @@ DEAL_II_NAMESPACE_OPEN
 
 
 template <int rank, int dim, typename Number>
-TensorFunction<rank, dim, Number>::TensorFunction(const Number initial_time)
-  : FunctionTime<Number>(initial_time)
+TensorFunction<rank, dim, Number>::TensorFunction(
+  const typename numbers::NumberTraits<Number>::real_type initial_time)
+  : FunctionTime<typename numbers::NumberTraits<Number>::real_type>(
+      initial_time)
 {}
 
 
@@ -83,8 +85,8 @@ TensorFunction<rank, dim, Number>::gradient_list(
 
 template <int rank, int dim, typename Number>
 ConstantTensorFunction<rank, dim, Number>::ConstantTensorFunction(
-  const Tensor<rank, dim, Number> &value,
-  const Number                     initial_time)
+  const Tensor<rank, dim, Number> &                       value,
+  const typename numbers::NumberTraits<Number>::real_type initial_time)
   : TensorFunction<rank, dim, Number>(initial_time)
   , _value(value)
 {}
@@ -147,7 +149,7 @@ ConstantTensorFunction<rank, dim, Number>::gradient_list(
 
 template <int rank, int dim, typename Number>
 ZeroTensorFunction<rank, dim, Number>::ZeroTensorFunction(
-  const Number initial_time)
+  const typename numbers::NumberTraits<Number>::real_type initial_time)
   : ConstantTensorFunction<rank, dim, Number>(
       dealii::Tensor<rank, dim, Number>(),
       initial_time)
index 7b1b7e9bd81047f0e58afa51ec44230b945b9ada..29e7c8b5e9bcc5955acc03c4d12001536b6bbf3e 100644 (file)
 
 #include "../tests.h"
 
-using namespace dealii;
-
 template <typename Number>
-class MyFunctionTime : dealii::FunctionTime<Number>
+class MyFunctionTime : public dealii::FunctionTime<Number>
 {};
 
 template <int dim, typename Number>
-class MyFunction : dealii::Function<dim, Number>
+class MyFunction : public dealii::Function<dim, Number>
 {};
 
 template <int rank, int dim, typename Number>
-class MyTensorFunction : dealii::TensorFunction<rank, dim, Number>
+class MyTensorFunction : public dealii::TensorFunction<rank, dim, Number>
 {};
 
 int
@@ -42,11 +40,30 @@ main()
 {
   initlog();
 
+  MyFunctionTime<std::complex<double>> function_time_1;
+  static_assert(std::is_same<decltype(function_time_1.get_time()),
+                             std::complex<double>>::value,
+                "Wrong return type for FunctionTime<double>::get_time().");
+  MyFunctionTime<std::complex<float>> function_time_2;
+  static_assert(std::is_same<decltype(function_time_2.get_time()),
+                             std::complex<float>>::value,
+                "Wrong return type for FunctionTime<float>::get_time().");
+
   MyFunction<1, std::complex<double>> function_1;
-  MyFunction<2, std::complex<float>>  function_2;
+  static_assert(std::is_same<decltype(function_1.get_time()), double>::value,
+                "Wrong return type for Function<double>::get_time().");
+  MyFunction<2, std::complex<float>> function_2;
+  static_assert(std::is_same<decltype(function_2.get_time()), float>::value,
+                "Wrong return type for Function<float>::get_time().");
 
   MyTensorFunction<1, 1, std::complex<double>> tensor_function_1;
-  MyTensorFunction<2, 2, std::complex<float>>  tensor_function_2;
+  static_assert(
+    std::is_same<decltype(tensor_function_1.get_time()), double>::value,
+    "Wrong return type for TensorFunction<double>::get_time().");
+  MyTensorFunction<2, 2, std::complex<float>> tensor_function_2;
+  static_assert(
+    std::is_same<decltype(tensor_function_2.get_time()), float>::value,
+    "Wrong return type for TensorFunction<float>::get_time().");
 
   deallog << "OK" << std::endl;
 }

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.