vspline 1.1.0
Generic C++11 Code for Uniform B-Splines
|
Namespaces | |
namespace | detail |
Classes | |
struct | allocator_traits |
vspline creates vigra::MultiArrays of vectorized types. As long as the vectorized types are Vc::SimdArray or vspline::simd_type, using std::allocator is fine, but when using other types, using a specific allocator may be necessary. Currently this is never the case, but I have the lookup of allocator type from this traits class in place if it should become necessary. More... | |
struct | allocator_traits< hwy_simd_type< T, N > > |
struct | allocator_traits< vc_simd_type< T, N > > |
struct | amplify_type |
amplify_type amplifies it's input with a factor. If the data are multi-channel, the factor is multi-channel as well and the channels are amplified by the corresponding elements of the factor. I added this class to make work with integer-valued splines more comfortable - if these splines are prefiltered with 'boost', the effect of the boost has to be reversed at some point, and amplify_type does just that when you use 1/boost as the 'factor'. More... | |
struct | basis_functor |
basis_functor is an object producing the b-spline basis function value for given arguments, or optionally a derivative of the basis function. While basis_functor can produce single basis function values for single arguments, it can also produce a set of basis function values for a given 'delta'. This set is a unit-spaced sampling of the basis function sampled at n + delta for all n E N. Such samplings are used to evaluate b-splines; they constitute the set of weights which have to be applied to a set of b-spline coefficients to form the weighted sum which is the spline's value at a given position. More... | |
struct | bf_grok_type |
If there are several differently-typed basis functors to be combined in a multi_bf_type object, we can erase their type, just like grok_type does for vspline::unary_functors. grokking a basis functor may cost a little bit of performance but it makes the code to handle multi_bf_types simple: instead of having to cope for several, potentially differently-typed per-axis functors there is only one type - which may be a bf_grok_type if the need arises to put differently-typed basis functors into the multi_bf_type. With this mechanism, the code to build evaluators can be kept simple (handling only one uniform type of basis functor used for all axes) and still use different basis functors. More... | |
struct | bracer |
class bracer encodes the entire bracing process. Note that contrary to my initial implementation, class bracer is now used exclusively for populating the frame around a core area of data. It has no code to determine which size a brace/frame should have. This is now determined in class bspline, see especially class bspline's methods get_left_brace_size(), get_right_brace_size() and setup_metrics(). More... | |
struct | broadcast |
struct broadcast is a mixin providing an 'eval' method to a functor which can process vectorized arguments. This mixin is inherited by the functor missing that capability, using CRTP. Because here, in the providing class, nothing is known (or, even, knowable) about the functor, we need to pass additional template arguments to establish the usual vspline unary functor frame of reference, namely in_type, out_type, vsize etc. The resulting 'vectorized' eval may not be efficient: it has to build individual 'in_type' values from the vectorized input, process them with the derived functor's eval routine, then insert the resulting out_type in the vectorized output. But it's a quick way of getting vectorized evaluation capability without writing the vector code. This is particularly useful when the functor's unvectorized eval() is complex (like, calling into legacy code or even into opaque binary) and 'proper' vectorization is hard to do. And with a bit of luck, the optimizer 'recognizes' what's going on and produces SIMD code anyway. Note that the derived class needs a using declaration for the vectorized eval overload inherited from this base class - see broadcast_type (below) for an example of using this mixin. More... | |
struct | broadcast_type |
struct | bspl_prefilter |
class to provide b-spline prefiltering, using 'iir_filter' above. The actual filter object has to interface with the data handling routine ('present', see filter.h). So this class functions as an adapter, combining the code needed to set up adequate buffers and creation of the actual IIR filter itself. The interface to the data handling routine is provided by inheriting from class buffer_handling More... | |
struct | bspline |
class bspline now builds on class bspline_base, adding coefficient storage, while bspline_base provides metadata handling. This separation makes it easier to generate classes which use the same metadata scheme. One concrete example is a class to type-erase a spline (not yet part of the library) which abstracts a spline, hiding the type of it's coefficients. Such a class can inherit from bspline_base and initialize the base class in the c'tor from a given bspline object, resulting in a uniform interface to the metadata. class bspline takes an additional template argument: the value type. This is the type of the coefficients stored in the spline, and it will typically be a single fundamental type or a small aggregate - like a vigra::TinyVector. vspline uses vigra's ExpandElementResult mechanism to inquire for a value type's elementary type and size, which makes it easy to adapt to new value types because the mechanism is traits-based. More... | |
struct | bspline_base |
struct bspline is the object in vspline holding b-spline coefficients. In a way, the b-spline 'is' it's coefficients, since it is totally determined by them - while, of course, the 'actual' spline is an n-dimensional curve. So, even if this is a bit sloppy, I often refer to the coefficients as 'the spline', and have named struct bspline so even if it just holds the coefficients. More... | |
struct | bspline_evaluator_tag |
tag class used to identify all vspline::evaluator instantiations More... | |
class | buffer_handling |
buffer_handling provides services needed for interfacing with a buffer of simdized/goading data. The init() routine receives two views: one to a buffer accepting incoming data, and one to a buffer providing results. Currently, all filters used in vspline operate in-place, but the two-argument form leaves room to manoevre. get() and put() receive 'bundle' arguments which are used to transfer incoming data to the view defined in in_window, and to transfer result data from the view defined in out_window back to target memory. More... | |
struct | bundle |
class 'bundle' holds all information needed to access a set of vsize 1D subarrays of an nD array. This is the data structure we use to tell the buffering and unbuffering code which data we want it to put into the buffer or distribute back out. The buffer itself holds the data in compact form, ready for vector code to access them at maximum speed. More... | |
struct | callable |
mixin 'callable' is used with CRTP: it serves as additional base to unary functors which are meant to provide operator() and takes the derived class as it's first template argument, followed be the argument types and vectorization width, so that the parameter and return type for operator() and - if vsize is greater than 1 - it's vectorized overload can be produced. This formulation has the advantage of not having to rely on the 'out_type_of' mechanism I was using before and provides precisely the operator() overload(s) which are appropriate. More... | |
struct | chain_type |
class chain_type is a helper class to pass one unary functor's result as argument to another one. We rely on T1 and T2 to provide a few of the standard types used in unary functors. Typically, T1 and T2 will both be vspline::unary_functors, but the type requirements could also be fulfilled 'manually'. More... | |
struct | clamp_gate |
clamp gate clamps out-of-bounds values. clamp_gate takes four arguments: the lower and upper limit of the gate, and the values which are returned if the input is outside the range: 'lfix' if it is below 'lower' and 'ufix' if it is above 'upper' More... | |
struct | convolve |
class convolve provides the combination of class fir_filter above with a vector-friendly buffer. Calling code provides information about what should be buffered, the data are sucked into the buffer, filtered, and moved back from there. The operation is orchestrated by the code in filter.h, which is also used to 'wield' the b-spline prefilter. Both operations are sufficiently similar to share the wielding code. More... | |
struct | dimension_mismatch |
dimension-mismatch is thrown if two arrays have different dimensions which should have the same dimensions. More... | |
struct | domain_type |
class domain is a coordinate transformation functor. It provides a handy way to translate an arbitrary range of incoming coordinates to an arbitrary range of outgoing coordinates. This is done with a linear translation function. if the source range is [s0,s1] and the target range is [t0,t1], the translation function s->t is: More... | |
class | evaluator |
class evaluator encodes evaluation of a spline-like object. This is a generalization of b-spline evaluation, which is the default if no other basis functions are specified. Technically, a vspline::evaluator is a vspline::unary_functor, which has the specific capability of evaluating a specific spline-like object. This makes it a candidate to be passed to the functions in transform.h, like remap() and transform(), and it also makes it suitable for vspline's functional constructs like chaining, mapping, etc. More... | |
struct | extrapolator |
struct extrapolator is a helper class providing extrapolated values for a 1D buffer indexed with possibly out-of-range indices. The extrapolated value is returned by value. boundary conditions PERIODIC , MIRROR , REFLECT, NATURAL and CONSTANT are currently supported. An extrapolator is set up by passing the boundary condition code (see common.h) and a const reference to the 1D data set, coded as a 1D vigra::MultiArrayView. The view has to refer to valid data for the time the extrapolator is in use. Now the extrapolator object can be indexed with arbitrary indices, and it will return extrapolated values. The indexing is done with operator() rather than operator[] to mark the semantic difference. Note how buffers with size 1 are treated specially for some boundary conditions: here we simply return the value at index 0. More... | |
struct | fir_filter |
class fir_filter provides the 'solve' routine which convolves a 1D signal with selectable extrapolation. Here, the convolution kernel is applied to the incoming signal and the result is written to the specified output location. Note that this operation can be done in-place, but input and output may also be different. While most of the time this routine will be invoked by class convolve (below), it is also directly used by the specialized code for 1D filtering. Note how we conveniently inherit from the specs class. This also enables us to use an instance of fir_filter or class convolve as specs argument to create further filters with the same arguments. More... | |
struct | fir_filter_specs |
fir_filter_specs holds the parameters for a filter performing a convolution along a single axis. In vspline, the place where the specifications for a filter are fixed and the place where it is finally created are far apart: the filter is created in the separate worker threads. So this structure serves as a vehicle to transport the arguments. Note the specification of 'headroom': this allows for non-symmetrical and even kernels. When applying the kernel to obtain output[i], the kernel is applied to input [ i - headroom ] , ... , input [ i - headroom + ksize - 1 ] More... | |
struct | flip |
flip functor produces it's input with component order reversed. This can be used to deal with situations where coordinates in the 'wrong' order have to be fed to a functor expecting the opposite order and should be a fast way of doing so, since the compiler can likely optimize it well. I added this class to provide simple handling of incoming NumPy coordinates, which are normally in reverse order of vigra coordinates More... | |
struct | grok_type |
class grok_type is a helper class wrapping a vspline::unary_functor so that it's type becomes opaque - a technique called 'type erasure', here applied to vspline::unary_functors with their specific capability of providing both vectorized and unvectorized operation in one common object. More... | |
struct | grok_type< IN, OUT, 1 > |
specialization of grok_type for _vsize == 1 this is the only possible specialization if vectorization is not used. here we don't use _v_ev but only the unvectorized evaluation. More... | |
struct | homogeneous_mbf_type |
homogeneous_mbf_type can be used for cases where all basis functors are the same. The evaluation code uses operator[] to pick the functor for each axis, so here we merely override operator[] to always yield a const reference to the same basis functor. More... | |
class | iir_filter |
class iir_filter implements an n-pole forward/backward recursive filter to be used for b-spline prefiltering. It inherits from the 'specs' class for easy initialization. More... | |
struct | iir_filter_specs |
structure to hold specifications for an iir_filter object. This set of parameters has to be passed through from the calling code through the multithreading code to the worker threads where the filter objects are finally constructed. Rather than passing the parameters via some variadic mechanism, it's more concise and expressive to contain them in a structure and pass that around. The filter itself inherits its specification type, and if the code knows the handler's type, it can derive the spec type. This way the argument passing can be formalized, allowing for uniform handling of several different filter types with the same code. Here we have the concrete parameter set needed for b-spline prefiltering. We'll pass one set of 'specs' per axis; it contains: More... | |
struct | invalid_scalar |
struct | map_functor |
finally we define class mapper which is initialized with a set of gate objects (of arbitrary type) which are applied to each component of an incoming nD coordinate in turn. The trickery with the variadic template argument list is necessary, because we want to be able to combine arbitrary gate types (which have distinct types) to make the mapper as efficient as possible. the only requirement for a gate type is that it has to provide the necessary eval() functions. More... | |
struct | mirror_gate |
mirror gate 'folds' coordinates into the range. From the infinite number of mirror images resulting from mirroring the input on the bounds, the only one inside the range is picked as the result. When using this gate type with splines with MIRROR boundary conditions, if the shape of the core for the axis in question is M, _lower would be passed 0 and _upper M-1. For splines with REFLECT boundary conditions, we'd pass -0.5 to _lower and M-0.5 to upper, since here we mirror 'between bounds' and the defined range is wider. More... | |
struct | multi_bf_type |
When several basis functors have to be passed to an evaluator, it's okay to pass a container type like a std::array or std::vector. All of these basis functors have to be of the same type, but using the method given above ('grokking' basis functors) it's possible to 'slot in' differently-typed basis functors - the functionality required by the evaluation code is provided by the 'grokked' functors and their inner workings remain opaque. More... | |
struct | multi_bf_type< vspline::basis_functor< math_type >, ndim > |
For b-spline processing, we use a multi_bf_type of b-spline basis functors (class vspline::basis_functor). We can't use the unspecialized template for this basis functor, because it takes the derivative specification, which is specified per-axis, so we need to 'pick it out' from the derivative specification and pass the per-axis value to the per-axis basis function c'tors. So here, instead of merely using the index_sequence to produce a sequence of basis function c'tor calls and ignoring the indices, we use the indices to pick out the per-axis derivative specs. More... | |
struct | not_implemented |
for interfaces which need specific implementations we use: More... | |
struct | not_supported |
exception which is thrown if an opertion is requested which vspline does not support More... | |
struct | numeric_overflow |
exception which is thrown when assigning an rvalue which is larger than what the lvalue can hold More... | |
struct | out_of_bounds |
out_of_bounds is thrown by mapping mode REJECT for out-of-bounds coordinates this exception is left without a message, it only has a very specific application, and there it may be thrown often, so we don't want anything slowing it down. More... | |
struct | pass_gate |
class pass_gate passes it's input to it's output unmodified. More... | |
struct | periodic_gate |
the periodic mapping also folds the incoming value into the allowed range. The resulting value will be ( N * period ) from the input value and inside the range, period being upper - lower. For splines done with PERIODIC boundary conditions, if the shape of the core for this axis is M, we'd pass 0 to _lower and M to _upper. More... | |
struct | reject_gate |
reject_gate throws vspline::out_of_bounds for invalid coordinates More... | |
struct | shape_mismatch |
shape mismatch is the exception which is thrown if the shapes of an input array and an output array do not match. More... | |
struct | simd_allocator |
struct | simd_traits |
traits class simd_traits provides three traits: More... | |
struct | simd_type |
class template simd_type provides a fixed-size container type for small sets of fundamentals which are stored in a C vector. The type offers arithmetic capabilities which are implemented by using loops over the elements in the vector, expecting that the compiler will autovectorize these loops into 'proper' SIMD code. The interface of this type is modelled to be compatible with Vc's SimdArray. Unfortunately, Vc::SimdArray requires additional template arguments, so at times it's difficult to use the two types instead of each other. The interface compatibility does not mean that the arithmetic will produce the same results - this is intended but neither tested nor enforced. More... | |
struct | sink_functor |
sink_functor is used for functors without an output - e.g. reductors which are used for analytic purposes on data sets. They use the same system of input types, but omit the output types. More... | |
struct | sink_functor_tag |
while 'normal' unary_functors are all derived from unary_functor_tag, sink functors will be derived from sink_functor_tag. More... | |
struct | unary_functor |
class unary_functor provides a functor object which offers a system of types for concrete unary functors derived from it. If vectorization isn't used, this is trivial, but with vectorization in use, we get vectorized types derived from plain IN and OUT via query of vspline::vector_traits. More... | |
struct | unary_functor_tag |
we derive all vspline::unary_functors from this empty class, to have a common base type for all of them. This enables us to easily check if a type is a vspline::unary_functor without having to wrangle with unary_functor's template arguments. More... | |
struct | vc_simd_type |
class template vc_simd_type provides a fixed-size SIMD type. This implementation of vspline::vc_simd_type uses Vc::SimdArray The 'acrobatics' may seem futile - why inherit privately from Vc::SimdArray, then code a class template which does essentially the same? There are several reasons: first, the wrapper class results in a common interface shared with the other SIMD implementations, second, there are some added members which can't be 'put into' Vc::SimdArray from the outside. And, third, the template signature is uniform, avoiding Vc::SimdArray's two additional template arguments. More... | |
struct | vector_traits |
with the definition of 'simd_traits', we can proceed to implement 'vector_traits': struct vector_traits is a traits class fixing the types used for vectorized code in vspline. These types go beyond mere vectors of fundamentals: most of the time, the data vspline has to process are not fundamentals, but what I call 'xel' data: pixels, voxels, stereo sound samples, etc. - so, small aggregates of a fundamental type. vector_traits defines how fundamentals and 'xel' data are to be vectorized. with the types defined by vector_traits, a system of type names is introduced which uses a set of patterns: More... | |
struct | vector_traits< T, _vsize, typename std::enable_if< vspline::is_element_expandable< T > ::value > ::type > |
specialization of vector_traits for 'element-expandable' types. These types are recognized by vigra's ExpandElementResult mechanism, resulting in the formation of a 'vectorized' version of the type. These data are what I call 'xel' data. As explained above, vectorization is horizontal, so if T is, say, a pixel of three floats, the type generated here will be a TinyVector of three vectors of vsize floats. More... | |
struct | yield_type |
at times we require reading access to an nD array at given coordinates, as a functor which, receiving the coordinates, produces the values from the array. In the scalar case, this is trivial: if the coordinate is integral, we have a simple indexed access, and if it is real, we can use std::round to produce a nearby discrete coordinate. But for the vectorized case we need a bit more effort: We need to translate the access with a vector of coordinates into a gather operation. We start out with a generalized template class 'yield-type': More... | |
struct | yield_type< crd_t, data_t, _vsize, typename std::enable_if< ! crd_integral< crd_t > ::value > ::type > |
Next, we specialize for real coordinates. More... | |
struct | yield_type< crd_t, data_t, _vsize, typename std::enable_if< crd_integral< crd_t > ::value > ::type > |
Typedefs | |
template<class spline_type , typename rc_type > | |
using | bspl_coordinate_type = vspline::canonical_type< rc_type, spline_type::dimension > |
using declaration for a coordinate suitable for bspline, given elementary type rc_type. This produces the elementary type itself if the spline is 1D, a TinyVector of rc_type otherwise. More... | |
template<class spline_type > | |
using | bspl_value_type = typename spline_type::value_type |
using declaration for a bspline's value type More... | |
typedef long double | xlf_type |
template<typename T > | |
using | ET = typename vigra::ExpandElementResult< T > ::type |
using definition for the 'elementary type' of a type via vigra's ExpandElementResult mechanism. Since this is a frequently used idiom in vspline but quite a mouthful, here's an abbreviation: More... | |
template<typename T > | |
using | EN = typename std::integral_constant< int, vigra::ExpandElementResult< T > ::size > ::type |
produce a std::integral_constant from the size obtained from vigra's ExpandElementResult mechanism More... | |
template<class T > | |
using | is_element_expandable = typename std::integral_constant< bool, ! std::is_same< typename vigra::ExpandElementResult< T > ::type, vigra::UnsuitableTypeForExpandElements > ::value > ::type |
is_element_expandable tests if a type T is known to vigra's ExpandElementResult mechanism. If this is so, the type is considered 'element-expandable'. More... | |
template<class T1 , class T2 > | |
using | promote_type = typename vigra::PromoteTraits< T1, T2 > ::Promote |
template<class T1 , class T2 > | |
using | common_math_type = typename vigra::NumericTraits< promote_type< T1, T2 > > ::RealPromote |
template<class T1 , class T2 > | |
using | promote_ele_type = typename vigra::ExpandElementResult< promote_type< T1, T2 > > ::type |
template<class T1 , class T2 > | |
using | common_math_ele_type = typename vigra::NumericTraits< promote_ele_type< T1, T2 > > ::RealPromote |
template<typename T , int N = EN < T > :: value> | |
using | canonical_type = typename std::conditional< N==1, ET< T >, vigra::TinyVector< ET< T >, N > > ::type |
produce the 'canonical type' for a given type T. If T is single-channel, this will be the elementary type itself, otherwise it's a TinyVector of the elementary type. optionally takes the number of elements the resulting type should have, to allow construction from a fundamental and a number of channels. More... | |
template<typename T , int N = EN < T > :: value> | |
using | synthetic_type = vigra::TinyVector< ET< T >, N > |
template<class T , size_t sz = 1> | |
using | scalar = typename std::conditional< sz==1, T, invalid_scalar< T, sz > > ::type |
definition of a scalar with the same template argument list as a simdized type, to use 'scalar' in the same syntactic slot More... | |
template<typename T , size_t SZ> | |
using | tiny_vector = typename vigra::TinyVector< T, static_cast< int >(SZ) > |
template<typename coordinate_type , typename value_type > | |
using | default_math_type = typename vigra::NumericTraits< typename vigra::PromoteTraits< typename vigra::ExpandElementResult< coordinate_type > ::type, typename vigra::ExpandElementResult< value_type > ::type > ::Promote > ::RealPromote |
template<typename _coordinate_type , typename _trg_type , typename _mbf_type , size_t _vsize = vspline::vector_traits < _trg_type > :: size, typename _math_ele_type = default_math_type < _coordinate_type , _trg_type >, typename _cf_type = _trg_type> | |
using | abf_evaluator = evaluator< _coordinate_type, _trg_type, _vsize, -1, _math_ele_type, _cf_type, _mbf_type > |
alias template to make the declaration of non-bspline evaluators easier. Note that we fix the 'specialze' template parameter at -1 (no specialization) and pass the mbf type argument in third position. 'abf' stands for 'alternative basis function'. One example of such an alternative basis function is vspline::area_basis_functor - see basis.h. More... | |
template<typename T , std::size_t N> | |
using | hwy_simd_type = HWY_NAMESPACE::hwy_simd_type< T, N > |
template<typename T > | |
using | atomic = std::atomic< T > |
template<typename T > | |
using | is_scalar = typename std::integral_constant< bool, std::is_fundamental< T > ::value||std::is_same< T, bool > ::value > ::type |
template<unsigned int dimension> | |
using | bcv_type = vigra::TinyVector< vspline::bc_code, dimension > |
a type for a set of boundary condition codes, one per axis More... | |
template<unsigned int dimension, typename rc_ele_type = float> | |
using | grid_spec = vigra::TinyVector< vigra::MultiArrayView< 1, rc_ele_type >, dimension > |
template<typename T > | |
using | crd_integral = typename std::conditional< std::is_fundamental< T >::value, std::is_integral< T >, std::is_integral< typename T::value_type > > ::type |
First, we specialize for integral coordinates. More... | |
template<typename T , size_t N> | |
using | simdized_type = typename vector_traits< T, N > ::type |
this alias is used as a shorthand to pick the vectorized type for a given type T and a size N from 'vector_traits': More... | |
Enumerations | |
enum | bc_code { MIRROR , PERIODIC , REFLECT , NATURAL , CONSTANT , ZEROPAD , GUESS , INVALID } |
This enumeration is used for codes connected to boundary conditions. There are two aspects to boundary conditions: During prefiltering, the initial causal and anticausal coefficients have to be calculated in a way specific to the chosen boundary conditions. Bracing also needs these codes to pick the appropriate extrapolation code to extrapolate the coefficients beyond the core array. More... | |
Functions | |
template<typename ic_t , typename rc_t > | |
void | odd_split (rc_t v, ic_t &iv, rc_t &fv) |
coordinates are split into an integral part and a remainder. this is used for weight generation, and also for calculating the basis function value. The split is done differently for odd and even splines. More... | |
template<typename ic_t , typename rc_t , int N> | |
void | odd_split (vigra::TinyVector< rc_t, N > v, vigra::TinyVector< ic_t, N > &iv, vigra::TinyVector< rc_t, N > &fv) |
template<typename ic_t , typename rc_t > | |
void | even_split (rc_t v, ic_t &iv, rc_t &fv) |
for even splines, the integral part is obtained by rounding. when the result of rounding is subtracted from the original coordinate, a value between -0.5 and 0.5 is obtained which is used for weight generation. More... | |
template<typename ic_t , typename rc_t , int N> | |
void | even_split (vigra::TinyVector< rc_t, N > v, vigra::TinyVector< ic_t, N > &iv, vigra::TinyVector< rc_t, N > &fv) |
template<class real_type > | |
real_type | bspline_basis_2 (int x2, int degree, int derivative=0) |
bspline_basis_2 yields the value of the b-spline basis function at multiples of 1/2, for which vspline has precomputed values. Instead of passing real values x which are multiples of 1/2, this routine takes the doubled argument, so instead of calling it with x = 0.5, you call it with x2 = 1. More... | |
template<class target_type > | |
void | calculate_weight_matrix (vigra::MultiArray< 2, target_type > &res) |
a 'weight matrix' is used to calculate a set of weights for a given remainder part of a real coordinate. The 'weight matrix' is multipled with a vector containing the power series of the given remainder, yielding a set of weights to apply to a window of b-spline coefficients. More... | |
template<class target_type > | |
void | get_kernel (const int °ree, vigra::MultiArrayView< 1, target_type > &kernel, const bool &odd=true) |
this function deposits the reconstruction kernel in the array 'kernel'. This kernel can be used to convolve a set of coefficients, to obtain the original signal. This is a convenience function which merely picks the right values from the precomputed values in precomputed_basis_function_values. if 'odd' is passed false, the result is an even kernel. This kernel can't be used for reconstruction (.5 phase shift), but it's handy to get values half a unit step from the knot points. More... | |
template<class real_type > | |
real_type | cdb_bspline_basis (real_type x, int degree, int derivative=0) |
Implementation of the Cox-de Boor recursion formula to calculate the value of the bspline basis function for arbitrary real x. This code was taken from vigra but modified to take the spline degree as a parameter. Since this routine uses recursion, it's usefulness is limited to smaller degrees. More... | |
template<typename real_type > | |
real_type | gaussian_bspline_basis_approximation (real_type x, int degree) |
Gaussian approximation to B-spline basis function. This routine approximates the basis function of degree spline_degree for real x. I checked for all degrees up to 45. The partition of unity quality of the resulting reconstruction filter is okay for larger degrees, the cumulated error over the covered interval is quite low. Still, as the basis function is never actually evaluated in vspline (whenever it's needed, it is needed at n * 1/2 and we have precomputed values for that) there is not much point in having this function around. I leave the code in for now. More... | |
template<class coordinate_type , size_t _vsize = vspline::vector_traits<coordinate_type>::vsize> | |
vspline::domain_type< coordinate_type, _vsize > | domain (const coordinate_type &in_low, const coordinate_type &in_high, const coordinate_type &out_low, const coordinate_type &out_high) |
factory function to create a domain_type type object from the desired lower and upper fix point for incoming coordinates and the lower and upper fix point for outgoing coordinates. the resulting functor maps incoming coordinates in the range of [in_low,in_high] to coordinates in the range of [out_low,out_high] More... | |
template<class coordinate_type , class spline_type , size_t _vsize = vspline::vector_traits<coordinate_type>::vsize> | |
vspline::domain_type< coordinate_type, _vsize > | domain (const coordinate_type &in_low, const coordinate_type &in_high, const spline_type &bspl) |
factory function to create a domain_type type object from the desired lower and upper reference point for incoming coordinates and a vspline::bspline object providing the lower and upper reference for outgoing coordinates the resulting functor maps incoming coordinates in the range of [ in_low , in_high ] to coordinates in the range of [ bspl.lower_limit() , bspl.upper_limit() ] More... | |
template<class coordinate_type , class spline_type , size_t _vsize = vspline::vector_traits<coordinate_type>::vsize> | |
vspline::domain_type< coordinate_type, _vsize > | domain (const spline_type &bspl, const coordinate_type &out_low, const coordinate_type &out_high) |
factory function to create a domain_type type object from a vspline::bspline object providing the lower and upper reference for incomoing coordinates and the desired lower and upper reference point for outgoing coordinates the resulting functor maps incoming coordinates in the range of [ bspl.lower_limit() , bspl.upper_limit() ] to coordinates in the range of [ out_low , out_high ] More... | |
template<class coordinate_type , class spline_type , size_t _vsize = vspline::vector_traits<coordinate_type>::vsize> | |
vspline::domain_type< coordinate_type, _vsize > | domain (const spline_type &bspl_in, const spline_type &bspl_out) |
factory function to create a domain_type type object from a vspline::bspline object providing the lower and upper reference for incomoing coordinates and a vspline::bspline object providing the lower and upper reference for outgoing coordinates the resulting functor maps incoming coordinates in the range of [ bspl_in.lower_limit() , bspl_in.upper_limit() ] to outgoing coordinates in the range of [ bspl_out.lower_limit() , bspl_out.upper_limit() ] More... | |
template<typename _grokkee_t , typename _delta_et = float, typename _target_et = _delta_et, std::size_t _vsize = vspline::vector_traits < _delta_et > :: size> | |
bf_grok_type< _delta_et, _target_et, _vsize > | bf_grok (const _grokkee_t &grokkee) |
template<class spline_type , typename rc_type = typename vigra::NumericTraits < typename spline_type::ele_type > :: RealPromote, size_t _vsize = vspline::vector_traits < typename spline_type::value_type > :: size, typename math_ele_type = default_math_type < typename spline_type::value_type , rc_type >, typename result_type = typename spline_type::value_type> | |
vspline::grok_type< bspl_coordinate_type< spline_type, rc_type >, result_type, _vsize > | make_evaluator (const spline_type &bspl, vigra::TinyVector< int, spline_type::dimension > dspec=vigra::TinyVector< int, spline_type::dimension >(0), int shift=0) |
make_evaluator is a factory function, producing a functor which provides access to an evaluator object. Evaluation using the resulting object is not intrinsically safe, it's the user's responsibility not to pass coordinates which are outside the spline's defined range. If you need safe access, see 'make_safe_evaluator' below. 'Not safe' in this context means that evaluation at out-of-bounds locations may result in a memory fault or produce wrong or undefined results. Note that vspline's bspline objects are set up as to allow evaluation at the lower and upper limit of the spline and will also tolerate values 'just outside' the bounds to guard against quantization errors. see vspline::bspline for details. More... | |
template<class spline_type , typename rc_type = typename vigra::NumericTraits < typename spline_type::ele_type > :: RealPromote, size_t _vsize = vspline::vector_traits < typename spline_type::value_type > :: size, typename math_ele_type = default_math_type < typename spline_type::value_type , rc_type >, typename result_type = typename spline_type::value_type> | |
vspline::grok_type< bspl_coordinate_type< spline_type, rc_type >, result_type, _vsize > | make_safe_evaluator (const spline_type &bspl, vigra::TinyVector< int, spline_type::dimension > dspec=vigra::TinyVector< int, spline_type::dimension >(0), int shift=0) |
make_safe_evaluator is a factory function, producing a functor which provides safe access to an evaluator object. This functor will map incoming coordinates into the spline's defined range, as given by the spline with it's lower_limit and upper_limit methods, honoring the bspline objects's boundary conditions. So if, for example, the spline is periodic, all incoming coordinates are valid and will be mapped to the first period. Note the use of lower_limit and upper_limit. These values also depend on the spline's boundary conditions, please see class vspline::bspline for details. If there is no way to meaningfully fold a coordinate into the defined range, the coordinate is clamped to the nearest limit. More... | |
template<class stype , class ttype , size_t vsize> | |
void | move (const bundle< stype, vsize > &src, ttype *trg, std::ptrdiff_t trg_stride) |
move bundle data to compact memory More... | |
template<class stype , class ttype , size_t vsize> | |
void | move (const bundle< stype, vsize > &src, ttype *trg, std::ptrdiff_t trg_stride, int ni) |
template<class stype , class ttype , size_t vsize> | |
void | move (const stype *src, std::ptrdiff_t src_stride, const bundle< ttype, vsize > &trg) |
move data from compact memory to bundle More... | |
template<class stype , class ttype , size_t vsize> | |
void | move (const stype *src, std::ptrdiff_t src_stride, const bundle< ttype, vsize > &trg, int ni) |
template<typename in_type , typename out_type , unsigned int D, class filter_type , typename ... types> | |
void | filter (const vigra::MultiArrayView< D, in_type > &input, vigra::MultiArrayView< D, out_type > &output, types ... args) |
vspline::filter is the common entry point for filter operations in vspline. This routine does not yet do any processing, it's purpose is to convert it's arguments to 'canonical' format and then call the actual filter code in namespace detail. It also determines the type used for arithmetic operations. The type specification for input and output assures that only arrays with the same dimensionality are accepted, and a static assertion makes sure the number of channels match. canonical form means that input and output value type are either fundamental (for single-channel data) or TinyVectors of a fundamental data type. This way, the input and output is presented in a neutral form, ignoring all specifics of T1 and T2, which may be TinyVectors, complex, RGBValue etc. More... | |
template<unsigned int dimension, typename in_value_type , typename out_value_type , typename math_ele_type = typename vspline::common_math_ele_type < in_value_type , out_value_type >, size_t vsize = vspline::vector_traits < math_ele_type > :: size> | |
void | forward_backward_recursive_filter (const vigra::MultiArrayView< dimension, in_value_type > &input, vigra::MultiArrayView< dimension, out_value_type > &output, vigra::TinyVector< bc_code, static_cast< int >(dimension) > bcv, std::vector< vspline::xlf_type > pv, xlf_type tolerance=-1, xlf_type boost=xlf_type(1), int axis=-1, int njobs=default_njobs) |
forward_backward_recursive_filter applies one or more pairs of simple recursive filters to the input. Each pair consists of: More... | |
template<unsigned int dimension, typename in_value_type , typename out_value_type , typename math_ele_type = typename vspline::common_math_ele_type < in_value_type , out_value_type >, size_t vsize = vspline::vector_traits < math_ele_type > :: size> | |
void | convolution_filter (const vigra::MultiArrayView< dimension, in_value_type > &input, vigra::MultiArrayView< dimension, out_value_type > &output, vigra::TinyVector< bc_code, static_cast< int >(dimension) > bcv, std::vector< vspline::xlf_type > kv, int headroom, int axis=-1, int njobs=default_njobs) |
convolution_filter implements convolution of the input with a fixed-size convolution kernel. Note that vspline does not follow the DSP convention of using the kernel's coefficients in reverse order. The standard is to calculate sum ( ck * u(n-k) ), vspline uses sum ( ck * u(n+k-h) ) where 'h' is the 'headroom' of the kernel - the number of coefficients which are applied 'to 'past' values, so that for for a kernel with three coefficients and headroom 1, the sum at position n would be c0 * u(n-1) + c1 * u(n) + c2 * u(n+1), and for a kernel with four coefficients and headroom 2, you'd get c0 * u(n-2) + c1 * u(n-1) + c2 * u(n) + c3 * u(n+1) If you use a 'normal' kernel in vspline, you must reverse it (unless it's symmetric, of course), and you must state the 'headroom': if your kernel is odd-sized, this will normally be half the kernel size (integer division!), producing the phase-correct result, and with an even kernel you can't get the phase right, and you have to make up your mind which way the phase should shift. The 'headroom' notation leaves you free to pick your choice. The input, output, bcv, axis and njobs parameters are as in the routine above. Here what's needed here is 'kv', a std::vector holding the filter's coefficients, and the 'headroom', usually kv.size() / 2. More... | |
template<typename rc_type , size_t _vsize = vspline::vector_traits < rc_type > :: size> | |
vspline::pass_gate< rc_type, _vsize > | pass () |
factory function to create a pass_gate type functor More... | |
template<typename rc_type , size_t _vsize = vspline::vector_traits < rc_type > :: size> | |
vspline::reject_gate< rc_type, _vsize > | reject (rc_type lower, rc_type upper) |
factory function to create a reject_gate type functor given a lower and upper limit for the allowed range. More... | |
template<typename rc_type , size_t _vsize = vspline::vector_traits < rc_type > :: size> | |
vspline::clamp_gate< rc_type, _vsize > | clamp (rc_type lower, rc_type upper, rc_type lfix, rc_type rfix) |
factory function to create a clamp_gate type functor given a lower and upper limit for the allowed range, and, optionally, the values to use if incoming coordinates are out-of-range More... | |
template<typename rc_type , size_t _vsize = vspline::vector_traits < rc_type > :: size> | |
vspline::clamp_gate< rc_type, _vsize > | clamp (rc_type lower, rc_type upper) |
template<typename rc_v > | |
rc_v | v_fmod (rc_v lhs, const typename rc_v::value_type &rhs) |
vectorized fmod function using std::trunc, which is fast, but checking the result to make sure it's always <= rhs. More... | |
template<typename rc_type , size_t _vsize = vspline::vector_traits < rc_type > :: size> | |
vspline::mirror_gate< rc_type, _vsize > | mirror (rc_type lower, rc_type upper) |
factory function to create a mirror_gate type functor given a lower and upper limit for the allowed range. More... | |
template<typename rc_type , size_t _vsize = vspline::vector_traits < rc_type > :: size> | |
vspline::periodic_gate< rc_type, _vsize > | periodic (rc_type lower, rc_type upper) |
factory function to create a periodic_gate type functor given a lower and upper limit for the allowed range. More... | |
template<typename nd_rc_type , size_t _vsize = vspline::vector_traits < nd_rc_type > :: size, class ... gate_types> | |
vspline::map_functor< nd_rc_type, _vsize, gate_types... > | mapper (gate_types ... args) |
factory function to create a mapper type functor given a set of gate_type objects. Please see vspline::make_safe_evaluator for code to automatically create a mapper object suitable for a specific vspline::bspline. More... | |
template<typename index_t > | |
bool | fetch_descending (vspline::atomic< index_t > &source, index_t &index) |
fetch_descending fetches the next index from an atomic, counting down. Indexes will range from one less than the value the atomic was initialized with, down to zero. If all indexes were distributed already, false is returned, true otherwise. Like the other fetch_XXX routines, the return of a boolean makes these functions good candidates to be used as conditional in a loop. More... | |
template<typename index_t > | |
bool | fetch_ascending (vspline::atomic< index_t > &source, const index_t &total, index_t &index) |
fetch_ascending counts up from zero to total-1, which is more efficient if the indexes are used to address memory. This is due to the side effects of accessing memory: if memory is accessed at an address x, the OS will typically fetch a chunk of data starting at or shortly before x. If the next fetch requires data just after x, there is a good chance that they are already in cache. More... | |
template<typename index_t > | |
bool | fetch_range_descending (vspline::atomic< index_t > &source, const index_t &count, index_t &low, index_t &high) |
fetch_range_ascending fetches the beginning and end of a range of indexes (in iterator semantic, low <= i < high) from a vspline::atomic which has been initialized with the total number of indexes that are to be processed. If the vspline::atomic, when accessed, already holds a value of or below zero, fetch_index_range returns false and leaves low and high unchanged. Otherwise it returns true and low and high will be set to the values gleaned from the atomic, raising 'low 'to zero if it would come out below. this function (and the next) enable calling code to process batches of adjacent indexes without any index artistry: the code is perfectly general and simple, and the use of the atomic and fetch_sub garantees that each fetch provides a distinct batch of indexes. More... | |
template<typename index_t > | |
bool | fetch_range_ascending (vspline::atomic< index_t > &source, const index_t &count, const index_t &total, index_t &low, index_t &high) |
fetch_range_ascending also uses an atomic initialized to the total number of indexes to be distributed, but the successive ranges are handed out in ascending order, which is more efficient if the indexes are used to address memory. More... | |
template<bool dummy = true> | |
int | multithread (std::function< void() > payload, std::size_t nr_workers=default_njobs) |
multithread uses a thread pool of worker threads to perform a multithreaded operation. It receives a functor (a single-threaded function used for all individual tasks), and, optionally, the desired number of worker instances to be used. These tasks are wrapped with a wrapper which takes care of signalling when the last task has completed. More... | |
template<unsigned int dimension, typename in_value_type , typename out_value_type , typename math_ele_type > | |
void | amplify (const vigra::MultiArrayView< dimension, in_value_type > &input, vigra::MultiArrayView< dimension, out_value_type > &output, math_ele_type boost=1, int njobs=vspline::default_njobs) |
amplify is used to copy input to output, optionally applying 'boost' in the process. If the operation is in-place and 'boost' is 1, 'amplify' returns prematurely. More... | |
template<unsigned int dimension, typename in_value_type , typename out_value_type , typename math_ele_type = typename vspline::common_math_ele_type < in_value_type , out_value_type >, size_t vsize = vspline::vector_traits < math_ele_type > :: size> | |
void | prefilter (const vigra::MultiArrayView< dimension, in_value_type > &input, vigra::MultiArrayView< dimension, out_value_type > &output, vigra::TinyVector< bc_code, static_cast< int >(dimension) > bcv, int degree, xlf_type tolerance=std::numeric_limits< math_ele_type > ::epsilon(), xlf_type boost=xlf_type(1), int njobs=default_njobs) |
'prefilter' handles b-spline prefiltering for the whole range of acceptable input and output. It combines two bodies of code to achieve this goal: More... | |
template<typename T , std::size_t vsize> | |
bool | any_of (simd_type< T, vsize > arg) |
template<typename T , std::size_t vsize> | |
bool | all_of (simd_type< T, vsize > arg) |
template<typename T , std::size_t vsize> | |
bool | none_of (simd_type< T, vsize > arg) |
template<typename unary_functor_type , unsigned int dimension> | |
void | transform (const unary_functor_type &functor, const vigra::MultiArrayView< dimension, typename unary_functor_type::in_type > &input, vigra::MultiArrayView< dimension, typename unary_functor_type::out_type > &output, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
implementation of two-array transform using wielding::coupled_wield. More... | |
template<typename sink_functor_type , unsigned int dimension> | |
void | reduce (const sink_functor_type &functor, const vigra::MultiArrayView< dimension, typename sink_functor_type::in_type > &input, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
reduce function processing by value. sink_functor_type is a functor taking values (which are loaded from an array). The reduce function, by delegation to wielding::value_reduce, performs the feeding part by parcelling the data, the functor receives either single in_type or vectorized data. The functor has to be a construct capable of accumulating data with several threads in parallel. A typical example would be a functor cumulating to some internal structure and 'offloading' to a mutex-protected pooling target during it's destructor. The wielding code creates a copy of the functor which is passed here for each worker thread and guarantees that this copy's destructor is executed before the call to reduce returns, so the user can be assured of having collected the partial results of all worker threads in the pooling target. The functor needs a copy constructor to communicate the pooling target to the copies made for the individual threads, see the index-based reduce function below this one for an illustration of the simple and elegant mechanism. Note that sink_functor will normally be a type inheriting from vspline::sink_functor, a close relative of vspline::unary_functor and also defined in unary_functor.h. This inheritance does not produce any functionality, just the standard vspline functor argument type system, which this code relies on. Of course the sink functor can be coded 'manually' as well - inheriting form vspline::sink_functor is merely convenient. More... | |
template<class unary_functor_type > | |
void | transform (const unary_functor_type &functor, vigra::MultiArrayView< unary_functor_type::dim_in, typename unary_functor_type::out_type > &output, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
reduce function processing by coordinate. sink_functor_type is a functor receiving integral coordinates, with no fixed meaning, but usually pertaining to an array. The reduce function feeds all coordinates encompassed by the shape, either in vectorized or single form. Here's an example of a functor for this function: More... | |
template<typename sink_functor_type > | |
void | reduce (const sink_functor_type &functor, typename std::conditional<(sink_functor_type::dim_in==1), int, vigra::TinyVector< int, sink_functor_type::dim_in > > ::type shape, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
while the function above fills an array by calling a functor with the target coordinates, this function accumulates the results by means of a sink_functor. The effect is as if an index-based transform was executed first, filling an intermediate array, followed by a call to reduce taking the intermediate array as input. The combination of both operations 'cuts out' the intermediate. Since this function does not have a target array, the shape has to be communicated explicitly. The functor will be invoked for every possible coordinate within the shape's scope. This is done by the function 'index_reduce' in wielding.h the type of the argument 'shape' is derived from the sink_functor's 'dim_in' value, to make the code 'respect the singular': if the data are 1D, reduce will accept a single integer as 'shape', whereas nD data require a TinyVector of int. For suitable sink_functors, see the remarks and example given above in the comments for the array-base overload of 'reduce'. More... | |
template<typename unary_functor_type , unsigned int dimension> | |
void | apply (const unary_functor_type &ev, vigra::MultiArrayView< dimension, typename unary_functor_type::out_type > &output, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
we code 'apply' as a special variant of 'transform' where the output is also used as input, so the effect is to feed the unary functor each 'output' value in turn, let it process it and store the result back to the same location. While this looks like a rather roundabout way of performing an apply, it has the advantage of using the same type of functor (namely one with const input and writable output), rather than a different functor type which modifies it's argument in-place. While, at this level, using such a functor looks like a feasible idea, It would require specialized code 'further down the line' when complex functors are built with vspline's functional programming tools: the 'apply-capable' functors would need to read the output values first and write them back after anyway, resulting in the same sequence of loads and stores as we get with the current 'fake apply' implementation. From the direct delegation to the two-array overload of 'transform', you can see that this is merely syntactic sugar. More... | |
template<typename original_type , typename result_type , typename coordinate_type , unsigned int cf_dimension, unsigned int trg_dimension, int bcv_dimension = cf_dimension> | |
void | remap (const vigra::MultiArrayView< cf_dimension, original_type > &input, const vigra::MultiArrayView< trg_dimension, coordinate_type > &coordinates, vigra::MultiArrayView< trg_dimension, result_type > &output, bcv_type< bcv_dimension > bcv=bcv_type< bcv_dimension >(MIRROR), int degree=3, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
Implementation of 'classic' remap, which directly takes an array of values and remaps it, internally creating a b-spline of given order just for the purpose. This is used for one-shot remaps where the spline isn't reused, and specific to b-splines, since the functor used is a b-spline evaluator. The spline defaults to a cubic b-spline with mirroring on the bounds. More... | |
template<typename functor_type > | |
void | transform (grid_spec< functor_type::dim_in, typename functor_type::in_ele_type > &grid, const functor_type &functor, vigra::MultiArrayView< functor_type::dim_in, typename functor_type::out_type > &result, int njobs, vspline::atomic< bool > *p_cancel, std::false_type) |
This overload of 'transform' takes a vspline::grid_spec as it's first parameter. This object defines a grid where each grid point is made up from components taken from the per-axis grid values. Let's say the grid_spec holds. More... | |
template<typename functor_type > | |
void | transform (grid_spec< 1, typename functor_type::in_ele_type > &grid, const functor_type &functor, vigra::MultiArrayView< 1, typename functor_type::out_type > &result, int njobs, vspline::atomic< bool > *p_cancel, std::true_type) |
overload of grid-based transform for 1D grids. Obviously, this is a case for using 'transform' directly taking the single set of per-axis grid values as input, and this overload is here only for completeness' sake - 'normal' user code will use an array-based transform straight away. More... | |
template<typename functor_type > | |
void | transform (grid_spec< functor_type::dim_in, typename functor_type::in_ele_type > &grid, const functor_type &functor, vigra::MultiArrayView< functor_type::dim_in, typename functor_type::out_type > &result, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
template<typename functor_type > | |
void | grid_eval (grid_spec< functor_type::dim_in, typename functor_type::in_ele_type > &grid, const functor_type &functor, vigra::MultiArrayView< functor_type::dim_in, typename functor_type::out_type > &result, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
for backward compatibility, we provide 'grid_eval'. This is deprecated and will eventually go away, new code should use 'transform' above. More... | |
template<typename functor_type > | |
void | gen_grid_eval (grid_spec< functor_type::dim_in, typename functor_type::in_ele_type > &grid, const functor_type &functor, vigra::MultiArrayView< functor_type::dim_in, typename functor_type::out_type > &result, int njobs=vspline::default_njobs, vspline::atomic< bool > *p_cancel=0) |
for backward compatibility, we provide 'gen_grid_eval'. This is deprecated and will eventually go away, new code should use 'transform' above. More... | |
template<unsigned int dimension, typename value_type , typename math_ele_type = typename vigra::NumericTraits < typename vigra::ExpandElementResult < value_type > :: type > :: RealPromote, size_t vsize = vspline::vector_traits < math_ele_type > :: size> | |
void | restore (const vspline::bspline< value_type, dimension > &bspl, vigra::MultiArrayView< dimension, value_type > &target) |
restore restores the original data from the b-spline coefficients. This is done efficiently using a separable convolution, the kernel is simply a unit-spaced sampling of the basis function. Since the filter uses internal buffering, using this routine in-place is safe - meaning that 'target' may be bspl.core itself. math_type, the data type for performing the actual maths on the buffered data, and the type the data are converted to when they are placed into the buffer, can be chosen, but I could not detect any real benefits from using anything but the default, which is to leave the data in their 'native' type. More... | |
template<int dimension, typename value_type , typename math_ele_type = typename vigra::NumericTraits < typename vigra::ExpandElementResult < value_type > :: type > :: RealPromote, size_t vsize = vspline::vector_traits<math_ele_type>::size> | |
void | restore (vspline::bspline< value_type, dimension > &bspl) |
overload of 'restore' writing the result of the operation back to the array which is passed in. This looks like an in-place operation, but the data are in fact moved to a buffer stripe by stripe, then the arithmetic is done on the buffer and finally the buffer is written back. This is repeated for each dimension of the array. More... | |
template<class IN , class OUT > | |
std::function< void(const IN &, OUT &) > | eval_wrap (std::function< OUT(const IN &) > f) |
eval_wrap is a helper function template, wrapping an 'ordinary' function which returns some value given some input in a void function taking input as const reference and writing output to a reference, which is the signature used for evaluation in vspline::unary_functors. More... | |
template<class T1 , class T2 > | |
vspline::chain_type< T1, T2 > | chain (const T1 &t1, const T2 &t2) |
chain is a factory function yielding the result of chaining two unary_functors. More... | |
template<typename T1 , typename T2 , typename enable = typename std::enable_if < std::is_base_of < vspline::unary_functor_tag < T2::vsize > , T1 > :: value && std::is_base_of < vspline::unary_functor_tag < T1::vsize > , T2 > :: value > :: type> | |
vspline::chain_type< T1, T2 > | operator+ (const T1 &t1, const T2 &t2) |
using operator overloading, we can exploit operator+'s semantics to chain several unary functors. We need to specifically enable this for types derived from unary_functor_tag to avoid a catch-all situation. More... | |
template<class grokkee_type > | |
vspline::grok_type< typename grokkee_type::in_type, typename grokkee_type::out_type, grokkee_type::vsize > | grok (const grokkee_type &grokkee) |
grok() is the corresponding factory function, wrapping grokkee in a vspline::grok_type. More... | |
template<typename T , typename U > | |
void | assign (T &t, const U &u) |
template<typename T , typename U , int N> | |
void | assign (vigra::TinyVector< T, N > &t, const vigra::TinyVector< U, N > &u) |
template<typename VT1 , typename PT , typename VT2 > | |
void | assign_if (VT1 &target, const PT &predicate, const VT2 &source) |
template<typename T > | |
void | assign_if (T &target, const bool &predicate, const T &source) |
Variables | |
const std::string | bc_name [] |
bc_name is for diagnostic output of bc codes More... | |
const int | ncores = vspline_threadpool::ncores |
const int | default_njobs = vspline_threadpool::default_njobs |
using vspline::abf_evaluator = typedef evaluator < _coordinate_type , _trg_type , _vsize , -1 , _math_ele_type , _cf_type , _mbf_type > |
alias template to make the declaration of non-bspline evaluators easier. Note that we fix the 'specialze' template parameter at -1 (no specialization) and pass the mbf type argument in third position. 'abf' stands for 'alternative basis function'. One example of such an alternative basis function is vspline::area_basis_functor - see basis.h.
using vspline::atomic = typedef std::atomic < T > |
Definition at line 224 of file multithread.h.
using vspline::bcv_type = typedef vigra::TinyVector < vspline::bc_code , dimension > |
a type for a set of boundary condition codes, one per axis
Definition at line 554 of file transform.h.
using vspline::bspl_coordinate_type = typedef vspline::canonical_type < rc_type , spline_type::dimension > |
using vspline::bspl_value_type = typedef typename spline_type::value_type |
using vspline::canonical_type = typedef typename std::conditional < N == 1 , ET < T > , vigra::TinyVector < ET < T > , N > > :: type |
produce the 'canonical type' for a given type T. If T is single-channel, this will be the elementary type itself, otherwise it's a TinyVector of the elementary type. optionally takes the number of elements the resulting type should have, to allow construction from a fundamental and a number of channels.
using vspline::common_math_ele_type = typedef typename vigra::NumericTraits < promote_ele_type < T1 , T2 > > :: RealPromote |
using vspline::common_math_type = typedef typename vigra::NumericTraits < promote_type < T1 , T2 > > :: RealPromote |
using vspline::crd_integral = typedef typename std::conditional < std::is_fundamental<T>::value , std::is_integral<T> , std::is_integral<typename T::value_type> > :: type |
First, we specialize for integral coordinates.
Definition at line 1063 of file unary_functor.h.
using vspline::default_math_type = typedef typename vigra::NumericTraits < typename vigra::PromoteTraits < typename vigra::ExpandElementResult < coordinate_type > :: type , typename vigra::ExpandElementResult < value_type > :: type > :: Promote > :: RealPromote |
using vspline::EN = typedef typename std::integral_constant < int , vigra::ExpandElementResult < T > :: size > :: type |
using vspline::ET = typedef typename vigra::ExpandElementResult < T > :: type |
using vspline::grid_spec = typedef vigra::TinyVector < vigra::MultiArrayView < 1 , rc_ele_type > , dimension > |
Definition at line 697 of file transform.h.
using vspline::hwy_simd_type = typedef HWY_NAMESPACE::hwy_simd_type < T , N > |
Definition at line 2057 of file hwy_simd_type.h.
using vspline::is_element_expandable = typedef typename std::integral_constant < bool , ! std::is_same < typename vigra::ExpandElementResult < T > :: type , vigra::UnsuitableTypeForExpandElements > :: value > :: type |
using vspline::is_scalar = typedef typename std::integral_constant < bool , std::is_fundamental < T > ::value || std::is_same < T , bool > :: value > :: type |
Definition at line 128 of file simd_type.h.
using vspline::promote_ele_type = typedef typename vigra::ExpandElementResult < promote_type < T1 , T2 > > :: type |
using vspline::promote_type = typedef typename vigra::PromoteTraits < T1 , T2 > :: Promote |
using vspline::scalar = typedef typename std::conditional < sz == 1 , T , invalid_scalar < T , sz > > :: type |
using vspline::simdized_type = typedef typename vector_traits < T , N > :: type |
this alias is used as a shorthand to pick the vectorized type for a given type T and a size N from 'vector_traits':
using vspline::synthetic_type = typedef vigra::TinyVector < ET < T > , N > |
using vspline::tiny_vector = typedef typename vigra::TinyVector < T , static_cast < int > ( SZ ) > |
typedef long double vspline::xlf_type |
enum vspline::bc_code |
This enumeration is used for codes connected to boundary conditions. There are two aspects to boundary conditions: During prefiltering, the initial causal and anticausal coefficients have to be calculated in a way specific to the chosen boundary conditions. Bracing also needs these codes to pick the appropriate extrapolation code to extrapolate the coefficients beyond the core array.
Enumerator | |
---|---|
MIRROR | |
PERIODIC | |
REFLECT | |
NATURAL | |
CONSTANT | |
ZEROPAD | |
GUESS | |
INVALID |
Definition at line 919 of file simd_type.h.
void vspline::amplify | ( | const vigra::MultiArrayView< dimension, in_value_type > & | input, |
vigra::MultiArrayView< dimension, out_value_type > & | output, | ||
math_ele_type | boost = 1 , |
||
int | njobs = vspline::default_njobs |
||
) |
amplify is used to copy input to output, optionally applying 'boost' in the process. If the operation is in-place and 'boost' is 1, 'amplify' returns prematurely.
Definition at line 994 of file prefilter.h.
Definition at line 910 of file simd_type.h.
void vspline::apply | ( | const unary_functor_type & | ev, |
vigra::MultiArrayView< dimension, typename unary_functor_type::out_type > & | output, | ||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
we code 'apply' as a special variant of 'transform' where the output is also used as input, so the effect is to feed the unary functor each 'output' value in turn, let it process it and store the result back to the same location. While this looks like a rather roundabout way of performing an apply, it has the advantage of using the same type of functor (namely one with const input and writable output), rather than a different functor type which modifies it's argument in-place. While, at this level, using such a functor looks like a feasible idea, It would require specialized code 'further down the line' when complex functors are built with vspline's functional programming tools: the 'apply-capable' functors would need to read the output values first and write them back after anyway, resulting in the same sequence of loads and stores as we get with the current 'fake apply' implementation. From the direct delegation to the two-array overload of 'transform', you can see that this is merely syntactic sugar.
Definition at line 531 of file transform.h.
void vspline::assign | ( | T & | t, |
const U & | u | ||
) |
void vspline::assign | ( | vigra::TinyVector< T, N > & | t, |
const vigra::TinyVector< U, N > & | u | ||
) |
void vspline::assign_if | ( | T & | target, |
const bool & | predicate, | ||
const T & | source | ||
) |
void vspline::assign_if | ( | VT1 & | target, |
const PT & | predicate, | ||
const VT2 & | source | ||
) |
bf_grok_type< _delta_et, _target_et, _vsize > vspline::bf_grok | ( | const _grokkee_t & | grokkee | ) |
real_type vspline::bspline_basis_2 | ( | int | x2, |
int | degree, | ||
int | derivative = 0 |
||
) |
bspline_basis_2 yields the value of the b-spline basis function at multiples of 1/2, for which vspline has precomputed values. Instead of passing real values x which are multiples of 1/2, this routine takes the doubled argument, so instead of calling it with x = 0.5, you call it with x2 = 1.
this is a helper routine to offer convenient access to vspline's precomputed basis function values, and it also handles the special case of degree 0, x = -1/2, where the b-spline basis function is not symmetric and yields 1.
Inside vspline, this routine is only used by calculate_weight_matrix. User code should rarely need it, but I hold on to it as a separate entity. The code also handles derivatives. This is done by recursion, which is potentially very slow (for high derivatives), so this routine should only be used for 'moderate' derivative values. For fast access to the basis function's derivatives, use of vspline::basis_functor is recommended instead.
void vspline::calculate_weight_matrix | ( | vigra::MultiArray< 2, target_type > & | res | ) |
a 'weight matrix' is used to calculate a set of weights for a given remainder part of a real coordinate. The 'weight matrix' is multipled with a vector containing the power series of the given remainder, yielding a set of weights to apply to a window of b-spline coefficients.
The routine 'calculate_weight_matrix' originated from vigra, but I rewrote it to avoid calculating values of derivatives of the basis function by using recursion. The content of the 'weight matrix' is now calculated directly with a forward iteration starting from precomputed basis function values, the derivatives needed are formed by repeatedly differencing these values. The recursive formulation in vigra makes sense, since the degree of the spline is a template argument in vigra and the recursive formulation can be evaluated at compile-time, allowing for certain optimizations. But in vspline, spline degree is a run-time variable, and vspline offers calculation of splines of degrees up to (currently) 45, which aren't feasibly calculable by recursion, with a recursion calling itself twice for every step: this would take a very long time or exceed the system's capacity. Analyzing the recursive implementation, it can be seen that it produces a great many redundant calculations, which soon exceed reasonable limits. With vigra's maximal spline degree the load is just about manageable, but beyond that (degree 19 or so at the time of this writing) it's clearly no-go.
The forward iteration is reasonably fast, even for high spline degrees, while my previous implementation did slow down very noticeably from, say, degree 20, making it unusable for high spline degrees. I really only noticed the problem after raising the maximal degree vspline can handle, following the rewrite of bootstrap.cc using arbitrary-precision maths.
A weight matrix is also used by 'basis_functor', in a similar way to how it's used for weight generation.
real_type vspline::cdb_bspline_basis | ( | real_type | x, |
int | degree, | ||
int | derivative = 0 |
||
) |
Implementation of the Cox-de Boor recursion formula to calculate the value of the bspline basis function for arbitrary real x. This code was taken from vigra but modified to take the spline degree as a parameter. Since this routine uses recursion, it's usefulness is limited to smaller degrees.
This routine operates in real and calculates the basis function value for arbitrary real x, but it suffers from cumulating errors, especially when the recursion is deep, so the results are not uniformly precise.
This code is expensive for higher spline orders because the routine calls itself twice recursively, so the performance is 2^N with the spline's degree. Luckily there are ways around using this routine at all
I leave this code in here for reference purposes - it's good to have another route to the basis function values, see self_test.cc.
vspline::chain_type< T1, T2 > vspline::chain | ( | const T1 & | t1, |
const T2 & | t2 | ||
) |
chain is a factory function yielding the result of chaining two unary_functors.
Definition at line 588 of file unary_functor.h.
vspline::clamp_gate< rc_type, _vsize > vspline::clamp | ( | rc_type | lower, |
rc_type | upper | ||
) |
vspline::clamp_gate< rc_type, _vsize > vspline::clamp | ( | rc_type | lower, |
rc_type | upper, | ||
rc_type | lfix, | ||
rc_type | rfix | ||
) |
factory function to create a clamp_gate type functor given a lower and upper limit for the allowed range, and, optionally, the values to use if incoming coordinates are out-of-range
void vspline::convolution_filter | ( | const vigra::MultiArrayView< dimension, in_value_type > & | input, |
vigra::MultiArrayView< dimension, out_value_type > & | output, | ||
vigra::TinyVector< bc_code, static_cast< int >(dimension) > | bcv, | ||
std::vector< vspline::xlf_type > | kv, | ||
int | headroom, | ||
int | axis = -1 , |
||
int | njobs = default_njobs |
||
) |
convolution_filter implements convolution of the input with a fixed-size convolution kernel. Note that vspline does not follow the DSP convention of using the kernel's coefficients in reverse order. The standard is to calculate sum ( ck * u(n-k) ), vspline uses sum ( ck * u(n+k-h) ) where 'h' is the 'headroom' of the kernel - the number of coefficients which are applied 'to 'past' values, so that for for a kernel with three coefficients and headroom 1, the sum at position n would be c0 * u(n-1) + c1 * u(n) + c2 * u(n+1), and for a kernel with four coefficients and headroom 2, you'd get c0 * u(n-2) + c1 * u(n-1) + c2 * u(n) + c3 * u(n+1) If you use a 'normal' kernel in vspline, you must reverse it (unless it's symmetric, of course), and you must state the 'headroom': if your kernel is odd-sized, this will normally be half the kernel size (integer division!), producing the phase-correct result, and with an even kernel you can't get the phase right, and you have to make up your mind which way the phase should shift. The 'headroom' notation leaves you free to pick your choice. The input, output, bcv, axis and njobs parameters are as in the routine above. Here what's needed here is 'kv', a std::vector holding the filter's coefficients, and the 'headroom', usually kv.size() / 2.
example: let's say you have data in 'image' and want to convolve in-place with [ 1 , 2 , 1 ], using mirror boundary conditions. Then call:
vspline::convolution_filter ( image , image , { vspline::MIRROR , vspline::MIRROR } , { 1 , 2 , 1 } , 1 ) ;
Definition at line 271 of file general_filter.h.
vspline::domain_type< coordinate_type, _vsize > vspline::domain | ( | const coordinate_type & | in_low, |
const coordinate_type & | in_high, | ||
const coordinate_type & | out_low, | ||
const coordinate_type & | out_high | ||
) |
factory function to create a domain_type type object from the desired lower and upper fix point for incoming coordinates and the lower and upper fix point for outgoing coordinates. the resulting functor maps incoming coordinates in the range of [in_low,in_high] to coordinates in the range of [out_low,out_high]
vspline::domain_type< coordinate_type, _vsize > vspline::domain | ( | const coordinate_type & | in_low, |
const coordinate_type & | in_high, | ||
const spline_type & | bspl | ||
) |
factory function to create a domain_type type object from the desired lower and upper reference point for incoming coordinates and a vspline::bspline object providing the lower and upper reference for outgoing coordinates the resulting functor maps incoming coordinates in the range of [ in_low , in_high ] to coordinates in the range of [ bspl.lower_limit() , bspl.upper_limit() ]
vspline::domain_type< coordinate_type, _vsize > vspline::domain | ( | const spline_type & | bspl, |
const coordinate_type & | out_low, | ||
const coordinate_type & | out_high | ||
) |
factory function to create a domain_type type object from a vspline::bspline object providing the lower and upper reference for incomoing coordinates and the desired lower and upper reference point for outgoing coordinates the resulting functor maps incoming coordinates in the range of [ bspl.lower_limit() , bspl.upper_limit() ] to coordinates in the range of [ out_low , out_high ]
vspline::domain_type< coordinate_type, _vsize > vspline::domain | ( | const spline_type & | bspl_in, |
const spline_type & | bspl_out | ||
) |
factory function to create a domain_type type object from a vspline::bspline object providing the lower and upper reference for incomoing coordinates and a vspline::bspline object providing the lower and upper reference for outgoing coordinates the resulting functor maps incoming coordinates in the range of [ bspl_in.lower_limit() , bspl_in.upper_limit() ] to outgoing coordinates in the range of [ bspl_out.lower_limit() , bspl_out.upper_limit() ]
std::function< void(const IN &, OUT &) > vspline::eval_wrap | ( | std::function< OUT(const IN &) > | f | ) |
eval_wrap is a helper function template, wrapping an 'ordinary' function which returns some value given some input in a void function taking input as const reference and writing output to a reference, which is the signature used for evaluation in vspline::unary_functors.
Definition at line 337 of file unary_functor.h.
void vspline::even_split | ( | rc_t | v, |
ic_t & | iv, | ||
rc_t & | fv | ||
) |
void vspline::even_split | ( | vigra::TinyVector< rc_t, N > | v, |
vigra::TinyVector< ic_t, N > & | iv, | ||
vigra::TinyVector< rc_t, N > & | fv | ||
) |
bool vspline::fetch_ascending | ( | vspline::atomic< index_t > & | source, |
const index_t & | total, | ||
index_t & | index | ||
) |
fetch_ascending counts up from zero to total-1, which is more efficient if the indexes are used to address memory. This is due to the side effects of accessing memory: if memory is accessed at an address x, the OS will typically fetch a chunk of data starting at or shortly before x. If the next fetch requires data just after x, there is a good chance that they are already in cache.
Definition at line 284 of file multithread.h.
bool vspline::fetch_descending | ( | vspline::atomic< index_t > & | source, |
index_t & | index | ||
) |
fetch_descending fetches the next index from an atomic, counting down. Indexes will range from one less than the value the atomic was initialized with, down to zero. If all indexes were distributed already, false is returned, true otherwise. Like the other fetch_XXX routines, the return of a boolean makes these functions good candidates to be used as conditional in a loop.
Definition at line 264 of file multithread.h.
bool vspline::fetch_range_ascending | ( | vspline::atomic< index_t > & | source, |
const index_t & | count, | ||
const index_t & | total, | ||
index_t & | low, | ||
index_t & | high | ||
) |
fetch_range_ascending also uses an atomic initialized to the total number of indexes to be distributed, but the successive ranges are handed out in ascending order, which is more efficient if the indexes are used to address memory.
Definition at line 336 of file multithread.h.
bool vspline::fetch_range_descending | ( | vspline::atomic< index_t > & | source, |
const index_t & | count, | ||
index_t & | low, | ||
index_t & | high | ||
) |
fetch_range_ascending fetches the beginning and end of a range of indexes (in iterator semantic, low <= i < high) from a vspline::atomic which has been initialized with the total number of indexes that are to be processed. If the vspline::atomic, when accessed, already holds a value of or below zero, fetch_index_range returns false and leaves low and high unchanged. Otherwise it returns true and low and high will be set to the values gleaned from the atomic, raising 'low 'to zero if it would come out below. this function (and the next) enable calling code to process batches of adjacent indexes without any index artistry: the code is perfectly general and simple, and the use of the atomic and fetch_sub garantees that each fetch provides a distinct batch of indexes.
Definition at line 311 of file multithread.h.
void vspline::filter | ( | const vigra::MultiArrayView< D, in_type > & | input, |
vigra::MultiArrayView< D, out_type > & | output, | ||
types ... | args | ||
) |
vspline::filter is the common entry point for filter operations in vspline. This routine does not yet do any processing, it's purpose is to convert it's arguments to 'canonical' format and then call the actual filter code in namespace detail. It also determines the type used for arithmetic operations. The type specification for input and output assures that only arrays with the same dimensionality are accepted, and a static assertion makes sure the number of channels match. canonical form means that input and output value type are either fundamental (for single-channel data) or TinyVectors of a fundamental data type. This way, the input and output is presented in a neutral form, ignoring all specifics of T1 and T2, which may be TinyVectors, complex, RGBValue etc.
void vspline::forward_backward_recursive_filter | ( | const vigra::MultiArrayView< dimension, in_value_type > & | input, |
vigra::MultiArrayView< dimension, out_value_type > & | output, | ||
vigra::TinyVector< bc_code, static_cast< int >(dimension) > | bcv, | ||
std::vector< vspline::xlf_type > | pv, | ||
xlf_type | tolerance = -1 , |
||
xlf_type | boost = xlf_type ( 1 ) , |
||
int | axis = -1 , |
||
int | njobs = default_njobs |
||
) |
forward_backward_recursive_filter applies one or more pairs of simple recursive filters to the input. Each pair consists of:
Definition at line 132 of file general_filter.h.
real_type vspline::gaussian_bspline_basis_approximation | ( | real_type | x, |
int | degree | ||
) |
Gaussian approximation to B-spline basis function. This routine approximates the basis function of degree spline_degree for real x. I checked for all degrees up to 45. The partition of unity quality of the resulting reconstruction filter is okay for larger degrees, the cumulated error over the covered interval is quite low. Still, as the basis function is never actually evaluated in vspline (whenever it's needed, it is needed at n * 1/2 and we have precomputed values for that) there is not much point in having this function around. I leave the code in for now.
void vspline::gen_grid_eval | ( | grid_spec< functor_type::dim_in, typename functor_type::in_ele_type > & | grid, |
const functor_type & | functor, | ||
vigra::MultiArrayView< functor_type::dim_in, typename functor_type::out_type > & | result, | ||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
for backward compatibility, we provide 'gen_grid_eval'. This is deprecated and will eventually go away, new code should use 'transform' above.
Definition at line 1376 of file transform.h.
void vspline::get_kernel | ( | const int & | degree, |
vigra::MultiArrayView< 1, target_type > & | kernel, | ||
const bool & | odd = true |
||
) |
this function deposits the reconstruction kernel in the array 'kernel'. This kernel can be used to convolve a set of coefficients, to obtain the original signal. This is a convenience function which merely picks the right values from the precomputed values in precomputed_basis_function_values. if 'odd' is passed false, the result is an even kernel. This kernel can't be used for reconstruction (.5 phase shift), but it's handy to get values half a unit step from the knot points.
void vspline::grid_eval | ( | grid_spec< functor_type::dim_in, typename functor_type::in_ele_type > & | grid, |
const functor_type & | functor, | ||
vigra::MultiArrayView< functor_type::dim_in, typename functor_type::out_type > & | result, | ||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
for backward compatibility, we provide 'grid_eval'. This is deprecated and will eventually go away, new code should use 'transform' above.
Definition at line 1360 of file transform.h.
vspline::grok_type< typename grokkee_type::in_type, typename grokkee_type::out_type, grokkee_type::vsize > vspline::grok | ( | const grokkee_type & | grokkee | ) |
grok() is the corresponding factory function, wrapping grokkee in a vspline::grok_type.
Definition at line 871 of file unary_functor.h.
vspline::grok_type< bspl_coordinate_type< spline_type, rc_type >, result_type, _vsize > vspline::make_evaluator | ( | const spline_type & | bspl, |
vigra::TinyVector< int, spline_type::dimension > | dspec = vigra::TinyVector<int,spline_type::dimension> ( 0 ) , |
||
int | shift = 0 |
||
) |
make_evaluator is a factory function, producing a functor which provides access to an evaluator object. Evaluation using the resulting object is not intrinsically safe, it's the user's responsibility not to pass coordinates which are outside the spline's defined range. If you need safe access, see 'make_safe_evaluator' below. 'Not safe' in this context means that evaluation at out-of-bounds locations may result in a memory fault or produce wrong or undefined results. Note that vspline's bspline objects are set up as to allow evaluation at the lower and upper limit of the spline and will also tolerate values 'just outside' the bounds to guard against quantization errors. see vspline::bspline for details.
The evaluator will be specialized to the spline's degree: degree 0 splines will be evaluated with nearest neighbour, degree 1 splines with linear interpolation, all other splines will use general b-spline evaluation.
This function returns the evaluator wrapped in an object which hides it's type. This object only 'knows' what coordinates it can take and what values it will produce. The extra level of indirection may cost a bit of performance, but having a common type simplifies handling. The wrapped evaluator also provides operator().
So, if you have a vspline::bspline object 'bspl', you can use this factory function like this:
auto ev = make_evaluator ( bspl ) ; typedef typename decltype(ev)::in_type coordinate_type ; coordinate_type c ; auto result = ev ( c ) ;
make_evaluator requires one template argument: spline_type, the type of the vspline::bspline object you want to have evaluated. Optionally, you can specify the elementary type for coordinates (use either float or double) and the vectorization width. The latter will only have an effect if vectorization is used and the spline's data type can be vectorized. Per default, the vectorization width will be inferred from the spline's value_type by querying vspline::vector_traits, which tries to provide a 'good' choice. Note that a specific evaluator will only be capable of processing vectorized coordinates of precisely the _vsize it has been created with. A recent addition to the template argument list is 'math_ele_type', which allows you to pick a different type for internal processing than the default. The default is a real type 'appropriate' to the data in the spline.
Note that the object created by this factory function will only work properly if you evaluate coordinates of the specific 'rc_type' used. If you create it with the default rc_type, which is float (and usually sufficiently precise for a coordinate), you can't evaluate double precision coordinates with it.
On top of the bspline object, you can optionally pass a derivative specification and a shift value, which are simply passed through to the evaluator's constructor, see there for the meaning of these optional parameters.
While the declaration of this function looks frightfully complex, using it is simple: in most cases it's simply
auto ev = make_evaluator ( bspl ) ;
For an explanation of the template arguments, please see make_safe_evaluator() below, which takes the same template args.
vspline::grok_type< bspl_coordinate_type< spline_type, rc_type >, result_type, _vsize > vspline::make_safe_evaluator | ( | const spline_type & | bspl, |
vigra::TinyVector< int, spline_type::dimension > | dspec = vigra::TinyVector<int,spline_type::dimension> ( 0 ) , |
||
int | shift = 0 |
||
) |
make_safe_evaluator is a factory function, producing a functor which provides safe access to an evaluator object. This functor will map incoming coordinates into the spline's defined range, as given by the spline with it's lower_limit and upper_limit methods, honoring the bspline objects's boundary conditions. So if, for example, the spline is periodic, all incoming coordinates are valid and will be mapped to the first period. Note the use of lower_limit and upper_limit. These values also depend on the spline's boundary conditions, please see class vspline::bspline for details. If there is no way to meaningfully fold a coordinate into the defined range, the coordinate is clamped to the nearest limit.
The evaluator will be specialized to the spline's degree: degree 0 splines will be evaluated with nearest neighbour, degree 1 splines with linear interpolation, all other splines will use general b-spline evaluation.
This function returns the functor wrapped in an object which hides it's type. This object only 'knows' what coordinates it can take and what values it will produce. The extra level of indirection costs a bit of performance, but having a common type simplifies handling: the type returned by this function only depends on the spline's data type, the coordinate type and the vectorization width.
Also note that the object created by this factory function will only work properly if you evaluate coordinates of the specific 'rc_type' used. If you create it with the default rc_type, which is float (and usually sufficiently precise for a coordinate), you can't evaluate double precision coordinates with it.
On top of the bspline object, you can optionally pass a derivative specification and a shift value, which are simply passed through to the evlauator's constructor, see there for the meaning of these optional parameters.
While the declaration of this function looks frightfully complex, using it is simple: in most cases it's simply
auto ev = make_safe_evaluator ( bspl ) ;
The first template argument, spline_type, is the type of a vspline::bspline object. This template argument has no default, since it determines the dimensionality and the coefficient type. But since the first argument to this factory function is of this type, spline_type can be fixed via ATD, so it can be omitted.
The second template argument, rc_type, can be used to pick a different elementary type for the coordinates the evaluator will accept. In most cases the default, float, will be sufficient.
The next template argument, _vsize, fixes the vectorization width. Per default, this will be what vspline deems appropriate for the spline's coefficient type.
math_ele_type can be used to specify a different fundamental type to be used for arithemtic operations during evaluation. The default used here is a real type of at least the precision of the coordinates or the spline's coefficients, but you may want to raise precision here, for example by passing 'double' while your data are all float.
Finally you can specify a result type. Per default the result will be of the same type as the spline's coefficients, but you may want a different value here - a typical example would be a spline with integral coefficients, where you might prefer to get the result in, say, float to avoid quantization errors on the conversion from the 'internal' result (which is in math_type) to the output.
vspline::map_functor< nd_rc_type, _vsize, gate_types... > vspline::mapper | ( | gate_types ... | args | ) |
factory function to create a mapper type functor given a set of gate_type objects. Please see vspline::make_safe_evaluator for code to automatically create a mapper object suitable for a specific vspline::bspline.
vspline::mirror_gate< rc_type, _vsize > vspline::mirror | ( | rc_type | lower, |
rc_type | upper | ||
) |
factory function to create a mirror_gate type functor given a lower and upper limit for the allowed range.
int vspline::multithread | ( | std::function< void() > | payload, |
std::size_t | nr_workers = default_njobs |
||
) |
multithread uses a thread pool of worker threads to perform a multithreaded operation. It receives a functor (a single-threaded function used for all individual tasks), and, optionally, the desired number of worker instances to be used. These tasks are wrapped with a wrapper which takes care of signalling when the last task has completed.
This, in a way, is the purest implementation of 'multithread': this implementation is not involved with any property of the jobs at hand, it merely forwards the arguments for the payload function to a specified number of workers. The code invoking 'multithread' has to set up whatever scheme it deems appropriate to convey, via the arguments, what the worker threads should do. with this mechanism in place, the calling code has the very efficient option of setting up a vspline::atomic holding some sort of 'joblet number' and passing a reference to this atomic to the payload code. The worker(s) executing payload code all get equal load until the job numbers are exhausted, which is when they terminate one by one. Since there is no inter-thread communication during the active phase, there is no signalling overhead at all, which allows fine granularity. The fine granularity ensures little tail-end idling (when the caller has to wait for the last worker to finish) and also makes it possible to choose some sort of partitioning which insinuates itself from the structure of the data at hand, rather than some preconceived parcel of the total job. When there are many more joblet numbers than workers, intermittent inactivity of a worker simply makes it consume fewer job numbers, rather than delaying it on a way to a preset goal.
Another effect is that, if job numbers relate to memory worked on, (think of lines of an image) - all activity is focussed in a narrow band of the memory, because the currently processed job numbers are all usually in sequence (unless the OS throws a spanner in the works by halting threads). This may or may not help - in some situations having several threads access adjacent memory locations may make it harder for the system to synchronize access. But payload code is free to use any interpretation of job numbers anyway, so that's an issue on the payload side.
Since the number of threads in the pool is static, requesting more workers than there are threads is futile (but still works). Requesting fewer may be useful to have some task 'on the back burner' while some critical task receives more workers to complete faster.
As I've pointed out in thread_pool.h, it seems beneficial (at least for vspline) to have a good deal more threads than physical cores. See there for reasons why this may be so.
Last but not least: the code is very simple :)
Why is multithread a template? So that several TUs can #include vspline without linker errors.
Definition at line 412 of file multithread.h.
Definition at line 928 of file simd_type.h.
void vspline::odd_split | ( | rc_t | v, |
ic_t & | iv, | ||
rc_t & | fv | ||
) |
coordinates are split into an integral part and a remainder. this is used for weight generation, and also for calculating the basis function value. The split is done differently for odd and even splines.
Note how the initial call to std::floor produces a real type, which is used to subtract from 'v', yielding the remainder in 'fv'. Only after having used this real representation of the integral part, it is cast to an integral type by assigning it to 'iv'. This is the most efficient route, better than producing an integral-typed integral part directly and subtracting that from 'v', which would require another conversion. Technically, one might consider this 'split' as a remainder division by 1.
void vspline::odd_split | ( | vigra::TinyVector< rc_t, N > | v, |
vigra::TinyVector< ic_t, N > & | iv, | ||
vigra::TinyVector< rc_t, N > & | fv | ||
) |
vspline::chain_type< T1, T2 > vspline::operator+ | ( | const T1 & | t1, |
const T2 & | t2 | ||
) |
using operator overloading, we can exploit operator+'s semantics to chain several unary functors. We need to specifically enable this for types derived from unary_functor_tag to avoid a catch-all situation.
Definition at line 613 of file unary_functor.h.
vspline::pass_gate< rc_type, _vsize > vspline::pass | ( | ) |
vspline::periodic_gate< rc_type, _vsize > vspline::periodic | ( | rc_type | lower, |
rc_type | upper | ||
) |
factory function to create a periodic_gate type functor given a lower and upper limit for the allowed range.
void vspline::prefilter | ( | const vigra::MultiArrayView< dimension, in_value_type > & | input, |
vigra::MultiArrayView< dimension, out_value_type > & | output, | ||
vigra::TinyVector< bc_code, static_cast< int >(dimension) > | bcv, | ||
int | degree, | ||
xlf_type | tolerance = std::numeric_limits < math_ele_type > :: epsilon() , |
||
xlf_type | boost = xlf_type ( 1 ) , |
||
int | njobs = default_njobs |
||
) |
'prefilter' handles b-spline prefiltering for the whole range of acceptable input and output. It combines two bodies of code to achieve this goal:
Note that vsize , the vectorization width, can be passed explicitly. If Vc is in use and math_ele_type can be used with hardware vectorization, the arithmetic will be done with Vc::SimdArrays of the given size. Otherwise 'goading' will be used: the data are presented in TinyVectors of vsize math_ele_type, hoping that the compiler may autovectorize the operation.
Definition at line 1082 of file prefilter.h.
void vspline::reduce | ( | const sink_functor_type & | functor, |
const vigra::MultiArrayView< dimension, typename sink_functor_type::in_type > & | input, | ||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
reduce function processing by value. sink_functor_type is a functor taking values (which are loaded from an array). The reduce function, by delegation to wielding::value_reduce, performs the feeding part by parcelling the data, the functor receives either single in_type or vectorized data. The functor has to be a construct capable of accumulating data with several threads in parallel. A typical example would be a functor cumulating to some internal structure and 'offloading' to a mutex-protected pooling target during it's destructor. The wielding code creates a copy of the functor which is passed here for each worker thread and guarantees that this copy's destructor is executed before the call to reduce returns, so the user can be assured of having collected the partial results of all worker threads in the pooling target. The functor needs a copy constructor to communicate the pooling target to the copies made for the individual threads, see the index-based reduce function below this one for an illustration of the simple and elegant mechanism. Note that sink_functor will normally be a type inheriting from vspline::sink_functor, a close relative of vspline::unary_functor and also defined in unary_functor.h. This inheritance does not produce any functionality, just the standard vspline functor argument type system, which this code relies on. Of course the sink functor can be coded 'manually' as well - inheriting form vspline::sink_functor is merely convenient.
Definition at line 286 of file transform.h.
void vspline::reduce | ( | const sink_functor_type & | functor, |
typename std::conditional<(sink_functor_type::dim_in==1), int, vigra::TinyVector< int, sink_functor_type::dim_in > > ::type | shape, | ||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
while the function above fills an array by calling a functor with the target coordinates, this function accumulates the results by means of a sink_functor. The effect is as if an index-based transform was executed first, filling an intermediate array, followed by a call to reduce taking the intermediate array as input. The combination of both operations 'cuts out' the intermediate. Since this function does not have a target array, the shape has to be communicated explicitly. The functor will be invoked for every possible coordinate within the shape's scope. This is done by the function 'index_reduce' in wielding.h the type of the argument 'shape' is derived from the sink_functor's 'dim_in' value, to make the code 'respect the singular': if the data are 1D, reduce will accept a single integer as 'shape', whereas nD data require a TinyVector of int. For suitable sink_functors, see the remarks and example given above in the comments for the array-base overload of 'reduce'.
Definition at line 480 of file transform.h.
vspline::reject_gate< rc_type, _vsize > vspline::reject | ( | rc_type | lower, |
rc_type | upper | ||
) |
factory function to create a reject_gate type functor given a lower and upper limit for the allowed range.
void vspline::remap | ( | const vigra::MultiArrayView< cf_dimension, original_type > & | input, |
const vigra::MultiArrayView< trg_dimension, coordinate_type > & | coordinates, | ||
vigra::MultiArrayView< trg_dimension, result_type > & | output, | ||
bcv_type< bcv_dimension > | bcv = bcv_type < bcv_dimension > ( MIRROR ) , |
||
int | degree = 3 , |
||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
Implementation of 'classic' remap, which directly takes an array of values and remaps it, internally creating a b-spline of given order just for the purpose. This is used for one-shot remaps where the spline isn't reused, and specific to b-splines, since the functor used is a b-spline evaluator. The spline defaults to a cubic b-spline with mirroring on the bounds.
So here we have the 'classic' remap, where the input array holds coordinates and the functor used is actually an interpolator. Since this is merely a special case of using transform(), we delegate to transform() once we have the evaluator.
The template arguments are chosen to allow the user to call 'remap' without template arguments; the template arguments can be found by ATD by looking at the MultiArrayViews passed in.
optionally, remap takes a set of boundary condition values and a spline degree, to allow creation of splines for specific use cases beyond the default. I refrain from extending the argument list further; user code with more specific requirements will have to create an evaluator and use 'transform'.
Note that remap can be called without template arguments, the types will be inferred by ATD from the arguments passed in.
See the remarks on the two trailing parameters given with the two-array overload of transform.
Definition at line 600 of file transform.h.
void vspline::restore | ( | const vspline::bspline< value_type, dimension > & | bspl, |
vigra::MultiArrayView< dimension, value_type > & | target | ||
) |
restore restores the original data from the b-spline coefficients. This is done efficiently using a separable convolution, the kernel is simply a unit-spaced sampling of the basis function. Since the filter uses internal buffering, using this routine in-place is safe - meaning that 'target' may be bspl.core itself. math_type, the data type for performing the actual maths on the buffered data, and the type the data are converted to when they are placed into the buffer, can be chosen, but I could not detect any real benefits from using anything but the default, which is to leave the data in their 'native' type.
an alternative way to restore is running an index-based transform with an evaluator for the spline. This is much less efficient, but the effect is the same:
auto ev = vspline::make_evaluator ( bspl ) ; vspline::transform ( ev , target ) ;
Note that vsize, the vectorization width, can be passed explicitly. If Vc is in use and math_ele_type can be used with hardware vectorization, the arithmetic will be done with Vc::SimdArrays of the given size. Otherwise 'goading' will be used: the data are presented in simd_type of vsize math_ele_type, hoping that the compiler may autovectorize the operation.
'math_ele_type', the type used for arithmetic inside the filter, defaults to the vigra RealPromote type of value_type's elementary. This ensures appropriate treatment of integral-valued splines.
Definition at line 1429 of file transform.h.
void vspline::restore | ( | vspline::bspline< value_type, dimension > & | bspl | ) |
overload of 'restore' writing the result of the operation back to the array which is passed in. This looks like an in-place operation, but the data are in fact moved to a buffer stripe by stripe, then the arithmetic is done on the buffer and finally the buffer is written back. This is repeated for each dimension of the array.
Definition at line 1524 of file transform.h.
void vspline::transform | ( | const unary_functor_type & | functor, |
const vigra::MultiArrayView< dimension, typename unary_functor_type::in_type > & | input, | ||
vigra::MultiArrayView< dimension, typename unary_functor_type::out_type > & | output, | ||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
implementation of two-array transform using wielding::coupled_wield.
'array-based' transform takes two template arguments:
this overload of transform takes three parameters:
transform can be used without template arguments, they will be inferred by ATD from the arguments.
trailing the argument list, all transform overloads have two further parameters with default values, which are only used for special purposes:
Definition at line 211 of file transform.h.
void vspline::transform | ( | const unary_functor_type & | functor, |
vigra::MultiArrayView< unary_functor_type::dim_in, typename unary_functor_type::out_type > & | output, | ||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
reduce function processing by coordinate. sink_functor_type is a functor receiving integral coordinates, with no fixed meaning, but usually pertaining to an array. The reduce function feeds all coordinates encompassed by the shape, either in vectorized or single form. Here's an example of a functor for this function:
implementation of index-based transform using wielding::index_wield
this overload of transform() is very similar to the first one, but instead of picking input from an array, it feeds the discrete coordinates of the successive places data should be rendered to to the unary_functor_type object.
This sounds complicated, but is really quite simple. Let's assume you have a 2X3 output array to fill with data. When this array is passed to transform, the functor will be called with every coordinate pair in turn, and the result the functor produces is written to the output array. So for the example given, with 'ev' being the functor, we have this set of operations:
output [ ( 0 , 0 ) ] = ev ( ( 0 , 0 ) ) ;
output [ ( 1 , 0 ) ] = ev ( ( 1 , 0 ) ) ;
output [ ( 2 , 0 ) ] = ev ( ( 2 , 0 ) ) ;
output [ ( 0 , 1 ) ] = ev ( ( 0 , 1 ) ) ;
output [ ( 1 , 1 ) ] = ev ( ( 1 , 1 ) ) ;
output [ ( 2 , 1 ) ] = ev ( ( 2 , 1 ) ) ;
this transform overload takes one template argument:
this transform overload takes two parameters:
Please note that vspline holds with vigra's coordinate handling convention, which puts the fastest-changing index first. In a 2D, image processing, context, this is the column index, or the x coordinate. C and C++ do instead put this index last when using multidimensional array access code.
transform can be used without template arguments, they will be inferred by ATD from the arguments.
See the remarks on the two trailing parameters given with the two-array overload of transform.
Definition at line 433 of file transform.h.
void vspline::transform | ( | grid_spec< 1, typename functor_type::in_ele_type > & | grid, |
const functor_type & | functor, | ||
vigra::MultiArrayView< 1, typename functor_type::out_type > & | result, | ||
int | njobs, | ||
vspline::atomic< bool > * | p_cancel, | ||
std::true_type | |||
) |
overload of grid-based transform for 1D grids. Obviously, this is a case for using 'transform' directly taking the single set of per-axis grid values as input, and this overload is here only for completeness' sake - 'normal' user code will use an array-based transform straight away.
Definition at line 1323 of file transform.h.
void vspline::transform | ( | grid_spec< functor_type::dim_in, typename functor_type::in_ele_type > & | grid, |
const functor_type & | functor, | ||
vigra::MultiArrayView< functor_type::dim_in, typename functor_type::out_type > & | result, | ||
int | njobs, | ||
vspline::atomic< bool > * | p_cancel, | ||
std::false_type | |||
) |
This overload of 'transform' takes a vspline::grid_spec as it's first parameter. This object defines a grid where each grid point is made up from components taken from the per-axis grid values. Let's say the grid_spec holds.
{ 1 , 2 } and { 3 , 4 },
then we get grid points
( 1 , 3 ) , ( 2 , 3 ) , ( 1 , 4 ) , ( 2 , 4 )
The grid points are taken as input for the functor, and the results are stored in 'result'. So with the simple example above, result would have to be a 2D array and after the transform it would contain
{ { f ( ( 1 , 3 ) ) , f ( ( 2 , 3 ) ) } ,
{ f ( ( 1 , 4 ) ) , f ( ( 2 , 4 ) ) } }
The grid specification obviously takes up less space than the whole grid would occupy, and if the functor is a b-spline evaluator, additional performance gain is possible by precalculating the evaluation weights. The distinction is made here by dispatching to two specific overloads in namespace detail, which provide the implementation. Note that only functors of type vspline::evaluator will qualify as b-spline evaluators. The functors produced by the factory functions 'make_evaluator' and 'make_safe_evaluator' will not qualify; they are of vspline::grok_type and even though they have a b-spline evaluator 'somewhere inside' this can't be accessed due to the type erasure.
See the remarks on the two trailing parameters given with the two-array overload of transform.
Definition at line 1294 of file transform.h.
void vspline::transform | ( | grid_spec< functor_type::dim_in, typename functor_type::in_ele_type > & | grid, |
const functor_type & | functor, | ||
vigra::MultiArrayView< functor_type::dim_in, typename functor_type::out_type > & | result, | ||
int | njobs = vspline::default_njobs , |
||
vspline::atomic< bool > * | p_cancel = 0 |
||
) |
Definition at line 1341 of file transform.h.
rc_v vspline::v_fmod | ( | rc_v | lhs, |
const typename rc_v::value_type & | rhs | ||
) |
const std::string vspline::bc_name[] |
const int vspline::default_njobs = vspline_threadpool::default_njobs |
Definition at line 220 of file multithread.h.
const int vspline::ncores = vspline_threadpool::ncores |
Definition at line 219 of file multithread.h.