TMB Documentation  v1.9.11
Classes | Enumerations | Functions
TMBad Namespace Reference

Automatic differentiation library designed for TMB. More...

Classes

struct  ad_plain_index
 Construct ad_plain from index. More...
 
struct  adaptive
 Enable weak comparison operators of an ad type. More...
 
struct  ADFun
 Automatic differentiation function object. More...
 
struct  ADFun_packed
 Container of ADFun object with packed input and output. More...
 
struct  Args
 Argument class to handle array access of operator inputs and outputs. More...
 
struct  AtomOp
 Generic checkpoint operator. More...
 
struct  autopar
 Split a computational graph using a simple heuristic. More...
 
struct  clique
 Utilility class for sequential_reduction. More...
 
struct  control
 User control parameters for R's integrate. More...
 
struct  Decomp2
 Decomposition of computational graph. More...
 
struct  Decomp3
 Decomposition of computational graph. More...
 
struct  forbid_remap
 Forbid remappings if not consecutive. More...
 
struct  ForwardArgs
 Access input/output values during a forward pass. Write access granted for the output value only. More...
 
struct  global
 Struct defining the main AD context. More...
 
struct  graph
 Operator graph in compressed row storage. More...
 
struct  IndirectAccessor
 Provide inplace read access to value or derivative arrays. More...
 
struct  Integral
 Interface to R's adaptive integrate routine. More...
 
struct  intervals
 Union of closed intervals. More...
 
struct  multivariate_index
 Utilility class for sequential_reduction. More...
 
struct  mvIntegral
 Multivariate integral class. More...
 
struct  old_state
 Read the current tape's state and restore it on request. More...
 
struct  omp_shared_ptr
 Manage shared operator data across multiple threads. More...
 
struct  op_info
 Bitwise collection of selected operator flags. More...
 
struct  PackOp
 Pack (PackOp) or unpack (UnpkOp) n consecutive values on the tape. More...
 
struct  PackWrap
 Transform a functor to have packed input/output. More...
 
struct  ParalOp
 Parallel operator. More...
 
struct  ParametersChanged
 Default tester for retaping_derivative_table. More...
 
struct  period
 Representation of a period in a sequence. More...
 
struct  periodic
 Period analyzer. More...
 
struct  retaping_derivative_table
 Adaptive derivative table used by AtomOp More...
 
struct  ReverseArgs
 Access input/output values and derivatives during a reverse pass. Write access granted for the input derivative only. More...
 
struct  segment_ref
 Provide read/write access to an array segment. More...
 
struct  SegmentRef
 Representation of a specific contiguous set of values on a specific tape. More...
 
struct  sequential_reduction
 Sequential reduction algorithm. More...
 
struct  SpJacFun_config
 Configuration parameters for SpJacFun() More...
 
struct  sr_grid
 Utilility class for sequential_reduction. More...
 
struct  standard_derivative_table
 Fixed derivative table used by AtomOp More...
 
struct  StdWrap
 Interoperability with other vector classes. More...
 
struct  UnpkOp
 Pack (PackOp) or unpack (UnpkOp) n consecutive values on the tape. More...
 

Enumerations

enum  ArrayAccess
 Define segment_ref array to access inside ForwardArgs or ReverseArgs
 

Functions

global accumulation_tree_split (global glob, bool sum_=false)
 Split a computational graph by it's accumulation tree. More...
 
template<class Functor , class Test = ParametersChanged>
ADFun ADFun_retaping (Functor &F, const std::vector< ad_aug > &x, Test test=Test())
 Construct ADFun that automatically retapes. More...
 
template<class Functor , class Test >
ADFun_packed ADFun_retaping (Functor &F, const std::vector< ad_segment > &x, Test test)
 Construct ADFun that automatically retapes. More...
 
void aggregate (global &glob, int sign=1)
 Reduce a multivariate output function by summing its range components. More...
 
bool all_allow_remap (const global &glob)
 Test if all operators in the stack allow input remapping. More...
 
template<class Matrix >
global::ad_segment contiguousBlock (const Matrix &x)
 Request a contiguous block on the tape. More...
 
std::vector< Index > find_op_by_name (global &glob, const char *name)
 Find nodes by name.
 
template<class V >
void forceContiguous (V &x)
 Make contiguous ad vector. More...
 
std::vector< Index > get_accumulation_tree (global &glob, bool boundary=false)
 Get node indices of the accumulation tree or its boundary. More...
 
globalget_glob ()
 Get pointer to current global AD context (or NULL if no context is active).
 
template<class V >
getContiguous (const V &x)
 Get contiguous (deep) copy of this vector. More...
 
