Reference documentation for deal.II version GIT relicensing-426-g7976cfd195 2024-04-18 21:10:01+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
Differentiation::SD::BatchOptimizer< ReturnType > Class Template Reference

#include <deal.II/differentiation/sd/symengine_optimizer.h>

Public Member Functions

 BatchOptimizer ()
 
 BatchOptimizer (const enum OptimizerType &optimization_method, const enum OptimizationFlags &optimization_flags=OptimizationFlags::optimize_all)
 
 BatchOptimizer (const BatchOptimizer &other)
 
 BatchOptimizer (BatchOptimizer &&) noexcept=default
 
 ~BatchOptimizer ()=default
 
void copy_from (const BatchOptimizer &other)
 
template<typename Stream >
void print (Stream &stream, const bool print_cse=false) const
 
template<class Archive >
void save (Archive &archive, const unsigned int version) const
 
template<class Archive >
void load (Archive &archive, const unsigned int version)
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
Independent variables
void register_symbols (const types::substitution_map &substitution_map)
 
void register_symbols (const SymEngine::map_basic_basic &substitution_map)
 
void register_symbols (const types::symbol_vector &symbols)
 
void register_symbols (const SymEngine::vec_basic &symbols)
 
types::symbol_vector get_independent_symbols () const
 
std::size_t n_independent_variables () const
 
Dependent variables
void register_function (const Expression &function)
 
template<int rank, int dim>
void register_function (const Tensor< rank, dim, Expression > &function_tensor)
 
template<int rank, int dim>
void register_function (const SymmetricTensor< rank, dim, Expression > &function_tensor)
 
void register_functions (const types::symbol_vector &functions)
 
void register_functions (const SymEngine::vec_basic &functions)
 
template<typename T >
void register_functions (const std::vector< T > &functions)
 
template<typename T , typename... Args>
void register_functions (const T &functions, const Args &...other_functions)
 
const types::symbol_vectorget_dependent_functions () const
 
std::size_t n_dependent_variables () const
 
Optimization
void set_optimization_method (const enum OptimizerType &optimization_method, const enum OptimizationFlags &optimization_flags=OptimizationFlags::optimize_all)
 
enum OptimizerType optimization_method () const
 
enum OptimizationFlags optimization_flags () const
 
bool use_symbolic_CSE () const
 
void optimize ()
 
bool optimized () const
 
Symbol substitution
void substitute (const types::substitution_map &substitution_map) const
 
void substitute (const SymEngine::map_basic_basic &substitution_map) const
 
void substitute (const types::symbol_vector &symbols, const std::vector< ReturnType > &values) const
 
void substitute (const SymEngine::vec_basic &symbols, const std::vector< ReturnType > &values) const
 
bool values_substituted () const
 
Evaluation / data extraction
const std::vector< ReturnType > & evaluate () const
 
ReturnType evaluate (const Expression &func) const
 
std::vector< ReturnType > evaluate (const std::vector< Expression > &funcs) const
 
template<int rank, int dim>
Tensor< rank, dim, ReturnType > evaluate (const Tensor< rank, dim, Expression > &funcs) const
 
template<int rank, int dim>
SymmetricTensor< rank, dim, ReturnType > evaluate (const SymmetricTensor< rank, dim, Expression > &funcs) const
 
ReturnType extract (const Expression &func, const std::vector< ReturnType > &cached_evaluation) const
 
std::vector< ReturnType > extract (const std::vector< Expression > &funcs, const std::vector< ReturnType > &cached_evaluation) const
 
template<int rank, int dim>
Tensor< rank, dim, ReturnType > extract (const Tensor< rank, dim, Expression > &funcs, const std::vector< ReturnType > &cached_evaluation) const
 
template<int rank, int dim>
SymmetricTensor< rank, dim, ReturnType > extract (const SymmetricTensor< rank, dim, Expression > &funcs, const std::vector< ReturnType > &cached_evaluation) const
 

Private Types

using map_dependent_expression_to_vector_entry_t = std::map< SD::Expression, std::size_t, SD::types::internal::ExpressionKeyLess >
 

Private Member Functions

