--- /dev/null
+New: The SD::BatchOptimizer::copy_from() function can be used to duplicate the
+data from one batch optimizer instance to another.
+<br>
+(Jean-Paul Pelteret, 2021/11/27)
*/
~BatchOptimizer() = default;
+ /**
+ * Duplicate the data stored in an @p other BatchOptimizer instance.
+ *
+ * @note The optimized data and results from previous substitutions
+ * executed by the @p other optimizer instance are not copied over.
+ * It is therefore necessary to call optimize() before it is possible to
+ * substitute() values and evaluate() data. One may, however, still
+ * extract() values using @p this optimizer instance if those results are
+ * stored elsewhere.
+ */
+ void
+ copy_from(const BatchOptimizer &other);
+
/**
* Print some information on state of the internal data
* structures stored in the class.
+ template <typename ReturnType>
+ void
+ BatchOptimizer<ReturnType>::copy_from(
+ const BatchOptimizer<ReturnType> &other)
+ {
+ method = other.method;
+ flags = other.flags;
+ independent_variables_symbols = other.independent_variables_symbols;
+ dependent_variables_functions = other.dependent_variables_functions;
+ dependent_variables_output.clear();
+ map_dep_expr_vec_entry = other.map_dep_expr_vec_entry;
+ ready_for_value_extraction = false;
+ has_been_serialized = false;
+ }
+
+
+
template <typename ReturnType>
void
BatchOptimizer<ReturnType>::set_optimization_method(
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 the BatchOptimizer::copy_from() function works as expected.
+//
+// Here we use only dictionary substitution, and invoke no symbolic
+// optimizations.
+
+#include "../tests.h"
+
+#include "sd_common_tests/batch_optimizer_08.h"
+
+int
+main()
+{
+ 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>();
+}
--- /dev/null
+
+DEAL:Copy:Float::f: x*y
+DEAL:Copy:Float::g: x/y
+DEAL:Copy:Float::f(x,y): 20.00000000
+DEAL:Copy:Float::g(x,y): 5.000000000
+DEAL:Copy:Float::
+DEAL:Copy:Float::Before copy
+DEAL:Copy:Float::f: x1*y1
+DEAL:Copy:Float::g: x1/y1
+DEAL:Copy:Float::f(x,y): 60.00000000
+DEAL:Copy:Float::g(x,y): 6.666666508
+DEAL:Copy:Float::
+DEAL:Copy:Float::After copy
+DEAL:Copy:Float::f: x*y
+DEAL:Copy:Float::g: x/y
+DEAL:Copy:Float::f(x,y): 20.00000000
+DEAL:Copy:Float::g(x,y): 5.000000000
+DEAL:Copy:Float::
+DEAL:Copy:Float::OK
+DEAL:Copy:Double::f: x*y
+DEAL:Copy:Double::g: x/y
+DEAL:Copy:Double::f(x,y): 20.00000000
+DEAL:Copy:Double::g(x,y): 5.000000000
+DEAL:Copy:Double::
+DEAL:Copy:Double::Before copy
+DEAL:Copy:Double::f: x1*y1
+DEAL:Copy:Double::g: x1/y1
+DEAL:Copy:Double::f(x,y): 60.00000000
+DEAL:Copy:Double::g(x,y): 6.666666667
+DEAL:Copy:Double::
+DEAL:Copy:Double::After copy
+DEAL:Copy:Double::f: x*y
+DEAL:Copy:Double::g: x/y
+DEAL:Copy:Double::f(x,y): 20.00000000
+DEAL:Copy:Double::g(x,y): 5.000000000
+DEAL:Copy:Double::
+DEAL:Copy:Double::OK
+DEAL:Copy:Complex float::f: x*y
+DEAL:Copy:Complex float::g: x/y
+DEAL:Copy:Complex float::f(x,y): (20.00000000,0.000000000)
+DEAL:Copy:Complex float::g(x,y): (5.000000000,0.000000000)
+DEAL:Copy:Complex float::
+DEAL:Copy:Complex float::Before copy
+DEAL:Copy:Complex float::f: x1*y1
+DEAL:Copy:Complex float::g: x1/y1
+DEAL:Copy:Complex float::f(x,y): (60.00000000,0.000000000)
+DEAL:Copy:Complex float::g(x,y): (6.666666508,0.000000000)
+DEAL:Copy:Complex float::
+DEAL:Copy:Complex float::After copy
+DEAL:Copy:Complex float::f: x*y
+DEAL:Copy:Complex float::g: x/y
+DEAL:Copy:Complex float::f(x,y): (20.00000000,0.000000000)
+DEAL:Copy:Complex float::g(x,y): (5.000000000,0.000000000)
+DEAL:Copy:Complex float::
+DEAL:Copy:Complex float::OK
+DEAL:Copy:Complex double::f: x*y
+DEAL:Copy:Complex double::g: x/y
+DEAL:Copy:Complex double::f(x,y): (20.00000000,0.000000000)
+DEAL:Copy:Complex double::g(x,y): (5.000000000,0.000000000)
+DEAL:Copy:Complex double::
+DEAL:Copy:Complex double::Before copy
+DEAL:Copy:Complex double::f: x1*y1
+DEAL:Copy:Complex double::g: x1/y1
+DEAL:Copy:Complex double::f(x,y): (60.00000000,0.000000000)
+DEAL:Copy:Complex double::g(x,y): (6.666666667,0.000000000)
+DEAL:Copy:Complex double::
+DEAL:Copy:Complex double::After copy
+DEAL:Copy:Complex double::f: x*y
+DEAL:Copy:Complex double::g: x/y
+DEAL:Copy:Complex double::f(x,y): (20.00000000,0.000000000)
+DEAL:Copy:Complex double::g(x,y): (5.000000000,0.000000000)
+DEAL:Copy:Complex double::
+DEAL:Copy:Complex double::OK
+DEAL::OK
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 the BatchOptimizer::copy_from() function works as expected.
+//
+// Here we use only lambda substitution, and invoke no symbolic optimizations.
+
+#include "../tests.h"
+
+#include "sd_common_tests/batch_optimizer_08.h"
+
+int
+main()
+{
+ 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::lambda;
+ const enum SD::OptimizationFlags opt_flags =
+ SD::OptimizationFlags::optimize_default;
+
+ const int dim = 2;
+ run_tests<dim, opt_method, opt_flags>();
+}
--- /dev/null
+
+DEAL:Copy:Float::f: x*y
+DEAL:Copy:Float::g: x/y
+DEAL:Copy:Float::f(x,y): 20.00000000
+DEAL:Copy:Float::g(x,y): 5.000000000
+DEAL:Copy:Float::
+DEAL:Copy:Float::Before copy
+DEAL:Copy:Float::f: x1*y1
+DEAL:Copy:Float::g: x1/y1
+DEAL:Copy:Float::f(x,y): 60.00000000
+DEAL:Copy:Float::g(x,y): 6.666666508
+DEAL:Copy:Float::
+DEAL:Copy:Float::After copy
+DEAL:Copy:Float::f: x*y
+DEAL:Copy:Float::g: x/y
+DEAL:Copy:Float::f(x,y): 20.00000000
+DEAL:Copy:Float::g(x,y): 5.000000000
+DEAL:Copy:Float::
+DEAL:Copy:Float::OK
+DEAL:Copy:Double::f: x*y
+DEAL:Copy:Double::g: x/y
+DEAL:Copy:Double::f(x,y): 20.00000000
+DEAL:Copy:Double::g(x,y): 5.000000000
+DEAL:Copy:Double::
+DEAL:Copy:Double::Before copy
+DEAL:Copy:Double::f: x1*y1
+DEAL:Copy:Double::g: x1/y1
+DEAL:Copy:Double::f(x,y): 60.00000000
+DEAL:Copy:Double::g(x,y): 6.666666667
+DEAL:Copy:Double::
+DEAL:Copy:Double::After copy
+DEAL:Copy:Double::f: x*y
+DEAL:Copy:Double::g: x/y
+DEAL:Copy:Double::f(x,y): 20.00000000
+DEAL:Copy:Double::g(x,y): 5.000000000
+DEAL:Copy:Double::
+DEAL:Copy:Double::OK
+DEAL:Copy:Complex float::f: x*y
+DEAL:Copy:Complex float::g: x/y
+DEAL:Copy:Complex float::f(x,y): (20.00000000,0.000000000)
+DEAL:Copy:Complex float::g(x,y): (5.000000000,0.000000000)
+DEAL:Copy:Complex float::
+DEAL:Copy:Complex float::Before copy
+DEAL:Copy:Complex float::f: x1*y1
+DEAL:Copy:Complex float::g: x1/y1
+DEAL:Copy:Complex float::f(x,y): (60.00000000,0.000000000)
+DEAL:Copy:Complex float::g(x,y): (6.666666508,0.000000000)
+DEAL:Copy:Complex float::
+DEAL:Copy:Complex float::After copy
+DEAL:Copy:Complex float::f: x*y
+DEAL:Copy:Complex float::g: x/y
+DEAL:Copy:Complex float::f(x,y): (20.00000000,0.000000000)
+DEAL:Copy:Complex float::g(x,y): (5.000000000,0.000000000)
+DEAL:Copy:Complex float::
+DEAL:Copy:Complex float::OK
+DEAL:Copy:Complex double::f: x*y
+DEAL:Copy:Complex double::g: x/y
+DEAL:Copy:Complex double::f(x,y): (20.00000000,0.000000000)
+DEAL:Copy:Complex double::g(x,y): (5.000000000,0.000000000)
+DEAL:Copy:Complex double::
+DEAL:Copy:Complex double::Before copy
+DEAL:Copy:Complex double::f: x1*y1
+DEAL:Copy:Complex double::g: x1/y1
+DEAL:Copy:Complex double::f(x,y): (60.00000000,0.000000000)
+DEAL:Copy:Complex double::g(x,y): (6.666666667,0.000000000)
+DEAL:Copy:Complex double::
+DEAL:Copy:Complex double::After copy
+DEAL:Copy:Complex double::f: x*y
+DEAL:Copy:Complex double::g: x/y
+DEAL:Copy:Complex double::f(x,y): (20.00000000,0.000000000)
+DEAL:Copy:Complex double::g(x,y): (5.000000000,0.000000000)
+DEAL:Copy:Complex double::
+DEAL:Copy:Complex double::OK
+DEAL::OK
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 the BatchOptimizer::copy_from() function works as expected.
+//
+// Here we invoke the LLVM optimizer before symbolic evaluation takes place,
+// and we invoke no additional symbolic optimizations as well.
+
+#include "../tests.h"
+
+#include "sd_common_tests/batch_optimizer_08.h"
+
+int
+main()
+{
+ 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::llvm;
+ const enum SD::OptimizationFlags opt_flags =
+ SD::OptimizationFlags::optimize_default;
+
+ const int dim = 2;
+ run_tests<dim, opt_method, opt_flags>();
+}
--- /dev/null
+
+DEAL:Copy:Float::f: x*y
+DEAL:Copy:Float::g: x/y
+DEAL:Copy:Float::f(x,y): 20.00000000
+DEAL:Copy:Float::g(x,y): 5.000000000
+DEAL:Copy:Float::
+DEAL:Copy:Float::Before copy
+DEAL:Copy:Float::f: x1*y1
+DEAL:Copy:Float::g: x1/y1
+DEAL:Copy:Float::f(x,y): 60.00000000
+DEAL:Copy:Float::g(x,y): 6.666666508
+DEAL:Copy:Float::
+DEAL:Copy:Float::After copy
+DEAL:Copy:Float::f: x*y
+DEAL:Copy:Float::g: x/y
+DEAL:Copy:Float::f(x,y): 20.00000000
+DEAL:Copy:Float::g(x,y): 5.000000000
+DEAL:Copy:Float::
+DEAL:Copy:Float::OK
+DEAL:Copy:Double::f: x*y
+DEAL:Copy:Double::g: x/y
+DEAL:Copy:Double::f(x,y): 20.00000000
+DEAL:Copy:Double::g(x,y): 5.000000000
+DEAL:Copy:Double::
+DEAL:Copy:Double::Before copy
+DEAL:Copy:Double::f: x1*y1
+DEAL:Copy:Double::g: x1/y1
+DEAL:Copy:Double::f(x,y): 60.00000000
+DEAL:Copy:Double::g(x,y): 6.666666667
+DEAL:Copy:Double::
+DEAL:Copy:Double::After copy
+DEAL:Copy:Double::f: x*y
+DEAL:Copy:Double::g: x/y
+DEAL:Copy:Double::f(x,y): 20.00000000
+DEAL:Copy:Double::g(x,y): 5.000000000
+DEAL:Copy:Double::
+DEAL:Copy:Double::OK
+DEAL::OK
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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 the BatchOptimizer::copy_from() function works as expected.
+
+#include <deal.II/differentiation/sd.h>
+
+#include "../../tests.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_copy()
+{
+ namespace SD = Differentiation::SD;
+ using SD_number_t = SD::Expression;
+
+ // Define
+ const SD_number_t x("x"), y("y");
+ const SD_number_t f = x * y;
+ const SD_number_t g = x / y;
+
+ // Substitution map
+ SD::types::substitution_map sub_vals =
+ SD::make_substitution_map(std::make_pair(x, NumberType(10.0)),
+ std::make_pair(y, NumberType(2.0)));
+
+ // Optimize
+ SD::BatchOptimizer<NumberType> optimizer(opt_method, opt_flags);
+ optimizer.register_symbols(sub_vals); // Independent symbols
+ optimizer.register_functions(f, g);
+ optimizer.optimize();
+
+ // Substitute
+ optimizer.substitute(sub_vals);
+
+ // Evaluate
+ const NumberType val_f = optimizer.evaluate(f);
+ const NumberType val_g = optimizer.evaluate(g);
+ deallog << "f: " << f << std::endl
+ << "g: " << g << std::endl
+ << "f(x,y): " << val_f << std::endl
+ << "g(x,y): " << val_g << std::endl
+ << std::endl;
+
+ // Define something new
+ {
+ deallog << "Before copy" << std::endl;
+
+ const SD_number_t x1("x1"), y1("y1");
+ const SD_number_t f1 = x1 * y1;
+ const SD_number_t g1 = x1 / y1;
+
+ SD::add_to_substitution_map(sub_vals,
+ std::make_pair(x1, NumberType(20.0)),
+ std::make_pair(y1, NumberType(3.0)));
+
+ // Optimizer
+ SD::BatchOptimizer<NumberType> new_optimizer;
+ new_optimizer.register_symbols(sub_vals); // Independent symbols
+ new_optimizer.register_functions(f1, g1);
+ new_optimizer.optimize();
+
+ // Substitute
+ new_optimizer.substitute(sub_vals);
+
+ // Evaluate
+ const NumberType val_f1 = new_optimizer.evaluate(f1);
+ const NumberType val_g1 = new_optimizer.evaluate(g1);
+ deallog << "f: " << f1 << std::endl
+ << "g: " << g1 << std::endl
+ << "f(x,y): " << val_f1 << std::endl
+ << "g(x,y): " << val_g1 << std::endl
+ << std::endl;
+
+ // Perform copy
+ deallog << "After copy" << std::endl;
+ new_optimizer.copy_from(optimizer);
+ Assert(new_optimizer.optimized() == false,
+ ExcMessage("Expected new optimizer to be unoptimized."));
+
+ // Extract
+ const NumberType val_f = new_optimizer.extract(f, optimizer.evaluate());
+ const NumberType val_g = new_optimizer.extract(g, optimizer.evaluate());
+ deallog << "f: " << f << std::endl
+ << "g: " << g << std::endl
+ << "f(x,y): " << val_f << std::endl
+ << "g(x,y): " << val_g << std::endl
+ << std::endl;
+ }
+
+ 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("Copy");
+ {
+ // deallog.push("Integer");
+ // test_tensor<dim,int>(n_runs, timer);
+ // deallog.pop();
+
+ deallog.push("Float");
+ test_copy<dim, float, opt_method, opt_flags>();
+ deallog.pop();
+
+ deallog.push("Double");
+ test_copy<dim, double, opt_method, opt_flags>();
+ deallog.pop();
+
+ // The LLVM optimizer does not currently support complex numbers.
+ if (opt_method != SD::OptimizerType::llvm)
+ {
+ deallog.push("Complex float");
+ test_copy<dim, std::complex<float>, opt_method, opt_flags>();
+ deallog.pop();
+
+ deallog.push("Complex double");
+ test_copy<dim, std::complex<double>, opt_method, opt_flags>();
+ deallog.pop();
+ }
+ }
+ deallog.pop();
+
+ deallog << "OK" << std::endl;
+}