template<class Integrand >
Integrand::Scalar integrate (Integrand f, typename Integrand::Scalar a=-INFINITY, typename Integrand::Scalar b=INFINITY, control c=control())
 Integrate function over finite or infinite interval. More...
 
template<class T >
std::vector< T > invperm (const std::vector< T > &perm)
 Inverse permutation.
 
template<class V >
bool isContiguous (V &x)
 Is this ad vector available as a contiguous block on the tape? More...
 
template<class T >
std::vector< bool > lmatch (const std::vector< T > &x, const std::vector< T > &y)
 Match x vector in y vector. More...
 
template<class T , class I >
void make_space_inplace (std::vector< T > &x, std::vector< I > &i, T space=T(0))
 Make space for new elements in an existing vector. More...
 
template<class T >
std::vector< size_t > match (const std::vector< T > &x, const std::vector< T > &y)
 Match x vector in y vector.
 
vmatrix matmul (const vmatrix &x, const vmatrix &y)
 Multiply two matrices of ad variables.
 
dmatrix matmul (const dmatrix &x, const dmatrix &y)
 Multiply two matrices of scalar types.
 
template<bool XT, bool YT, bool ZT, bool UP>
void matmul (Map< const dmatrix > x, Map< const dmatrix > y, Map< dmatrix > z)
 
template<class Integrand >
mvIntegral0< Integrand > mvIntegrate (Integrand &f, control c=control())
 Multivariate integration. More...
 
bool operator< (const ad_aug &x, const ad_aug &y)
 Allow / disallow comparison operators for ad_aug More...
 
template<class T >
std::vector< size_t > order (std::vector< T > x)
 Get permutation that sorts a vector. More...
 
ad_segment pack (const ad_segment &x)
 Pack consecutive values on the tape.
 
size_t prod_int (const std::vector< size_t > &x)
 Integer product function.
 
std::vector< Index > remap_identical_sub_expressions (global &glob, std::vector< Index > inv_remap)
 Remap identical sub-expressions. More...
 
void reorder_depth_first (global &glob)
 Depth-first reordering of computational graph. More...
 
void reorder_graph (global &glob, std::vector< Index > inv_idx)
 Reorder computational graph such that selected independent variables come last. More...
 
void reorder_sub_expressions (global &glob)
 Re-order computational graph to make it more compressible. More...
 
void reorder_temporaries (global &glob)
 Re-order computational graph to make it more compressible. More...
 
std::vector< bool > reverse_boundary (global &glob, const std::vector< bool > &vars)
 Find boundary without using the graph.
 
template<class T >
void sort_inplace (std::vector< T > &x)
 Utility: sort inplace.
 
template<class T >
void sort_unique_inplace (std::vector< T > &x)
 Utility: sort unique inplace.
 
std::vector< periodsplit_period (global *glob, period p, size_t max_period_size)
 Helper. More...
 
template<class T >
std::vector< T > subset (const std::vector< T > &x, const std::vector< bool > &y)
 Vector subset by boolean mask.
 
template<class T , class I >
std::vector< T > subset (const std::vector< T > &x, const std::vector< I > &ind)
 Vector subset by index vector.
 
std::vector< Index > substitute (global &glob, const std::vector< Index > &seq, bool inv_tags=true, bool dep_tags=true)
 substitute node index sequence by independent variables More...
 
std::vector< Index > substitute (global &glob, const char *name, bool inv_tags=true, bool dep_tags=true)
 substitute by op_name
 
template<class T >
ad_segment unpack (const std::vector< T > &x, Index j)
 Unpack consecutive values on the tape.
 
ad_segment unpack (const ad_segment &x)
 Unpack consecutive values on the tape.
 
template<class T >
double value (T x)
 Namespace with utility functions for adaptive numerical integration. More...
 
template<class I >
std::vector< I > which (const std::vector< bool > &x)
 Convert logical vector to index vector.
 

Detailed Description

Automatic differentiation library designed for TMB.

See TMBad::ADFun.

Function Documentation

§ accumulation_tree_split()

global TMBad::accumulation_tree_split ( global  glob,
bool  sum_ = false 
)

Split a computational graph by it's accumulation tree.

This routine transforms a function object with univariate output into a function object with multiple outputs. Each new output represents a term of the original function, i.e. by summing all the terms one can recover the original function.

Note
FIXME: Handle the case Range() > 1. Select which component ?

Definition at line 3613 of file TMBad.cpp.

Referenced by TMBad::ADFun<>::marginal_gk(), TMBad::ADFun<>::marginal_sr(), and order().