bool is_valid_nonunique_dependent_variable (const SD::Expression &function) const
 
bool is_valid_nonunique_dependent_variable (const SymEngine::RCP< const SymEngine::Basic > &function) const
 
void register_scalar_function (const SD::Expression &function)
 
void register_vector_functions (const types::symbol_vector &functions)
 
void create_optimizer (std::unique_ptr< SymEngine::Visitor > &optimizer)
 
void substitute (const std::vector< ReturnType > &substitution_values) const
 

Private Attributes

enum OptimizerType method
 
enum OptimizationFlags flags
 
types::substitution_map independent_variables_symbols
 
types::symbol_vector dependent_variables_functions
 
std::vector< ReturnType > dependent_variables_output
 
map_dependent_expression_to_vector_entry_t map_dep_expr_vec_entry
 
std::unique_ptr< SymEngine::Visitor > optimizer
 
bool ready_for_value_extraction
 
bool has_been_serialized
 

Detailed Description

template<typename ReturnType>
class Differentiation::SD::BatchOptimizer< ReturnType >

A class that facilitates the optimization of symbol expressions.

This expression will be optimized by this class; that is to say that the code path taken to substitute the set of (independent) symbols into a collection of (dependent) symbolic functions will be optimized using a chosen approach.

This snippet of pseudo-code describes the general usage of this class:

// Define some independent variables
const Expression x("x");
const Expression y("y");
...
// Compute some symbolic expressions that are dependent on the
// independent variables. These could be, for example, scalar
// expressions or tensors of expressions.
const auto f = calculate_f(x, y, ...);
const auto g = calculate_g(x, y, ...);
...
// Now create a optimizer to evaluate the dependent functions.
// The numerical result will be of type double, and a "lambda" optimizer,
// which employs common subexpression elimination, will be used.
using ReturnType = double;
// Register symbols that represent independent variables...
optimizer.register_symbols(x, y, ...);
// ... and symbolic expressions that represent dependent functions.
optimizer.register_functions(f, g, ...);
// Now we determine an equivalent code path that will evaluate
// all of the dependent functions at once, but with less computational
// cost than when evaluating the symbolic expression directly.
optimizer.optimize(); // Note: This is an expensive call.
// Next we pass the optimizer the numeric values that we wish the
// independent variables to represent.
const auto substitution_map
= make_substitution_map({x, ...}, {y, ...}, ...);
// When making this next call, the call path used to (numerically)
// evaluate the dependent functions is quicker than dictionary
// substitution.
optimizer.substitute(substitution_map);
// Finally, we can get the numeric equivalent of the dependent functions
// from the optimizer.
const auto result_f = optimizer.evaluate(f);
const auto result_g = optimizer.evaluate(g);
std::unique_ptr< SymEngine::Visitor > optimizer
types::substitution_map make_substitution_map(const Expression &symbol, const Expression &value)

Since the call to optimize() may be quite costly, there are a few "best practices" that can be adopted in order to mitigate this cost as much as possible:

  1. Reuse a single instance of the class as much as possible. The most obvious way that this can be achieved would be to place an instance of this class in a centralized location where it can potentially be used by multiple calling functions and objects, if contextually possible.
  2. Another form of reuse would entail generalizing the dependent functions/expressions to be evaluated by the optimizer as much as possible. For example, material coefficients need not necessarily be hard-coded, and one generalized statement of a constitutive law could then be broadly used in other material subdomains governed by the same class of constitutive law, but with different constitutive parameters. The same principle applies if using symbolic expressions to describe boundary conditions, systems of linear equations, etc.
  3. When possible, consider using serialization to save and load the state of an optimizer that has already "optimized", i.e. it has been placed in a state where it is ready to evaluate expressions. With the exception of "lambda" optimization, all other forms of optimization permit checkpointing, meaning that the optimization could be done up front before executing the main body of code. It could also be used to duplicate an optimizer in an efficient manner, should multiple instances of the same optimizer be required.
Template Parameters
ReturnTypeThe number type that is to be returned after value substitution and evaluation. Floating point and complex numbers are currently supported.
Warning
This class is not thread-safe.
The LLVM optimizer does not yet support complex numbers. If this incompatible combination of ReturnType and optimization method are selected, then an error will be thrown at run time.

