]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a test 13133/head
authorJean-Paul Pelteret <jppelteret@gmail.com>
Mon, 27 Dec 2021 10:46:28 +0000 (11:46 +0100)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Mon, 27 Dec 2021 10:47:03 +0000 (11:47 +0100)
tests/symengine/batch_optimizer_09_1.cc [new file with mode: 0644]
tests/symengine/batch_optimizer_09_1.output [new file with mode: 0644]
tests/symengine/sd_common_tests/batch_optimizer_09.h [new file with mode: 0644]

diff --git a/tests/symengine/batch_optimizer_09_1.cc b/tests/symengine/batch_optimizer_09_1.cc
new file mode 100644 (file)
index 0000000..f6598bf
--- /dev/null
@@ -0,0 +1,45 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// Check that assertions for some illegal operations when using the
+// BatchOptimizer get triggered.
+//
+// Here we use only dictionary substitution, and invoke no symbolic
+// optimizations.
+
+#include "../tests.h"
+
+#include "sd_common_tests/batch_optimizer_09.h"
+
+int
+main()
+{
+  deal_II_exceptions::disable_abort_on_exception();
+
+  initlog();
+  deallog << std::setprecision(10);
+
+  // Show the difference between a SymEngine "value" and
+  // an evaluated, floating point number
+  // deallog << std::setprecision(3);
+
+  const enum SD::OptimizerType     opt_method = SD::OptimizerType::dictionary;
+  const enum SD::OptimizationFlags opt_flags =
+    SD::OptimizationFlags::optimize_all;
+
+  const int dim = 2;
+  run_tests<dim, opt_method, opt_flags>();
+}
diff --git a/tests/symengine/batch_optimizer_09_1.output b/tests/symengine/batch_optimizer_09_1.output
new file mode 100644 (file)
index 0000000..f927e52
--- /dev/null
@@ -0,0 +1,24 @@
+
+DEAL:Double:Evaluation::f: x*y
+DEAL:Double:Evaluation::G: (x + y)*T_00 (x + y)*T_01 (x + y)*T_10 (x + y)*T_11
+DEAL:Double:Evaluation::H: S_00/(x - y) S_01/(x - y) S_01/(x - y) S_11/(x - y)
+DEAL:Double:Evaluation::f(x,y): 20.00000000
+DEAL:Double:Evaluation::G(x,y,T): 48.00000000 48.00000000 60.00000000 84.00000000
+DEAL:Double:Evaluation::H(x,y,S): 0.6250000000 0.7500000000 0.7500000000 0.7500000000
+DEAL:Double:Evaluation::
+DEAL:Double:Extraction::f(x,y): 20.00000000
+DEAL:Double:Extraction::G(x,y,T): 48.00000000 48.00000000 60.00000000 84.00000000
+DEAL:Double:Extraction::H(x,y,S): 0.6250000000 0.7500000000 0.7500000000 0.7500000000
+DEAL:Double:Extraction::
+DEAL:Double:Invalid evaluation::Going to try evaluation without substitution...
+DEAL:Double:Invalid evaluation::Caught invalid evaluation before substitution: Scalar
+DEAL:Double:Invalid evaluation::Caught invalid evaluation before substitution: Tensor
+DEAL:Double:Invalid evaluation::Caught invalid evaluation before substitution: SymmetricTensor
+DEAL:Double:Invalid evaluation::Going to try substitution without optimization...
+DEAL:Double:Invalid evaluation::Caught invalid substitution before optimization.
+DEAL:Double:Invalid evaluation::Going to try evaluation without optimization...
+DEAL:Double:Invalid evaluation::Caught invalid evaluation before optimization: Scalar
+DEAL:Double:Invalid evaluation::Caught invalid evaluation before optimization: Tensor
+DEAL:Double:Invalid evaluation::Caught invalid evaluation before optimization: SymmetricTensor
+DEAL:Double::OK
+DEAL::OK
diff --git a/tests/symengine/sd_common_tests/batch_optimizer_09.h b/tests/symengine/sd_common_tests/batch_optimizer_09.h
new file mode 100644 (file)
index 0000000..75ccd32
--- /dev/null
@@ -0,0 +1,225 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2020 - 2021 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// Check that assertions for some illegal operations when using the
+// BatchOptimizer get triggered.
+
+#include <deal.II/differentiation/sd.h>
+
+#include "../../tests.h"
+
+#include "utilities.h"
+
+using namespace dealii;
+namespace SD = Differentiation::SD;
+
+template <int dim,
+          typename NumberType,
+          enum SD::OptimizerType     opt_method,
+          enum SD::OptimizationFlags opt_flags>
+void
+test_extract_evaluate()
+{
+  namespace SD           = Differentiation::SD;
+  using SD_number_t      = SD::Expression;
+  using SD_tensor_t      = Tensor<2, dim, SD_number_t>;
+  using SD_symm_tensor_t = SymmetricTensor<2, dim, SD_number_t>;
+
+  // Define
+  const SD_number_t      x("x"), y("y");
+  const SD_tensor_t      T = SD::make_tensor_of_symbols<2, dim>("T");
+  const SD_symm_tensor_t S = SD::make_symmetric_tensor_of_symbols<2, dim>("S");
+  const SD_number_t      f = x * y;
+  const SD_tensor_t      G = T * (x + y);
+  const SD_symm_tensor_t H = S / (x - y);
+
+  // Substitution map
+  const SD::types::substitution_map sub_vals = SD::make_substitution_map(
+    std::make_pair(x, NumberType(10.0)),
+    std::make_pair(y, NumberType(2.0)),
+    std::make_pair(T, make_tensor<dim>(NumberType(3.0))),
+    std::make_pair(S, make_symm_tensor<dim>(NumberType(4.0))));
+
+  // Optimize
+  SD::BatchOptimizer<NumberType> optimizer(opt_method, opt_flags);
+  optimizer.register_symbols(sub_vals);  // Independent symbols
+  optimizer.register_functions(f, G, H); // Dependent functions
+  optimizer.optimize();
+
+  // Substitute
+  optimizer.substitute(sub_vals);
+
+  // Evaluate
+  deallog.push("Evaluation");
+  {
+    const NumberType                          val_f = optimizer.evaluate(f);
+    const Tensor<2, dim, NumberType>          val_G = optimizer.evaluate(G);
+    const SymmetricTensor<2, dim, NumberType> val_H = optimizer.evaluate(H);
+    deallog << "f: " << f << std::endl
+            << "G: " << G << std::endl
+            << "H: " << H << std::endl
+            << "f(x,y): " << val_f << std::endl
+            << "G(x,y,T): " << val_G << std::endl
+            << "H(x,y,S): " << val_H << std::endl
+            << std::endl;
+  }
+  deallog.pop();
+
+  // Now we setup a second optimizer, which will only perform value extraction
+  SD::BatchOptimizer<NumberType> extraction_optimizer(opt_method, opt_flags);
+  extraction_optimizer.register_functions(f, G, H); // Dependent functions
+
+  // Extract
+  deallog.push("Extraction");
+  {
+    const auto values = optimizer.evaluate();
+
+    const NumberType extr_val_f = extraction_optimizer.extract(f, values);
+    const Tensor<2, dim, NumberType> extr_val_G =
+      extraction_optimizer.extract(G, values);
+    const SymmetricTensor<2, dim, NumberType> extr_val_H =
+      extraction_optimizer.extract(H, values);
+    deallog << "f(x,y): " << extr_val_f << std::endl
+            << "G(x,y,T): " << extr_val_G << std::endl
+            << "H(x,y,S): " << extr_val_H << std::endl
+            << std::endl;
+  }
+  deallog.pop();
+
+  // Attempt evaluation with extraction-specific optimiser. We expect this to
+  // fail since extraction_optimizer.optimize() hasn't been called.
+  deallog.push("Invalid evaluation");
+  {
+    extraction_optimizer.register_symbols(sub_vals);
+
+    deallog << "Going to try evaluation without substitution..." << std::endl;
+
+    try
+      {
+        const NumberType val_f = extraction_optimizer.evaluate(f);
+        deallog << "f(x,y): " << val_f << std::endl << std::endl;
+      }
+    catch (std::exception &)
+      {
+        deallog << "Caught invalid evaluation before substitution: Scalar"
+                << std::endl;
+      }
+
+    try
+      {
+        const Tensor<2, dim, NumberType> val_G =
+          extraction_optimizer.evaluate(G);
+        deallog << "G(x,y,T): " << val_G << std::endl << std::endl;
+      }
+    catch (std::exception &)
+      {
+        deallog << "Caught invalid evaluation before substitution: Tensor"
+                << std::endl;
+      }
+
+    try
+      {
+        const SymmetricTensor<2, dim, NumberType> val_H =
+          extraction_optimizer.evaluate(H);
+        deallog << "H(x,y,S): " << val_H << std::endl << std::endl;
+      }
+    catch (std::exception &)
+      {
+        deallog
+          << "Caught invalid evaluation before substitution: SymmetricTensor"
+          << std::endl;
+      }
+
+    deallog << "Going to try substitution without optimization..." << std::endl;
+
+    try
+      {
+        extraction_optimizer.substitute(sub_vals);
+      }
+    catch (std::exception &)
+      {
+        deallog << "Caught invalid substitution before optimization."
+                << std::endl;
+      }
+
+    deallog << "Going to try evaluation without optimization..." << std::endl;
+
+    try
+      {
+        const NumberType val_f = extraction_optimizer.evaluate(f);
+        deallog << "f(x,y): " << val_f << std::endl << std::endl;
+      }
+    catch (std::exception &)
+      {
+        deallog << "Caught invalid evaluation before optimization: Scalar"
+                << std::endl;
+      }
+
+    try
+      {
+        const Tensor<2, dim, NumberType> val_G =
+          extraction_optimizer.evaluate(G);
+        deallog << "G(x,y,T): " << val_G << std::endl << std::endl;
+      }
+    catch (std::exception &)
+      {
+        deallog << "Caught invalid evaluation before optimization: Tensor"
+                << std::endl;
+      }
+
+    try
+      {
+        const SymmetricTensor<2, dim, NumberType> val_H =
+          extraction_optimizer.evaluate(H);
+        deallog << "H(x,y,S): " << val_H << std::endl << std::endl;
+      }
+    catch (std::exception &)
+      {
+        deallog
+          << "Caught invalid evaluation before optimization: SymmetricTensor"
+          << std::endl;
+      }
+  }
+  deallog.pop();
+
+  deallog << "OK" << std::endl;
+}
+
+
+template <int                        dim,
+          enum SD::OptimizerType     opt_method,
+          enum SD::OptimizationFlags opt_flags>
+void
+run_tests()
+{
+  // Show the difference between a SymEngine "value" and
+  // an evaluated, floating point number
+  // deallog << std::setprecision(3);
+
+  deallog.push("Double");
+  try
+    {
+      test_extract_evaluate<dim, double, opt_method, opt_flags>();
+    }
+  catch (...)
+    {
+      deallog.pop();
+      deallog << "NOT OK" << std::endl;
+    }
+  deallog.pop();
+
+  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.