§ ADFun_retaping() [1/2]

template<class Functor , class Test = ParametersChanged>
ADFun TMBad::ADFun_retaping ( Functor &  F,
const std::vector< ad_aug > &  x,
Test  test = Test() 
)

Construct ADFun that automatically retapes.

By default, retaping takes place whenever the parameter vector changes compared to previous evaluation. However, this can be controlled in more detail by passing a custom tester object - see the prototype ParametersChanged.

Parameters
FFunctor to be taped
xEvaluation point
Template Parameters
TestClass to test if parameters have changed

Definition at line 1165 of file TMBad.hpp.

Referenced by sparse_matrix_exponential::expm_series< T >::operator()().

§ ADFun_retaping() [2/2]

template<class Functor , class Test >
ADFun_packed TMBad::ADFun_retaping ( Functor &  F,
const std::vector< ad_segment > &  x,
Test  test 
)

Construct ADFun that automatically retapes.

By default, retaping takes place whenever the parameter vector changes compared to previous evaluation. However, this can be controlled in more detail by passing a custom tester object - see the prototype ParametersChanged.

Parameters
FFunctor to be taped
xEvaluation point
Template Parameters
TestClass to test if parameters have changed

The resulting object is ADFun_packed, i.e. it has packed inputs and outputs. Such packed I/O can compactly represent e.g. matrices, vectors or other large objects with a consequtive memory layout.

Definition at line 1194 of file TMBad.hpp.

§ aggregate()

void TMBad::aggregate ( global glob,
int  sign = 1 
)

Reduce a multivariate output function by summing its range components.

On input the tape represents a function f:R^m->R^n. On output the tape represents a function f:R^m->R.

Parameters
globTape which is modified by this function

Definition at line 3653 of file TMBad.cpp.

Referenced by TMBad::autopar::extract(), TMBad::ADFun<>::marginal_gk(), TMBad::ADFun<>::marginal_sr(), and order().

§ all_allow_remap()

bool TMBad::all_allow_remap ( const global glob)

Test if all operators in the stack allow input remapping.

I.e. test if op_info.allow_remap = true for all operators.

Definition at line 4291 of file TMBad.cpp.

Referenced by reorder_graph().

§ contiguousBlock()

template<class Matrix >
global::ad_segment TMBad::contiguousBlock ( const Matrix &  x)

Request a contiguous block on the tape.

  1. Check if x already is on the tape and satisfies the storage requirement.
  2. If no invoke a deep copy of x to the tape and update x with the new tape addresses.
Returns
A reference to x as a contiguous block on the tape.
Note
The update step is critical as it ensures that a given matrix can be used several times without invoking a deep copy more than once.

Definition at line 23 of file ad_blas.hpp.

§ forceContiguous()

template<class V >
void TMBad::forceContiguous ( V &  x)

Make contiguous ad vector.

Template type 'V::value_type' can be

  • ad_plain
  • ad_aug
  • ad_adapt

Definition at line 3083 of file global.hpp.

Referenced by TMBad::sequential_reduction::tabulate().

§ get_accumulation_tree()

std::vector< Index > TMBad::get_accumulation_tree ( global glob,
bool  boundary = false 
)

Get node indices of the accumulation tree or its boundary.

The tree is defined as the maximal topologically closed sub graph T, containing only linear operators. In other words, the nodes in the tree are not allowed to produce variables that are inputs to any non linear operators.

Definition at line 3550 of file TMBad.cpp.

Referenced by accumulation_tree_split(), and order().

§ getContiguous()

template<class V >
V TMBad::getContiguous ( const V &  x)

Get contiguous (deep) copy of this vector.

Template type 'V::value_type' can be

  • ad_plain
  • ad_aug
  • ad_adapt

Definition at line 3071 of file global.hpp.

Referenced by forceContiguous().

§ integrate()

template<class Integrand >
Integrand::Scalar TMBad::integrate ( Integrand  f,
typename Integrand::Scalar  a = -INFINITY,
typename Integrand::Scalar  b = INFINITY,
control  c = control() 
)

Integrate function over finite or infinite interval.

Parameters
fUnivariate integrand (functor)
aLower integration limit. Default is negative infinity.
aUpper integration limit. Default is positive infinity.
cOptional control parameters.

Example:

template<class Float>
struct Gauss_t {
typedef Float Scalar;
Float a; // Parameter
// Evaluate integrand
Float operator(Float x) () {
Float ans = exp(- a*x*x);
return ans;
}
// Integrate wrt x
Float my_integrate() {
Float ans = integrate(*this);
return ans;
}
};