Definition at line 1434 of file symengine_optimizer.h.

Member Typedef Documentation

◆ map_dependent_expression_to_vector_entry_t

template<typename ReturnType >
using Differentiation::SD::BatchOptimizer< ReturnType >::map_dependent_expression_to_vector_entry_t = std::map<SD::Expression, std::size_t, SD::types::internal::ExpressionKeyLess>
private

A map type used to indicate which dependent variable is associated with which entry in the output vector.

Note
We use a custom comparator here because otherwise we can't use std::map::find; it is sensitive to the some data from the underlying SymEngine::Basic other than the value that it represents.

Definition at line 2064 of file symengine_optimizer.h.

Constructor & Destructor Documentation

◆ BatchOptimizer() [1/4]

template<typename ReturnType >
Differentiation::SD::BatchOptimizer< ReturnType >::BatchOptimizer ( )

Default constructor.

By default, dictionary substitution will be selected when this constructor is called. In order to select a specific optimization approach, a call to set_optimization_method() is necessary.

Definition at line 35 of file symengine_optimizer.cc.

◆ BatchOptimizer() [2/4]

template<typename ReturnType >
Differentiation::SD::BatchOptimizer< ReturnType >::BatchOptimizer ( const enum OptimizerType optimization_method,
const enum OptimizationFlags optimization_flags = OptimizationFlags::optimize_all 
)

Constructor.

Parameters
[in]optimization_methodThe optimization method that is to be employed.
[in]optimization_flagsThe optimization flags that indicate which expression manipulation mechanisms are to be employed.
Note
As the optimization method is fixed, a further call to set_optimization_method() is not necessary and will result in an error being thrown.
In the case that the optimization_method is not implemented for the required ReturnType, or the desired feature is not active, then an error will be thrown. Currently the LLVM optimization method is not compatible with complex numbers.

Definition at line 45 of file symengine_optimizer.cc.

◆ BatchOptimizer() [3/4]

template<typename ReturnType >
Differentiation::SD::BatchOptimizer< ReturnType >::BatchOptimizer ( const BatchOptimizer< ReturnType > &  other)

Copy constructor.

Note
The optimized data and results from previous substitutions executed by the other optimizer instance are not copied over. It is therefore necessary to re-optimize the data stored in this class, and it is possible to do so with a different optimization scheme.

Definition at line 56 of file symengine_optimizer.cc.

◆ BatchOptimizer() [4/4]

template<typename ReturnType >
Differentiation::SD::BatchOptimizer< ReturnType >::BatchOptimizer ( BatchOptimizer< ReturnType > &&  )
defaultnoexcept

Move constructor.

◆ ~BatchOptimizer()

template<typename ReturnType >
Differentiation::SD::BatchOptimizer< ReturnType >::~BatchOptimizer ( )
default

Destructor.

Member Function Documentation

◆ copy_from()

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::copy_from ( const BatchOptimizer< ReturnType > &  other)

Duplicate the data stored in an other BatchOptimizer instance.

Note
The optimized data and results from previous substitutions executed by the 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 this optimizer instance if those results are stored elsewhere.

Definition at line 72 of file symengine_optimizer.cc.

◆ print()

template<typename ReturnType >
template<typename Stream >
void Differentiation::SD::BatchOptimizer< ReturnType >::print ( Stream &  stream,
const bool  print_cse = false 
) const

Print some information on state of the internal data structures stored in the class.

Template Parameters
StreamThe type for the output stream.
Parameters
streamThe output stream to print to.
print_cseA flag to indicate whether or not all common subexpressions should be printed to the stream.

◆ save()

template<typename ReturnType >
template<class Archive >
void Differentiation::SD::BatchOptimizer< ReturnType >::save ( Archive &  archive,
const unsigned int  version 
) const

Write the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

This effectively saves the value stored into the archive with the given version number into this object.

◆ load()

template<typename ReturnType >
template<class Archive >
void Differentiation::SD::BatchOptimizer< ReturnType >::load ( Archive &  archive,
const unsigned int  version 
)

Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