Definition at line 1111 of file TMBad/integrate.hpp.

§ isContiguous()

template<class V >
bool TMBad::isContiguous ( V &  x)

Is this ad vector available as a contiguous block on the tape?

Template type 'V::value_type' can be

  • ad_plain
  • ad_aug
  • ad_adapt

Definition at line 3045 of file global.hpp.

Referenced by forceContiguous().

§ lmatch()

template<class T >
std::vector<bool> TMBad::lmatch ( const std::vector< T > &  x,
const std::vector< T > &  y 
)

Match x vector in y vector.

Returns
Logical vector of the same length as x
Note
Intended for small vectors only

Definition at line 19 of file graph_transform.hpp.

Referenced by TMBad::clique::get_stride(), and match().

§ make_space_inplace()

template<class T , class I >
void TMBad::make_space_inplace ( std::vector< T > &  x,
std::vector< I > &  i,
space = T(0) 
)

Make space for new elements in an existing vector.

Resize and move elements of a vector to make space for index i new members. Vector x and index vector i are modified such that

  • Index vector points at the same elements as before: x_new[i_new] == x_old[i_old] except that i_new is sorted while i_old might not be.
  • Undefined elements are at positions i-1
  • Removing the undefined elements you get the old vector

Definition at line 68 of file graph_transform.hpp.

Referenced by substitute().

§ matmul()

template<bool XT, bool YT, bool ZT, bool UP>
void TMBad::matmul ( Map< const dmatrix >  x,
Map< const dmatrix >  y,
Map< dmatrix >  z 
)

Expand all 16 combinations

Definition at line 92 of file ad_blas.hpp.

§ mvIntegrate()

template<class Integrand >
mvIntegral0<Integrand> TMBad::mvIntegrate ( Integrand &  f,
control  c = control() 
)

Multivariate integration.

Parameters
fMultivariate integrand (functor)
cOptional control parameters

Example:

template<class Float>
struct Gauss2D_t {
typedef Float Scalar;
Float a, b; // Parameters
Float x, y; // Integration variables
// Evaluate integrand (u1,u2)
Float operator() () {
Float ans = exp(- a*x*x - b*y*y);
return ans;
}
// Integrate wrt (x,y)
Float my_integrate() {
Float ans = mvIntegrate(*this).
wrt(x, -INFINITY, INFINITY).
wrt(y, -INFINITY, INFINITY) ();
return ans;
}
};

Definition at line 1199 of file TMBad/integrate.hpp.

§ operator<()

bool TMBad::operator< ( const ad_aug x,
const ad_aug y 
)

Allow / disallow comparison operators for ad_aug

It is recommended to disable comparison for the non-retaping type ad_aug, i.e. not including this file!

Definition at line 4493 of file TMBad.cpp.

§ order()

template<class T >
std::vector<size_t> TMBad::order ( std::vector< T >  x)

Get permutation that sorts a vector.

If x equals

{5, 4, 3, 5, 4, 3, 5, 2, 1}

Then order(x) is

{8, 7, 2, 5, 1, 4, 0, 3, 6}

Note
Order of duplicates within groups is preserved

Definition at line 117 of file graph_transform.hpp.

Referenced by TMBad::autopar::max_tree_depth(), and TMBad::ADFun<>::set_inv_positions().

§ remap_identical_sub_expressions()

std::vector< Index > TMBad::remap_identical_sub_expressions ( global glob,
std::vector< Index >  inv_remap 
)

Remap identical sub-expressions.

Description

Recall that all variables are stored sequentially on the tape and that each variable can be thought of as an 'expression'. The purpuse of the present function is to find out for any given variable whether it has already been calculated and, if so, find its first occurance. Finally, having constructed such a variable remapping table, we apply the remapping to all operator inputs, and return.

This routine is thus only useful when followed by a call to the method global::eliminate, which will remove all redundant variables.

Algorithm

The algorithm essentially consists of two forward passes:

  1. Hash step Calculate a table of likely valid remappings using global::hash_sweep and radix::first_occurance. We call the resulting vector 'remap'. It satisfies remap[i] <= i with equality signifying that the variable i must be kept as is (which is always valid).
  2. Proof step Assume by induction that remap[j] is valid for all j strictly less than i. We must decide if remap[i] is valid, that is if the expressions remap[i] and i are identical. The two expressions are identical if they are result of the same operator and if their inputs are identical expressions. Because these inputs have index smaller than i we know that remap is valid for them. In other words we simply check that remap[inputs(remap[i])] is the same as remap[inputs(i)]. If this is the case we have proved that the expressions are equal and can accept the remapping. Otherwise we reject the remapping by setting remap[i]=i. In any case remap[i] now contains a valid remapping which completes the induction step.