This effectively loads the value stored out of the archive with the given version number into this object. In doing so, the previous contents of this object are thrown away.

Note
When deserializing a symbolic expression, it is imperative that you first create or deserialize all of the symbolic variables used in the serialized expression.

◆ serialize()

template<typename ReturnType >
template<class Archive >
void Differentiation::SD::BatchOptimizer< ReturnType >::serialize ( Archive &  archive,
const unsigned int  version 
)

Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

This effectively saves or loads the value stored into/out of the archive with the given version number into this object. If deserializing data, then the previous contents of this object are thrown away.

Note
When deserializing a batch optimizer, it is imperative that you first create or deserialize all of the symbolic variables and symbolic functions used in the optimizer.
Complete serialization is not possible when the "lambda" optimization method is invoked. Although the registered symbols and dependent function expressions are stored, the optimization is itself not stored. It might, therefore, take some time for the deserialization when "lambda" optimization is used as the optimization step will be (automatically) performed once more.

◆ register_symbols() [1/4]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_symbols ( const types::substitution_map substitution_map)

Register a collection of symbols that represents an independent variable. These symbols are stored as the key to the substitution_map.

Definition at line 165 of file symengine_optimizer.cc.

◆ register_symbols() [2/4]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_symbols ( const SymEngine::map_basic_basic &  substitution_map)

Register a collection of symbols that represents an independent variable. These symbols are stored as the key to the substitution_map.

Definition at line 192 of file symengine_optimizer.cc.

◆ register_symbols() [3/4]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_symbols ( const types::symbol_vector symbols)

Register a collection of symbols that represents independent variables.

Warning
When using this function is no mechanism to check that the ordering of the later used substitution_values vector or map matches the internal ordering of the registered symbols. This function is therefore typically used in conjunction with the substitute() function that takes in a vector of values. With this pair of functions to the class interface, the management of symbol ordering is maintained by the user.

Definition at line 203 of file symengine_optimizer.cc.

◆ register_symbols() [4/4]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_symbols ( const SymEngine::vec_basic &  symbols)

Register a collection of symbols that represents independent variables.

Warning
When using this function is no mechanism to check that the ordering of the later used substitution_values vector or map matches the internal ordering of the registered symbols. This function is therefore typically used in conjunction with the substitute() function that takes in a vector of values. With this pair of functions to the class interface, the management of symbol ordering is maintained by the user.

Definition at line 224 of file symengine_optimizer.cc.

◆ get_independent_symbols()

template<typename ReturnType >
SD::types::symbol_vector Differentiation::SD::BatchOptimizer< ReturnType >::get_independent_symbols ( ) const

Return a vector of symbols that have been registered as independent variables.

Definition at line 235 of file symengine_optimizer.cc.

◆ n_independent_variables()

template<typename ReturnType >
std::size_t Differentiation::SD::BatchOptimizer< ReturnType >::n_independent_variables ( ) const

The number of independent variables that this optimizer will recognize. This is equal to the number of unique symbols passed to this class instance through the register_symbols() function.

Definition at line 244 of file symengine_optimizer.cc.

◆ register_function() [1/3]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_function ( const Expression function)

Register a scalar symbolic expression that represents a dependent variable.

Definition at line 253 of file symengine_optimizer.cc.

◆ register_function() [2/3]

template<typename ReturnType >
template<int rank, int dim>
void Differentiation::SD::BatchOptimizer< ReturnType >::register_function ( const Tensor< rank, dim, Expression > &  function_tensor)

Register a tensor of symbolic expressions that represents a dependent variable.

◆ register_function() [3/3]

template<typename ReturnType >
template<int rank, int dim>
void Differentiation::SD::BatchOptimizer< ReturnType >::register_function ( const SymmetricTensor< rank, dim, Expression > &  function_tensor)

Register a symmetric tensor of symbolic expressions that represents a dependent variable.

◆ register_functions() [1/4]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_functions ( const types::symbol_vector functions)

Register a collection of symbolic expressions that represent dependent variables.

Definition at line 266 of file symengine_optimizer.cc.

◆ register_functions() [2/4]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_functions ( const SymEngine::vec_basic &  functions)

Register a collection of symbolic expressions that represent multiple dependent variables.

Definition at line 280 of file symengine_optimizer.cc.

◆ register_functions() [3/4]

template<typename ReturnType >
template<typename T >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_functions ( const std::vector< T > &  functions)

Register a collection of symbolic expressions that represent multiple dependent variables.

Template Parameters
TA compatible type that may be used to represent a single dependent variable. This includes scalar Expressions, Tensors of Expressions and SymmetricTensors of Expressions.
Parameters
functionsA vector of symbolic dependent variables.

◆ register_functions() [4/4]

template<typename ReturnType >
template<typename T , typename... Args>
void Differentiation::SD::BatchOptimizer< ReturnType >::register_functions ( const T &  functions,
const Args &...  other_functions 
)

Register a collection of symbolic expressions that represent dependent variables.

Template Parameters
TA compatible type that may be used to represent a single dependent variable. This includes scalar Expressions, Tensors of Expressions, SymmetricTensors of Expressions, and std::vectors of Expressions.
ArgsA variadic template that represents a collection of any compatible symbolic dependent variable types.
Parameters
functionsOne or more symbolic dependent variables.
other_functionsAn arbitrary collection of symbolic dependent variables.

◆ get_dependent_functions()

template<typename ReturnType >
const SD::types::symbol_vector & Differentiation::SD::BatchOptimizer< ReturnType >::get_dependent_functions ( ) const

Return a vector of expressions that have been registered as dependent variables.

Definition at line 291 of file symengine_optimizer.cc.

◆ n_dependent_variables()

template<typename ReturnType >
std::size_t Differentiation::SD::BatchOptimizer< ReturnType >::n_dependent_variables ( ) const

The number of dependent symbolic expressions that this optimizer will optimize. This is equal to the number of unique symbolic functions / expressions passed to this class instance through the register_functions() method.

Definition at line 300 of file symengine_optimizer.cc.

◆ set_optimization_method()

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::set_optimization_method ( const enum OptimizerType optimization_method,
const enum OptimizationFlags optimization_flags = OptimizationFlags::optimize_all 
)

Select the optimization_method for the batch optimizer to employ, in conjunction with the optimization_flags.

It is required that the this class instance is not yet optimized, i.e. that the optimize() method has not yet been called.

Note
In the case that the method is not implemented for the required ReturnType, or the desired feature is not active, then a safe default will be selected.

Definition at line 89 of file symengine_optimizer.cc.

◆ optimization_method()

template<typename ReturnType >
enum OptimizerType Differentiation::SD::BatchOptimizer< ReturnType >::optimization_method ( ) const

Return the optimization method that has been selected for use.

Definition at line 112 of file symengine_optimizer.cc.

◆ optimization_flags()

template<typename ReturnType >
enum OptimizationFlags Differentiation::SD::BatchOptimizer< ReturnType >::optimization_flags ( ) const

Return the optimization flags that have been selected for use.

Definition at line 121 of file symengine_optimizer.cc.

◆ use_symbolic_CSE()

template<typename ReturnType >
bool Differentiation::SD::BatchOptimizer< ReturnType >::use_symbolic_CSE ( ) const

State whether the internal selection of optimization methods and flags will render an optimizer that uses common subexpression elimination (CSE).

Definition at line 130 of file symengine_optimizer.cc.

◆ optimize()

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::optimize ( )

Perform the optimization of all registered dependent functions using the registered symbols.

Note
This function, which should only be called once per instance of this class, finalizes the set of accepted independent symbols and dependent functions that are recognized and used by the optimizer.
This may be a time-consuming process, but if the class instance is retained throughout the course of a simulation (and both the independent and dependent variables that are associated with the class instance remain unchanged) then it need only be performed once. Serialization also offers the opportunity to reuse the already computed optimized evaluation call path.

Definition at line 317 of file symengine_optimizer.cc.

◆ optimized()

template<typename ReturnType >
bool Differentiation::SD::BatchOptimizer< ReturnType >::optimized ( ) const

Returns a flag which indicates whether the optimize() function has been called and the class is finalized.

Definition at line 139 of file symengine_optimizer.cc.