Applying the remap

To summarize the above algorithm, it gives as output a vector remap such that expressions i and remap[i] are identical. We can thus apply the remap to all operator inputs and get an equivalent computational graph.

However, special attention must be payed to operators with pointer inputs because they assume a contiguous memory layout. If an operator assumes a contiguous memory layout of the variables a:b we must require that the remapped variables remap[a:b] are also contiguous. This is done as a post rejection step by setting remap[a:b] := a:b for invalid remappings.

Parameters
globFunction object to be modified
all_allow_remapSkip extra check

Definition at line 4303 of file TMBad.cpp.

§ reorder_depth_first()

void TMBad::reorder_depth_first ( global glob)

Depth-first reordering of computational graph.

Application: Register optimization for GPU

Definition at line 595 of file TMBad.cpp.

Referenced by TMBad::periodic< T >::find_best_period(), and reorder_depth_first().

§ reorder_graph()

void TMBad::reorder_graph ( global glob,
std::vector< Index >  inv_idx 
)

Reorder computational graph such that selected independent variables come last.

Parameters
inv_idxSorted vector of independent variables.
Note
Nothing is done if any operators with pointer inputs are detected on the tape. (FIXME: Not optimal)

Definition at line 4456 of file TMBad.cpp.

Referenced by TMBad::ADFun<>::reorder().

§ reorder_sub_expressions()

void TMBad::reorder_sub_expressions ( global glob)

Re-order computational graph to make it more compressible.

The degree of compression offered by the periodic sequence analysis can be greatly improved by re-ordering the computational graph before the analysis. The re-ordering algorithm works as follows (see also remap_identical_sub_expressions).

  1. Hash step Two sub-expressions are considered 'weakly similar' if they are result of the same operator sequence without necessarily sharing any constants or independent variables. Assume we have a vector of hash codes 'h' - one for each sub-expression. Expressions with the same hash code are likely 'weakly similar'. Now remap:=first_occurance(h) gives a likely valid reordering of variables satisfying remap[i] <= i.
  2. Proof step Assume by induction that remap[j] is valid for all j strictly less than i. Test if remap[i] is also valid (i.e. max(remap[dependencies(v2o(i))]) < remap[i]). If yes keep it. If no reject it be setting remap[i]=i (which is always valid).
  3. Now remap provides a valid ordering of the graph using the permutation radix::order(remap)

Definition at line 537 of file TMBad.cpp.

Referenced by TMBad::periodic< T >::find_best_period(), and reorder_sub_expressions().

§ reorder_temporaries()

void TMBad::reorder_temporaries ( global glob)

Re-order computational graph to make it more compressible.

Potential re-ordering to be applied after reorder_sub_expressions. Pushes all temporaries forward. Here a temporary* is defined as a variable which is only used once as input to an operator.

Definition at line 567 of file TMBad.cpp.

Referenced by TMBad::periodic< T >::find_best_period(), and reorder_temporaries().

§ split_period()

std::vector< period > TMBad::split_period ( global glob,
period  p,
size_t  max_period_size 
)

Helper.

  • Given a operator period, i.e. (begin, size, rep)
  • Given the full input sequence (glob.inputs) Make a split of the operator period such that the inputs are periodic as well. The split is represented as a vector of periods. Algorithm:
  1. Find matrix of input diffs
  2. Row-wise, find all periods and mark the 'period boundaries', i.e. (begin, begin+size*rep) in a common boolean workspace (shared by all rows).
  3. Construct the period split. Each mark signifies the beginning of a new period.

FIXME: This function implementation is still a bit unclear

Definition at line 190 of file TMBad.cpp.

Referenced by TMBad::periodic< T >::find_best_period(), reorder_depth_first(), and split_period().

§ substitute()

std::vector< Index > TMBad::substitute ( global glob,
const std::vector< Index > &  seq,
bool  inv_tags = true,
bool  dep_tags = true 
)

substitute node index sequence by independent variables

Note
Currently all operators in sequence must be different from InvOp

Definition at line 3584 of file TMBad.cpp.

Referenced by accumulation_tree_split(), TMBad::ADFun<>::decompose(), order(), and substitute().

§ value()

template<class T >
double TMBad::value ( x)

Namespace with utility functions for adaptive numerical integration.

Interfaces to R's integrator that can be used with forward mode AD.

Definition at line 15 of file TMBad/integrate.hpp.

Referenced by TMBad::global::clear_array_subgraph().

License: GPL v2