◆ substitute() [1/5]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::substitute ( const types::substitution_map substitution_map) const

Perform batch substitution of all of the registered symbols into the registered functions. The result is cached and can be extracted by calls to evaluate().

Note
Calling substitute() again with a new substitution_map overwrites any previously computed results.

Definition at line 394 of file symengine_optimizer.cc.

◆ substitute() [2/5]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::substitute ( const SymEngine::map_basic_basic &  substitution_map) const

Perform batch substitution of all of the registered symbols into the registered functions. The result is cached and can be extracted by calls to evaluate().

Note
Calling substitute() again with a new substitution_map overwrites any previously computed results.

Definition at line 433 of file symengine_optimizer.cc.

◆ substitute() [3/5]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::substitute ( const types::symbol_vector symbols,
const std::vector< ReturnType > &  values 
) const

Perform batch substitution of all of the registered symbols into the registered functions. The result is cached and can be extracted by calls to evaluate(). It is expected that there is a 1-1 correspondence between each of the symbols and values.

Note
Calling substitute() again with a new set of values overwrites any previously computed results.

Definition at line 444 of file symengine_optimizer.cc.

◆ substitute() [4/5]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::substitute ( const SymEngine::vec_basic &  symbols,
const std::vector< ReturnType > &  values 
) const

Perform batch substitution of all of the registered symbols into the registered functions. The result is cached and can be extracted by calls to evaluate(). It is expected that there is a 1-1 correspondence between each of the symbols and values.

Note
Calling substitute() again with a new set of values overwrites any previously computed results.

Definition at line 458 of file symengine_optimizer.cc.

◆ values_substituted()

template<typename ReturnType >
bool Differentiation::SD::BatchOptimizer< ReturnType >::values_substituted ( ) const

Returns a flag to indicate whether the substitute() function has been called and if there are meaningful values that will be returned upon evaluation.

Definition at line 156 of file symengine_optimizer.cc.

◆ evaluate() [1/5]

template<typename ReturnType >
const std::vector< ReturnType > & Differentiation::SD::BatchOptimizer< ReturnType >::evaluate ( ) const

Returns the result of a value substitution into the optimized counterpart of all dependent functions. This function fetches all of those cached values. If these results are stored elsewhere, and you want to use this optimizer to extract components of the cached results then you can make use of the extract() functions listed below.

These values were computed by substituting a substitution_values map during substitute() call.

Note
In contrast to the other variants of this function, the order of the returned entries is an internal implementation detail, and cannot be guaranteed under all conditions. The entries in the returned vector are, in general, identical to the order in which the dependent expressions are originally registered. However, when registering tensors and symmetric tensors of expressions, these are "unrolled" and their components are individually registered. If it is necessary to control the order in which results appear in the returned vector, then the register_function() method that takes a Tensor or a SymmetricTensor as an argument should be avoided. Instead, the individual entries of these types of data should be registered one by one.

Definition at line 528 of file symengine_optimizer.cc.

◆ evaluate() [2/5]

template<typename ReturnType >
ReturnType Differentiation::SD::BatchOptimizer< ReturnType >::evaluate ( const Expression func) const

Returns the result of a value substitution into the optimized counterpart of func. This function fetches that one cached value.

This value was computed by substituting a substitution_values map during substitute() call.

Definition at line 636 of file symengine_optimizer.cc.

◆ evaluate() [3/5]

template<typename ReturnType >
std::vector< ReturnType > Differentiation::SD::BatchOptimizer< ReturnType >::evaluate ( const std::vector< Expression > &  funcs) const

Returns the result of a value substitution into the optimized counterpart of funcs. This function fetches that subset of cached values.

This value was computed by substituting a substitution_values map during substitute() call.

Definition at line 668 of file symengine_optimizer.cc.

◆ evaluate() [4/5]

template<typename ReturnType >
template<int rank, int dim>
Tensor< rank, dim, ReturnType > Differentiation::SD::BatchOptimizer< ReturnType >::evaluate ( const Tensor< rank, dim, Expression > &  funcs) const

Returns the result of a tensor value substitution into the optimized counterpart of funcs. This function fetches those cached tensor components.

This value was computed by substituting a substitution_values map during substitute() call.

◆ evaluate() [5/5]

template<typename ReturnType >
template<int rank, int dim>
SymmetricTensor< rank, dim, ReturnType > Differentiation::SD::BatchOptimizer< ReturnType >::evaluate ( const SymmetricTensor< rank, dim, Expression > &  funcs) const

Returns the result of a tensor value substitution into the optimized counterpart of funcs. This function fetches those cached symmetric tensor components.

This value was computed by substituting a substitution_values map during substitute() call.

◆ extract() [1/4]

template<typename ReturnType >
ReturnType Differentiation::SD::BatchOptimizer< ReturnType >::extract ( const Expression func,
const std::vector< ReturnType > &  cached_evaluation 
) const

Returns the result of a value substitution into the optimized counterpart of func. This function fetches that one value stored in the cached_evaluation vector. This vector is, most typically, first attained by a call to the evaluate() variant that takes no arguments.

Definition at line 543 of file symengine_optimizer.cc.

◆ extract() [2/4]

template<typename ReturnType >
std::vector< ReturnType > Differentiation::SD::BatchOptimizer< ReturnType >::extract ( const std::vector< Expression > &  funcs,
const std::vector< ReturnType > &  cached_evaluation 
) const

Returns the result of a value substitution into the optimized counterpart of funcs. This function fetches that subset of values stored in the cached_evaluation vector. This vector is, most typically, first attained by a call to the evaluate() variant that takes no arguments.

Definition at line 651 of file symengine_optimizer.cc.

◆ extract() [3/4]

template<typename ReturnType >
template<int rank, int dim>
Tensor< rank, dim, ReturnType > Differentiation::SD::BatchOptimizer< ReturnType >::extract ( const Tensor< rank, dim, Expression > &  funcs,
const std::vector< ReturnType > &  cached_evaluation 
) const

Returns the result of a tensor value substitution into the optimized counterpart of funcs. This function fetches those tensor components stored in the cached_evaluation vector. This vector is, most typically, first attained by a call to the evaluate() variant that takes no arguments.

◆ extract() [4/4]

template<typename ReturnType >
template<int rank, int dim>
SymmetricTensor< rank, dim, ReturnType > Differentiation::SD::BatchOptimizer< ReturnType >::extract ( const SymmetricTensor< rank, dim, Expression > &  funcs,
const std::vector< ReturnType > &  cached_evaluation 
) const

Returns the result of a tensor value substitution into the optimized counterpart of funcs. This function fetches those symmetric tensor components stored in the cached_evaluation vector. This vector is, most typically, first attained by a call to the evaluate() variant that takes no arguments.

◆ is_valid_nonunique_dependent_variable() [1/2]

template<typename ReturnType >
bool Differentiation::SD::BatchOptimizer< ReturnType >::is_valid_nonunique_dependent_variable ( const SD::Expression function) const
private

A check to see if a function is exactly equal to one of the logical results of a differentiation operation.

Definition at line 683 of file symengine_optimizer.cc.

◆ is_valid_nonunique_dependent_variable() [2/2]

template<typename ReturnType >
bool Differentiation::SD::BatchOptimizer< ReturnType >::is_valid_nonunique_dependent_variable ( const SymEngine::RCP< const SymEngine::Basic > &  function) const
private

A check to see if a function is exactly equal to one of the logical results of a differentiation operation.

Definition at line 693 of file symengine_optimizer.cc.

◆ register_scalar_function()

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_scalar_function ( const SD::Expression function)
private

Register a single symbol that represents a dependent variable.

Definition at line 726 of file symengine_optimizer.cc.

◆ register_vector_functions()

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::register_vector_functions ( const types::symbol_vector functions)
private

Register a collection of symbols that represent dependent variables.

Definition at line 754 of file symengine_optimizer.cc.

◆ create_optimizer()

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::create_optimizer ( std::unique_ptr< SymEngine::Visitor > &  optimizer)
private

Create an instance of the selected optimizer.

Definition at line 788 of file symengine_optimizer.cc.

◆ substitute() [5/5]

template<typename ReturnType >
void Differentiation::SD::BatchOptimizer< ReturnType >::substitute ( const std::vector< ReturnType > &  substitution_values) const
private

Perform batch substitution of all of the registered symbols into the registered functions. The result is cached and can be extracted by calls to evaluate().

Note
Calling substitute() again with a new set of substitution_values overwrites any previously computed results.
Warning
When using this function there is no mechanism to check that the ordering of the substitution_values vector matches the internal ordering of the registered symbols. This function is therefore typically used in conjunction with the register_symbols() function that takes in a vector of symbols. With this pair of functions to the class interface, the management of symbol ordering is maintained by the user.

Definition at line 471 of file symengine_optimizer.cc.

Member Data Documentation

◆ method

template<typename ReturnType >
enum OptimizerType Differentiation::SD::BatchOptimizer< ReturnType >::method
private

The optimization methods that is to be employed.

Definition at line 1998 of file symengine_optimizer.h.

◆ flags

template<typename ReturnType >
enum OptimizationFlags Differentiation::SD::BatchOptimizer< ReturnType >::flags
private

The optimization flags that indicate which expression manipulation mechanisms are to be employed.

Definition at line 2004 of file symengine_optimizer.h.

◆ independent_variables_symbols

template<typename ReturnType >
types::substitution_map Differentiation::SD::BatchOptimizer< ReturnType >::independent_variables_symbols
private

A map that represents the symbols that form the set of independent variables upon which optimized symbolic expressions are to be based.

Note
As the ordering of the input symbols is fixed at the time at which optimization is performed, we store all of the entries in a map to ensure that both we and the user never mistakenly swap the order of two or more symbols during evaluation.

Definition at line 2015 of file symengine_optimizer.h.

◆ dependent_variables_functions

template<typename ReturnType >
types::symbol_vector Differentiation::SD::BatchOptimizer< ReturnType >::dependent_variables_functions
private

A set of symbolic expressions to be optimized. It is required that the symbols on which these dependent functions are based are registered in the independent_variables_symbols map.

Definition at line 2022 of file symengine_optimizer.h.

◆ dependent_variables_output

template<typename ReturnType >
std::vector<ReturnType> Differentiation::SD::BatchOptimizer< ReturnType >::dependent_variables_output
mutableprivate

The output of substituting symbolic values with floating point values through the use of the optimizer.

It is necessary to use this intermediate storage mechanism to store the result of a substitution sweep as some optimizers work on all symbolic expressions in a single batch. In this way they can employ methods such as common subexpression elimination to minimize the number of terms evaluated across all symbolic functions.

This variable is marked as mutable. This facilitates the substitution functionality being used in logically constant get_* functions.

Definition at line 2054 of file symengine_optimizer.h.

◆ map_dep_expr_vec_entry

template<typename ReturnType >
map_dependent_expression_to_vector_entry_t Differentiation::SD::BatchOptimizer< ReturnType >::map_dep_expr_vec_entry
mutableprivate

A map indicating which dependent variable is associated with which entry in the output vector.

Definition at line 2073 of file symengine_optimizer.h.

◆ optimizer

template<typename ReturnType >
std::unique_ptr<SymEngine::Visitor> Differentiation::SD::BatchOptimizer< ReturnType >::optimizer
mutableprivate

A pointer to an instance of an optimizer that will be used to reformulate the substitution of symbolic expressions in a manner that is more efficient than plain dictionary-based approach.

Definition at line 2081 of file symengine_optimizer.h.

◆ ready_for_value_extraction

template<typename ReturnType >
bool Differentiation::SD::BatchOptimizer< ReturnType >::ready_for_value_extraction
mutableprivate

A flag to record whether or not substitution has taken place and values can now be extracted.

Note
This variable is marked as mutable. This facilitates the substitution functionality being used in logically constant get_* functions.

Definition at line 2091 of file symengine_optimizer.h.

◆ has_been_serialized

template<typename ReturnType >
bool Differentiation::SD::BatchOptimizer< ReturnType >::has_been_serialized
mutableprivate

A flag to record whether or not this class instance has been serialized in the past.

Definition at line 2097 of file symengine_optimizer.h.


The documentation for this class was generated from the following files: