Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
c2_function.hh
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 /**
27  * \file electromagnetic/TestEm7/include/c2_function.hh
28  * \brief Provides the headers for the general c2_function algebra which supports
29  * fast, flexible operations on piecewise-twice-differentiable functions
30  *
31  * \author Created by R. A. Weller and Marcus H. Mendenhall on 7/9/05.
32  * \author Copyright 2005 __Vanderbilt University__. All rights reserved.
33  *
34  * \version c2_function.hh 490 2012-04-10 19:05:40Z marcus
35  * \see \ref c2_factory "Factory Functions" for information on constructing things in here
36  */
37 
38 //
39 // $Id: c2_function.hh 66587 2012-12-21 11:06:44Z ihrivnac $
40 
41 #ifndef __has_c2_function_hh
42 #define __has_c2_function_hh 1
43 
44 // MSVC does not automatically define numerical constants such as M_PI without this.
45 // this came from the msdn website, so it should be right...
46 #ifdef _MSC_VER
47 #define _USE_MATH_DEFINES
48 #define c2_isnan _isnan
49 #define c2_isfinite _finite
50 #else
51 #define c2_isnan std::isnan
52 #define c2_isfinite std::isfinite
53 #endif
54 
55 #include <cmath>
56 #include <vector>
57 #include <utility>
58 #include <string>
59 #include <stdexcept>
60 #include <typeinfo>
61 #include <sstream>
62 #include <limits> // fails under gcc-4.3 without this here, was ok in c2_function.cc before
63 
64 /// \brief the exception class for c2_function operations.
65 class c2_exception : public std::exception {
66 public:
67  /// \brief construct the exception with an error message
68  /// \param msgcode the message
69  c2_exception(const char msgcode[]) : info(msgcode) { }
70  virtual ~c2_exception() throw() { }
71  /** Returns a C-style character string describing the general cause
72  * of the current error. */
73  virtual const char* what() const throw() { return info.c_str(); }
74 private:
75  std::string info;
76 };
77 
78 // put these forward references here, and with a bogus typename to make swig happy.
79 template <typename float_type> class c2_composed_function_p;
80 template <typename float_type> class c2_sum_p;
81 template <typename float_type> class c2_diff_p;
82 template <typename float_type> class c2_product_p;
83 template <typename float_type> class c2_ratio_p;
84 template <typename float_type> class c2_piecewise_function_p;
85 template <typename float_type> class c2_quadratic_p;
86 template <typename float_type> class c2_ptr;
87 /**
88  \defgroup abstract_classes Abstract Classes
89  \defgroup arithmetic_functions Arithmetic Functions
90  \defgroup math_functions Mathemetical Functions
91  \defgroup parametric_functions Parametric Families of Functions
92  \defgroup interpolators Interpolating Functions
93  \defgroup containers Functions which are containers for, or functions of, other functions
94  \defgroup factories Factory classes which reduce silly template typing
95  \defgroup transforms Classes which provide coordinate system transformations, wih derivatives
96 */
97 
98 /// \brief structure used to hold evaluated function data at a point.
99 ///
100 /// Contains all the information for the function at one point.
101 template <typename float_type> class c2_fblock
102 {
103 public:
104  /// \brief the abscissa
105  float_type x;
106  /// \brief the value of the function at \a x
107  float_type y;
108  /// \brief the derivative at \a x
109  float_type yp;
110  /// \brief the second derivative at \a x
111  float_type ypp;
112  /// flag, filled in by c2_function::fill_fblock(), indicating the derivative is NaN of Inf
113  bool ypbad;
114  /// flag, filled in by c2_function::fill_fblock(), indicating the second derivative is NaN of Inf
115  bool yppbad;
116 };
117 
118 /**
119  \brief the parent class for all c2_functions.
120  \ingroup abstract_classes
121  c2_functions know their value, first, and second derivative at almost every point.
122  They can be efficiently combined with binary operators, via c2_binary_function,
123  composed via c2_composed_function_,
124  have their roots found via find_root(),
125  and be adaptively integrated via partial_integrals() or integral().
126  They also can carry information with them about how to find 'interesting' points on the function.
127  This information is set with set_sampling_grid() and extracted with get_sampling_grid().
128 
129  Particularly important subclasses are the interpolating functions classes,
130  interpolating_function , lin_log_interpolating_function, log_lin_interpolating_function,
131  log_log_interpolating_function, and arrhenius_interpolating_function,
132  as well as the template functions
133  inverse_integrated_density_function().
134 
135  For a discussion of memory management, see \ref memory_management
136 
137  */
138 template <typename float_type=double> class c2_function {
139 public:
140  /// \brief get versioning information for the header file
141  /// \return the CVS Id string
142  const std::string cvs_header_vers() const { return
143  "c2_function.hh 490 2012-04-10 19:05:40Z marcus ";
144  }
145 
146  /// \brief get versioning information for the source file
147  /// \return the CVS Id string
148  const std::string cvs_file_vers() const ;
149 
150 public:
151  /// \brief destructor
152  virtual ~c2_function() {
154  if(root_info) delete root_info;
155  if(owner_count) {
156  std::ostringstream outstr;
157  outstr << "attempt to delete an object with non-zero ownership in class ";
158  outstr << typeid(*this).name() << std::endl;
159  throw c2_exception(outstr.str().c_str());
160  }
161  }
162 
163  /// \brief get the value and derivatives.
164  ///
165  /// There is required checking for null pointers on the derivatives,
166  /// and most implementations should operate faster if derivatives are not needed.
167  /// \param[in] x the point at which to evaluate the function
168  /// \param[out] yprime the first derivative (if pointer is non-null)
169  /// \param[out] yprime2 the second derivative (if pointer is non-null)
170  /// \return the value of the function
171  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception) =0 ; // { return 0; };
172 
173  /// \brief evaluate the function in the classic way, ignoring derivatives.
174  /// \param x the point at which to evaluate
175  /// \return the value of the function
176  inline float_type operator () (float_type x) const throw(c2_exception)
177  { return value_with_derivatives(x, (float_type *)0, (float_type *)0); }
178 
179  /// \brief get the value and derivatives.
180  ///
181  /// \param[in] x the point at which to evaluate the function
182  /// \param[out] yprime the first derivative (if pointer is non-null)
183  /// \param[out] yprime2 the second derivative (if pointer is non-null)
184  /// \return the value of the function
185  inline float_type operator () (float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
186  { return value_with_derivatives(x, yprime, yprime2); }
187 
188  /// \brief solve f(x)==value very efficiently, with explicit knowledge of derivatives of the function
189  ///
190  /// find_root solves by iterated inverse quadratic extrapolation for a solution to f(x)=y. It
191  /// includes checks against bad convergence, so it should never be able to fail. Unlike typical
192  /// secant method or fancier Brent's method finders, this does not depend in any strong wasy on the
193  /// brackets, unless the finder has to resort to successive approximations to close in on a root.
194  /// Often, it is possible to make the brackets equal to the domain of the function, if there is
195  /// any clue as to where the root lies, as given by the parameter \a start.
196  /// \param lower_bracket the lower bound for the search
197  /// \param upper_bracket the upper bound for the search. Function sign must be
198  /// opposite to that at \a lower_bracket
199  /// \param start starting value for the search
200  /// \param value the value of the function being sought (solves f(x) = \a value)
201  /// \param[out] error If pointer is zero, errors raise exception. Otherwise, returns error here.
202  /// \param[out] final_yprime If pointer is not zero, return derivative of function at root
203  /// \param[out] final_yprime2 If pointer is not zero, return second derivative of function at root
204  /// \return the position of the root.
205  /// \see \ref rootfinder_subsec "Root finding sample"
206  float_type find_root(float_type lower_bracket, float_type upper_bracket, float_type start,
207  float_type value, int *error=0,
208  float_type *final_yprime=0, float_type *final_yprime2=0 ) const throw(c2_exception) ; // solve f(x)=value
209 
210  /// \brief for points in xgrid, adaptively return Integral[f(x),{x,xgrid[i],xgrid[i+1]}] and return in vector, along with sum
211  ///
212  /// partial_integrals uses a method with an error O(dx**10) with full information from the derivatives,
213  /// and falls back to lower order methods if informed of incomplete derivatives.
214  /// It uses exact midpoint splitting of the intervals for recursion, resulting in no recomputation of the function
215  /// during recursive descent at previously computed points.
216  /// \param xgrid points between which to evaluate definite integrals.
217  /// \param partials if non-NULL, a vector in which to receive the partial integrals.
218  /// It will automatically be sized apprpropriately, if provided, to contain \a n - 1 elements where \a n is the length of \a xgrid
219  /// \param abs_tol the absolute error bound for each segment
220  /// \param rel_tol the fractional error bound for each segment.
221  /// If the error is smaller than either the relative or absolute tolerance, the integration step is finished.
222  /// \param derivs number of derivatives to trust, which sets the order of the integrator. The order
223  /// is 3*\a derivs + 4. \a derivs can be 0, 1, or 2.
224  /// \param adapt if true, use recursive adaptation, otherwise do simple evaluation on the grid provided
225  /// with no error checking.
226  /// \param extrapolate if true, use simple Richardson extrapolation on the final 2 steps to reduce the error.
227  /// \return sum of partial integrals, which is the definite integral from the first value in \a xgrid to the last.
228  float_type partial_integrals(std::vector<float_type> xgrid, std::vector<float_type> *partials = 0,
229  float_type abs_tol=1e-12, float_type rel_tol=1e-12, int derivs=2, bool adapt=true, bool extrapolate=true)
230  const throw(c2_exception);
231 
232  /// \brief a fully-automated integrator which uses the information provided by the get_sampling_grid() function
233  /// to figure out what to do.
234  ///
235  /// It returns the integral of the function over the domain requested
236  /// with error tolerances as specified. It is just a front-end to partial_integrals()
237  ///
238  /// \param amin lower bound of the domain for integration
239  /// \param amax upper bound of the domain for integration
240  /// \param partials if non-NULL, a vector in which to receive the partial integrals.
241  /// It will automatically be sized appropriately, if provided, to contain \a n - 1 elements where \a n is the length of \a xgrid
242  /// \param abs_tol the absolute error bound for each segment
243  /// \param rel_tol the fractional error bound for each segment.
244  /// If the error is smaller than either the relative or absolute tolerance, the integration step is finished.
245  /// \param derivs number of derivatives to trust, which sets the order of the integrator. The order
246  /// is 3*\a derivs + 4. \a derivs can be 0, 1, or 2.
247  /// \param adapt if true, use recursive adaptation, otherwise do simple evaluation on the grid provided
248  /// with no error checking.
249  /// \param extrapolate if true, use simple Richardson extrapolation on the final 2 steps to reduce the error.
250  /// \return sum of partial integrals, which is the definite integral from the first value in \a xgrid to the last.
251  float_type integral(float_type amin, float_type amax, std::vector<float_type> *partials = 0,
252  float_type abs_tol=1e-12, float_type rel_tol=1e-12, int derivs=2, bool adapt=true, bool extrapolate=true)
253  const throw(c2_exception);
254 
255  /// \brief create a c2_piecewise_function_p from c2_connector_function_p segments which
256  /// is a representation of the parent function to the specified accuracy, but maybe much cheaper to evaluate
257  ///
258  /// This method has three modes, depending on the \a derivs flag.
259  ///
260  /// If \a derivs is 2,
261  /// it computes a c2_piecewise_function_p representation of its parent function, which may be a much faster
262  /// function to use in codes if the parent function is expensive. If \a xvals and \a yvals are non-null,
263  /// it will also fill them in with the function values at each grid point the adaptive algorithm chooses.
264  ///
265  /// If \a derivs is 1, this does not create the connectors,
266  /// and returns an null pointer, but will fill in the \a xvals and \a yvals
267  /// vectors with values of the function at points such that the linear interpolation error between the points
268  /// is bounded by the tolerance values given. Because it uses derivative information from the function to manage the
269  /// error control, it is almost completely free of issues with missing periods of oscillatory functions,
270  /// even with no information provided in the sampling grid.
271  /// This is typically useful for sampling a function for plotting.
272  ///
273  /// If \a derivs is 0, this does something very like what it does if \a derivs = 1, but without derivatives.
274  /// Instead, to compute the intermediate value of the function for error control, it just uses
275  /// 3-point parabolic interpolation. This is useful amost exclusively for converting a non-c2_function,
276  /// with no derivatives, but wrapped in a c2_classic_function wrapper, into a table of values to seed an interpolating_function_p.
277  /// Note, however, that without derivatives, this is very susceptible to missing periods of oscillatory
278  /// functions, so it is important to set a sampling grid which isn't too much coarser than the typical oscillations.
279  ///
280  /// \note the \a sampling_grid of the returned function matches the \a sampling_grid of its parent.
281  /// \see \ref sample_function_for_plotting "Adaptive Sampling Examples"
282  /// \param amin lower bound of the domain for sampling
283  /// \param amax upper bound of the domain for sampling
284  /// \param abs_tol the absolute error bound for each segment
285  /// \param rel_tol the fractional error bound for each segment.
286  /// \param derivs if 0 or 1, return a useless function, but fill in the \a xvals and \a yvals vectors (if non-null).
287  /// Also, if 0 or 1, tolerances refer to linear interpolation, not high-order interpolation.
288  /// If 2, return a full piecewise collection of c2_connector_function_p segments. See discussion above.
289  /// \param [in,out] xvals vector of abscissas at which the function was actually sampled (if non-null)
290  /// \param [in,out] yvals vector of function values corresponding to \a xvals (if non-null)
291  /// \return a new, sampled representation, if \a derivs is 2. A null pointer if \a derivs is 0 or 1.
292  c2_piecewise_function_p<float_type> *adaptively_sample(float_type amin, float_type amax,
293  float_type abs_tol=1e-12, float_type rel_tol=1e-12,
294  int derivs=2, std::vector<float_type> *xvals=0, std::vector<float_type> *yvals=0) const throw(c2_exception);
295 
296  /// \brief return the lower bound of the domain for this function as set by set_domain()
297  inline float_type xmin() const { return fXMin; }
298  /// \brief return the upper bound of the domain for this function as set by set_domain()
299  inline float_type xmax() const { return fXMax; }
300  /// \brief set the domain for this function.
301  void set_domain(float_type amin, float_type amax) { fXMin=amin; fXMax=amax; }
302 
303  /// \brief this is a counter owned by the function but which can be used to monitor efficiency of algorithms.
304  ///
305  /// It is not maintained automatically in general! The root finder, integrator, and sampler do increment it.
306  /// \return number of evaluations logged since last reset.
307  size_t get_evaluations() const { return evaluations; }
308  /// \brief reset the counter
309  void reset_evaluations() const { evaluations=0; } // evaluations are 'invisible' to constant
310  /// \brief count evaluations
311  inline void increment_evaluations() const { evaluations++; }
312 
313  /// \brief check that a vector is monotonic, throw an exception if not, and return a flag if it is reversed
314  ///
315  /// \param data a vector of data points which are expected to be monotonic.
316  /// \param message an informative string to include in an exception if this throws c2_exception
317  /// \return true if in decreasing order, false if increasing
318  bool check_monotonicity(const std::vector<float_type> &data, const char message[]) const throw(c2_exception);
319 
320  /// \brief establish a grid of 'interesting' points on the function.
321  ///
322  /// The sampling grid describes a reasonable initial set of points to look at the function.
323  /// this should generally be set at a scale which is quite coarse, and sufficient for initializing
324  /// adaptive integration or possibly root bracketing. For sampling a function to build a new interpolating
325  /// function, one may want to refine this for accuracy. However, interpolating_functions themselves
326  /// return their original X grid by default, so refining the grid in this case might be a bad idea.
327  /// \param grid a vector of abscissas. The contents is copied into an internal vector, so the \a grid can be discarded after passingin.
328  virtual void set_sampling_grid(const std::vector<float_type> &grid) throw(c2_exception);
329 
330  /// \brief get the sampling grid, which may be a null pointer
331  /// \return pointer to the sampling grid
332  std::vector<float_type> *get_sampling_grid_pointer() const { return sampling_grid; }
333 
334  /// \brief return the grid of 'interesting' points along this function which lie in the region requested
335  ///
336  /// if a sampling grid is defined, work from there, otherwise return vector of (amin, amax)
337  /// \param amin the lower bound for which the function is to be sampled
338  /// \param amax the upper bound for which the function is to be sampled
339  /// \param [in,out] grid filled vector containing the samplng grid.
340  virtual void get_sampling_grid(float_type amin, float_type amax, std::vector<float_type> &grid) const ;
341 
342  /// \brief clean up endpoints on a grid of points
343  /// \param[in,out] result the sampling grid with excessively closely space endpoints removed.
344  /// The grid is modified in place.
345  void preen_sampling_grid(std::vector<float_type> *result) const;
346  /// \brief refine a grid by splitting each interval into more intervals
347  /// \param [in,out] grid the grid to refine in place
348  /// \param refinement the number of new steps for each old step
349  void refine_sampling_grid(std::vector<float_type> &grid, size_t refinement) const;
350 
351  /// \brief create a new c2_function from this one which is normalized on the interval
352  /// \param amin lower bound of the domain for integration
353  /// \param amax upper bound of the domain for integration
354  /// \param norm the desired integral for the function over the region
355  /// \return a new c2_function with the desired \a norm.
356  c2_function<float_type> &normalized_function(float_type amin, float_type amax, float_type norm=1.0) const throw(c2_exception);
357  /// \brief create a new c2_function from this one which is square-normalized on the interval
358  /// \param amin lower bound of the domain for integration
359  /// \param amax upper bound of the domain for integration
360  /// \param norm the desired integral for the function over the region
361  /// \return a new c2_function with the desired \a norm.
362  c2_function<float_type> &square_normalized_function(float_type amin, float_type amax, float_type norm=1.0) const throw(c2_exception);
363  /// \brief create a new c2_function from this one which is square-normalized with the provided \a weight on the interval
364  /// \param amin lower bound of the domain for integration
365  /// \param amax upper bound of the domain for integration
366  /// \param weight a c2_function providing the weight
367  /// \param norm the desired integral for the function over the region
368  /// \return a new c2_function with the desired \a norm.
370  float_type amin, float_type amax, const c2_function<float_type> &weight, float_type norm=1.0)
371  const throw(c2_exception);
372 
373  /// \brief factory function to create a c2_sum_p from a regular algebraic expression.
374  /// \param rhs the right-hand term of the sum
375  /// \return a new c2_function
376  c2_sum_p<float_type> &operator + (const c2_function<float_type> &rhs) const
377  { return *new c2_sum_p<float_type>(*this, rhs); }
378  /// \brief factory function to create a c2_diff_p from a regular algebraic expression.
379  /// \param rhs the right-hand term of the difference
380  /// \return a new c2_function
382  { return *new c2_diff_p<float_type>(*this, rhs); }
383  /// \brief factory function to create a c2_product_p from a regular algebraic expression.
384  /// \param rhs the right-hand term of the product
385  /// \return a new c2_function
387  { return *new c2_product_p<float_type>(*this, rhs); }
388  /// \brief factory function to create a c2_ratio_p from a regular algebraic expression.
389  /// \param rhs the right-hand term of the ratio (the denominator)
390  /// \return a new c2_function
392  { return *new c2_ratio_p<float_type>(*this, rhs); }
393  /// \brief compose this function outside another.
394  /// \param inner the inner function
395  /// \return the composed function
396  /// \anchor compose_operator
398  { return *new c2_composed_function_p<float_type>((*this), inner); }
399 
400  /// \brief Find out where a calculation ran into trouble, if it got a nan.
401  /// If the most recent computation did not return a nan, this is undefined.
402  /// \return \a x value of point at which something went wrong, if integrator (or otherwise) returned a nan.
403  float_type get_trouble_point() const { return bad_x_point; }
404 
405  /// \brief increment our reference count. Destruction is only legal if the count is zero.
406  void claim_ownership() const { owner_count++; }
407  /// \brief decrement our reference count. Do not destroy at zero.
408  /// \return final owner count, to check whether object should disappear.
410  if(!owner_count) {
411  std::ostringstream outstr;
412  outstr << "attempt to release ownership of an unowned function in class ";
413  outstr << typeid(*this).name() << std::endl;
414  throw c2_exception(outstr.str().c_str());
415  }
416  owner_count--;
417  return owner_count;
418  }
419  /// \brief decrement our reference count. If the count reaches zero, destroy ourself.
420  void release_ownership() const throw(c2_exception) {
421  if(!release_ownership_for_return()) delete this;
422  }
423  /// \brief get the reference count, mostly for debugging
424  /// \return the count
425  size_t count_owners() const { return owner_count; }
426 
427 protected:
429  no_overwrite_grid(false),
430  fXMin(src.fXMin), fXMax(src.fXMax), root_info(0), owner_count(0)
431  {} // copy constructor only copies domain, and is only for internal use
434  fXMin(-std::numeric_limits<float_type>::max()),
435  fXMax(std::numeric_limits<float_type>::max()), root_info(0), owner_count(0)
436  {} // prevent accidental naked construction (impossible any since this has pure virtual methods)
437 
438  // this should only be called very early on, by a constructor, before anyone else
439  // sets a sampling grid, or it will leak memory
440  virtual void set_sampling_grid_pointer(std::vector<float_type> &grid)
441  {
442  if (sampling_grid && !no_overwrite_grid) delete sampling_grid; // grid was ours, lose it.
444  }
445 
446  std::vector<float_type> * sampling_grid;
448 
449  float_type fXMin, fXMax;
450  mutable size_t evaluations;
451  /// \brief this point may be used to record where a calculation ran into trouble
452  mutable float_type bad_x_point;
453 public:
454  /// \brief fill in a c2_fblock<float_type>... a shortcut for the integrator & sampler
455  /// \param [in,out] fb the block to fill in with information
456  inline void fill_fblock(c2_fblock<float_type> &fb) const throw(c2_exception)
457  {
458  fb.y=value_with_derivatives(fb.x, &fb.yp, &fb.ypp);
459  fb.ypbad=c2_isnan(fb.yp) || !c2_isfinite(fb.yp);
460  fb.yppbad=c2_isnan(fb.ypp) || !c2_isfinite(fb.ypp);
462  }
463 
464 private:
465  /// \brief the data element for the internal recursion stack for the sampler and integrator
466  struct recur_item {
467  c2_fblock<float_type> f1; size_t depth;
468  float_type previous_estimate, abs_tol, step_sum;
469  bool done;
470  size_t f0index, f2index;
471  };
472 
473 
474  /// \brief structure used to pass information recursively in integrator.
475  ///
476  /// the \a abs_tol is scaled by a factor of two at each division.
477  /// Everything else is just passed down.
478  struct c2_integrate_recur {
479  c2_fblock<float_type> *f0, *f1;
480  float_type abs_tol, rel_tol, eps_scale, extrap_coef, extrap2, dx_tolerance, abs_tol_min;
481  std::vector< recur_item > *rb_stack;
482  int derivs;
483  bool adapt, extrapolate, inited;
484  };
485 
486  /// \brief structure used to pass information recursively in sampler.
487  ///
488  struct c2_sample_recur {
489  c2_fblock<float_type> *f0, *f1;
490  float_type abs_tol, rel_tol, dx_tolerance, abs_tol_min;
491  int derivs;
493  std::vector<float_type> *xvals, *yvals;
494  std::vector< recur_item > *rb_stack;
495  bool inited;
496  };
497 
498  /// \brief structure used to hold root bracketing information
499  ///
500  struct c2_root_info {
501  c2_fblock<float_type> lower, upper;
502  bool inited;
503  };
504 
505  /// \brief Carry out the recursive subdivision and integration.
506  ///
507  /// This passes information recursively through the \a recur block pointer
508  /// to allow very efficient recursion.
509  /// \param rb a pointer to the recur struct.
510  float_type integrate_step(struct c2_integrate_recur &rb) const throw(c2_exception);
511 
512  /// \brief Carry out the recursive subdivision for sampling.
513  ///
514  /// This passes information recursively through the \a recur block pointer
515  /// to allow very efficient recursion.
516  /// \param rb a pointer to the recur struct.
517  void sample_step(struct c2_sample_recur &rb) const throw(c2_exception);
518 
519  /// this carry a memory of the last root bracketing,
520  /// to avoid the necessity of evaluating the function on the brackets every time
521  /// if the brackets have not been changed.
522  /// it is declared as a pointer, since many c2_functions may never need one allocated
523  mutable struct c2_root_info *root_info;
524 
525  mutable size_t owner_count;
526 };
527 
528 /// \brief a container into which any conventional c-style function can be dropped,
529 /// to create a degenerate c2_function without derivatives.
530 /// Mostly useful for sampling into interpolating functions.
531 /// construct a reference to this with c2_classic_function()
532 /// \ingroup containers
533 /// The factory function c2_factory::classic_function() creates *new c2_classic_function_p()
534 template <typename float_type=double> class c2_classic_function_p : public c2_function<float_type> {
535 public:
536  /// \brief construct the container
537  /// \param c_func a pointer to a conventional c-style function
538  c2_classic_function_p(const float_type (*c_func)(float_type)) : c2_function<float_type>(), func(c_func) {}
539 
540  /// \copydoc c2_function::value_with_derivatives
541  /// Uses the internal function pointer set by set_function().
542  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
543  {
544  if(!func) throw c2_exception("c2_classic_function called with null function");
545  if(yprime) *yprime=0;
546  if(yprime2) *yprime2=0;
547  return func(x);
548  }
550 
551 protected:
552  /// \brief pointer to our function
553  const float_type (*func)(float_type);
554 };
555 
556 /// \brief create a container for a c2_function which handles the reference counting.
557 /// \ingroup containers
558 /// It is useful as a smart container to hold a c2_function and keep the reference count correct.
559 /// The recommended way for a class to store a c2_function which is handed in from the outside
560 /// is for it to have a c2_ptr member into which the passed-in function is stored.
561 /// This way, when the class instance is deleted, it will automatically dereference any function
562 /// which it was handed.
563 ///
564 /// This class contains a copy constructor and operator=, to make it fairly easy to make
565 /// a std::vector of these objects, and have it work as expected.
566 template <typename float_type> class c2_const_ptr {
567 public:
568  /// \brief construct the container with no function
569  c2_const_ptr() : func(0) {}
570  /// \brief construct the container with a pre-defined function
571  /// \param f the function to store
573  { this->set_function(&f); }
574  /// \brief copy constructor
575  /// \param src the container to copy
577  { this->set_function(src.get_ptr()); }
578  /// \brief fill the container with a new function, or clear it with a null pointer
579  /// \param f the function to store, releasing any previously held function
581  {
582  if(func) func->release_ownership();
583  func=f;
584  if(func) func->claim_ownership();
585  }
586 
587  /// \brief fill the container from another container
588  /// \param f the container to copy
590  { this->set_function(f.get_ptr()); return f; }
591  /// \brief fill the container with a function
592  /// \param f the function
594  { this->set_function(&f); return f; }
595  /// \brief release the function without destroying it, so it can be returned from a function
596  ///
597  /// This is usually the very last line of a function before the return statement, so that
598  /// any exceptions that happen during execution of the function will cause proper cleanup.
599  /// Once the function has been released from its container this way, it is an orhpaned object
600  /// until the caller claims it, so it could get lost if an exception happens.
602  {
603  if(func) func->release_ownership_for_return();
604  func=0;
605  }
606  /// \brief clear the function
607  ///
608  /// Any attempt to use this c2_plugin_function_p throws an exception if the saved function is cleared.
609  void unset_function(void) { this->set_function(0); }
610  /// \brief destructor
611  ~c2_const_ptr() { this->set_function(0); }
612 
613  /// \brief get a reference to our owned function
614  inline const c2_function<float_type> &get() const throw(c2_exception)
615  {
616  if(!func) throw c2_exception("c2_ptr accessed uninitialized");
617  return *func;
618  }
619  /// \brief get an unchecked pointer to our owned function
620  inline const c2_function<float_type> *get_ptr() const { return func; }
621  /// \brief get a checked pointer to our owned function
622  inline const c2_function<float_type> *operator -> () const
623  { return &get(); }
624  /// \brief check if we have a valid function
625  bool valid() const { return func != 0; }
626 
627  /// \brief type coercion operator which lets us use a pointer as if it were a const c2_function
628  operator const c2_function<float_type>& () const { return this->get(); }
629 
630  /// \brief convenience operator to make us look like a function
631  /// \param x the value at which to evaluate the contained function
632  /// \return the evaluated function
633  /// \note If you using this repeatedly, do const c2_function<float_type> &func=ptr;
634  /// and use func(x). Calling this operator wastes some time, since it checks the validity of the
635  /// pointer every time.
636  float_type operator()(float_type x) const throw(c2_exception) { return get()(x); }
637  /// \brief convenience operator to make us look like a function
638  /// \param x the value at which to evaluate the contained function
639  /// \param yprime the derivative
640  /// \param yprime2 the second derivative
641  /// \return the evaluated function
642  /// \note If you using this repeatedly, do const c2_function<float_type> &func=ptr;
643  /// and use func(x). Calling this operator wastes some time, since it checks the validity of the
644  /// pointer every time.
645  float_type operator()(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
646  { return get().value_with_derivatives(x, yprime, yprime2); }
647  /// \brief factory function to create a c2_sum_p from a regular algebraic expression.
648  /// \param rhs the right-hand term of the sum
649  /// \return a new c2_function
651  { return *new c2_sum_p<float_type>(get(), rhs); }
652  /// \brief factory function to create a c2_diff_p from a regular algebraic expression.
653  /// \param rhs the right-hand term of the difference
654  /// \return a new c2_function
656  { return *new c2_diff_p<float_type>(get(), rhs); }
657  /// \brief factory function to create a c2_product_p from a regular algebraic expression.
658  /// \param rhs the right-hand term of the product
659  /// \return a new c2_function
661  { return *new c2_product_p<float_type>(get(), rhs); }
662  /// \brief factory function to create a c2_ratio_p from a regular algebraic expression.
663  /// \param rhs the right-hand term of the ratio (the denominator)
664  /// \return a new c2_function
666  { return *new c2_ratio_p<float_type>(get(), rhs); }
667  /// \brief compose this function outside another.
668  /// \param inner the inner function
669  /// \return the composed function
671  { return *new c2_composed_function_p<float_type>(get(), inner); }
672 
673 protected:
675 };
676 
677 /// \brief create a container for a c2_function which handles the reference counting.
678 /// \ingroup containers
679 ///
680 /// \see c2_const_ptr and \ref memory_management "Use of c2_ptr for memory management"
681 template <typename float_type> class c2_ptr : public c2_const_ptr<float_type >
682 {
683 public:
684  /// \brief construct the container with no function
685  c2_ptr() : c2_const_ptr<float_type>() {}
686  /// \brief construct the container with a pre-defined function
687  /// \param f the function to store
689  c2_const_ptr<float_type>() { this->set_function(&f); }
690  /// \brief copy constructor
691  /// \param src the container to copy
692  c2_ptr(const c2_ptr<float_type> &src) :
693  c2_const_ptr<float_type>() { this->set_function(src.get_ptr()); }
694  /// \brief get a checked pointer to our owned function
695  inline c2_function<float_type> &get() const throw(c2_exception)
696  { return *const_cast<c2_function<float_type>*>(&c2_const_ptr<float_type>::get()); }
697  /// \brief get an unchecked pointer to our owned function
699  { return const_cast<c2_function<float_type>*>(this->func); }
700  /// \brief get a checked pointer to our owned function
702  { return &get(); }
703  /// \brief fill the container from another container
704  /// \param f the container to copy
706  { this->set_function(f.get_ptr()); return f; }
707  /// \brief fill the container with a function
708  /// \param f the function
710  { this->set_function(&f); return f; }
711 private:
712  /// \brief hidden non-const-safe version of operator=
713  void operator =(const c2_const_ptr<float_type> &) { }
714  /// \brief hidden non-const-safe version of operator=
715  void operator =(const c2_function<float_type> &) { }
716 };
717 
718 /// \brief create a non-generic container for a c2_function which handles the reference counting.
719 /// \ingroup containers
720 ///
721 /// \see c2_const_ptr and \ref memory_management "Use of c2_ptr for memory management"
722 ///
723 /// \note Overuse of this class will generate massive bloat. Use c2_ptr and c2_const_ptr if you don't _really_ need specific pointer types.
724 /// \see \ref memory_management "Use of c2_ptr for memory management"
725 template <typename float_type, template <typename> class c2_class > class c2_typed_ptr : public c2_const_ptr<float_type> {
726 public:
727  /// \brief construct the container with no function
728  c2_typed_ptr() : c2_ptr<float_type>() {}
729  /// \brief construct the container with a pre-defined function
730  /// \param f the function to store
731  c2_typed_ptr(c2_class<float_type> &f)
732  : c2_const_ptr<float_type>() { this->set_function(&f); }
733  /// \brief copy constructor
734  /// \param src the container to copy
736  : c2_const_ptr<float_type>() { this->set_function(src.get_ptr()); }
737 
738  /// \brief get a reference to our owned function
739  inline c2_class<float_type> &get() const throw(c2_exception)
740  {
741  return *static_cast<c2_class<float_type> *>(const_cast<c2_function<float_type>*>(&c2_const_ptr<float_type>::get()));
742  }
743  /// \brief get a checked pointer to our owned function
744  inline c2_class<float_type> *operator -> () const
745  { return &get(); }
746  /// \brief get an unchecked pointer to our owned function
747  inline c2_class<float_type> *get_ptr() const
748  { return static_cast<c2_class<float_type> *>(const_cast<c2_function<float_type>*>(this->func)); }
749  /// \brief type coercion operator which lets us use a pointer as if it were a c2_function
750  operator c2_class<float_type>&() const { return get(); }
751  /// \brief fill the container from another container
752  /// \param f the container to copy
754  { this->set_function(f.get_ptr()); }
755  /// \brief fill the container with a function
756  /// \param f the function
757  void operator =(c2_class<float_type> &f)
758  { this->set_function(&f); }
759 private:
760  /// \brief hidden downcasting version of operator=
761  void operator =(const c2_const_ptr<float_type> &) { }
762  /// \brief hidden downcasting version of operator=. Use an explicit dynamic_cast<c2_class<float_type>&>(f) if you need to try this.
763  void operator =(const c2_function<float_type> &) { }
764 };
765 
766 /// \brief a container into which any other c2_function can be dropped, to allow expressions
767 /// with replacable components.
768 /// \ingroup containers
769 ///It is useful for plugging different InterpolatingFunctions into a c2_function expression.
770 ///It saves a lot of effort in other places with casting away const declarations.
771 ///
772 /// It is also useful as a wrapper for a function if it is necessary to have a copy of a function
773 /// which has a different domain or sampling grid than the parent function. This can be
774 /// be used, for example, to patch badly-behaved functions with c2_piecewise_function_p by
775 /// taking the parent function, creating two plugins of it with domains on each side of the
776 /// nasty bit, and then inserting a nice function in the hole.
777 ///
778 /// This can also be used as a fancier c2_ptr which allows direct evaluation
779 /// instead of having to dereference the container first.
780 ///
781 /// The factory function c2_factory::plugin_function() creates *new c2_plugin_function_p()
782 template <typename float_type=double> class c2_plugin_function_p :
783  public c2_function<float_type> {
784 public:
785  /// \brief construct the container with no function
786  c2_plugin_function_p() : c2_function<float_type>(), func() {}
787  /// \brief construct the container with a pre-defined function
789  c2_function<float_type>(),func(f) { }
790  /// \brief fill the container with a new function, or clear it with a null pointer
791  /// and copy our domain
793  {
794  func.set_function(f);
795  if(f) this->set_domain(f->xmin(), f->xmax());
796  }
797  /// \copydoc c2_function::value_with_derivatives
798  /// Uses the internal function pointer set by set_function().
799  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
800  {
801  if(!func.valid()) throw c2_exception("c2_plugin_function_p called uninitialized");
802  return func->value_with_derivatives(x, yprime, yprime2);
803  }
804  /// \brief destructor
805  virtual ~c2_plugin_function_p() { }
806 
807  /// \brief clear our function
808  void unset_function() { func.unset_function(); }
809 
810  virtual void get_sampling_grid(float_type amin, float_type amax, std::vector<float_type> &grid) const
811  {
812  if(!func.valid()) throw c2_exception("c2_plugin_function_p called uninitialized");
813  if(this->sampling_grid) c2_function<float_type>::get_sampling_grid(amin, amax, grid);
814  else func->get_sampling_grid(amin, amax, grid);
815  }
816 protected:
818 };
819 
820 /// \brief a c2_plugin_function_p which promises not to fiddle with the plugged function.
821 /// \ingroup containers
822 ///
823 /// The factory function c2_factory::const_plugin_function() creates *new c2_const_plugin_function_p()
824 template <typename float_type=double> class c2_const_plugin_function_p : public c2_plugin_function_p<float_type> {
825 public:
826  /// \brief construct the container with no function
828  /// \brief construct the container with a pre-defined function
830  c2_plugin_function_p<float_type>() { this->set_function(&f); }
831  /// \brief fill the container with a new function, or clear it with a null pointer
834  /// \brief destructor
836 
837  /// \brief get a const reference to our owned function, for direct access
838  const c2_function<float_type> &get() const throw(c2_exception)
839  { return this->func.get(); }
840 };
841 
842 /// \brief Provides support for c2_function objects which are constructed from two other c2_function
843 /// objects.
844 /// \ingroup abstract_classes
845 template <typename float_type=double> class c2_binary_function : public c2_function<float_type> {
846 public:
847  /// \brief function to manage the binary operation, used by c2_binary_function::value_with_derivatives()
848  ///
849 
850  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw (c2_exception)
851  {
852  if(stub) throw c2_exception("attempt to evaluate a c2_binary_function stub");
853  return this->combine(*Left.get_ptr(), *Right.get_ptr(), x, yprime, yprime2);
854  }
855 
856  /// \brief destructor releases ownership of member functions
857  ///
858  virtual ~c2_binary_function() { }
859 
860 protected:
861  /// \brief construct the binary function
862  /// \param combiner pointer to the function which actualy knows how to execute the binary
863  /// \param left the c2_function to be used in the left side of the binary relation
864  /// \param right the c2_function to be used in the right side of the binary relation
866  float_type (*combiner)(const c2_function<float_type> &left, const c2_function<float_type> &right,
867  float_type x, float_type *yprime, float_type *yprime2),
868  const c2_function<float_type> &left, const c2_function<float_type> &right) :
869  c2_function<float_type>(), combine(combiner), Left(left), Right(right), stub(false)
870  {
871  this->set_domain(
872  (left.xmin() > right.xmin()) ? left.xmin() : right.xmin(),
873  (left.xmax() < right.xmax()) ? left.xmax() : right.xmax()
874  );
875  }
876 
877  /// \brief construct a 'stub' c2_binary_function, which provides access to the combine() function
878  /// \note Do not evaluate a 'stub' ever. It is only used so that combine() can be called
880  float_type (*combiner)(const c2_function<float_type> &left, const c2_function<float_type> &right,
881  float_type x, float_type *yprime, float_type *yprime2)
882  ) : c2_function<float_type>(), combine(combiner), Left(), Right(), stub(true) { }
883 
884 public:
885  float_type (* const combine)(const c2_function<float_type> &left, const c2_function<float_type> &right,
886  float_type x, float_type *yprime, float_type *yprime2);
887 
888 protected:
890  /// \brief if true, we don't own any functions, we are just a source of a combining function.
891  bool stub;
892 
893 };
894 
895 /// \brief Create a very lightweight method to return a scalar multiple of another function.
896 /// \ingroup containers \ingroup arithmetic_functions \ingroup parametric_functions
897 ///
898 /// The factory function c2_factory::scaled_function() creates *new c2_scaled_function_p
899 template <typename float_type=double> class c2_scaled_function_p : public c2_function<float_type> {
900 public:
901  /// \brief construct the function with its scale factor.
902  ///
903  /// \param outer the function to be scaled
904  /// \param scale the multiplicative scale factor
905  c2_scaled_function_p(const c2_function<float_type> &outer, float_type scale) :
906  c2_function<float_type>(), func(outer), yscale(scale) { }
907 
908  /// \brief set a new scale factor
909  /// \param scale the new factor
910  void reset(float_type scale) { yscale=scale; }
911 
912  /// \copydoc c2_function::value_with_derivatives
913  ///
914  /// provide our own value_with_derivatives which bypasses the combiner for quicker operation
915  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw (c2_exception)
916  {
917  float_type y=this->func->value_with_derivatives(x, yprime, yprime2);
918  if(yprime) (*yprime)*=yscale;
919  if(yprime2) (*yprime2)*=yscale;
920  return y*yscale;
921  }
922 
923 protected:
924  c2_scaled_function_p<float_type>() : func() {} // hide default constructor, since its use is almost always an error.
925  /// \brief the scaling factor for the function
927  float_type yscale;
928 };
929 
930 /// \brief A container into which any other c2_function can be dropped.
931 /// \ingroup containers
932 /// It allows a function to be pre-evaluated at a point, and used at multiple places in an expression
933 /// efficiently. If it is re-evaluated at the previous point, it returns the remembered values;
934 /// otherwise, it re-evauates the function at the new point.
935 ///
936 /// The factory function c2_factory::cached_function() creates *new c2_cached_function_p
937 template <typename float_type=double> class c2_cached_function_p : public c2_function<float_type> {
938 public:
939  /// \brief construct the container
940  ///
941  /// \param f the function to be cached
943  func(f), init(false) {}
944  /// \copydoc c2_function::value_with_derivatives
945  ///
946  /// Checks to see if the function is being re-evaluated at the previous point, and
947  /// returns remembered values if so.
948  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
949  {
950  if(!init || x != x0) {
951  y=this->func->value_with_derivatives(x, &yp, &ypp);
952  x0=x;
953  init=true;
954  }
955  if(yprime) *yprime=yp;
956  if(yprime2) *yprime2=ypp;
957  return y;
958  }
959 
960 protected:
961  c2_cached_function_p() : func() {} // hide default constructor, since its use is almost always an error.
963  mutable bool init;
964  mutable float_type x0, y, yp, ypp;
965 
966 };
967 
968 /// \brief Provides function composition (nesting)
969 /// \ingroup arithmetic_functions
970 /// This allows evaluation of \a f(g(x)) where \a f and \a g are c2_function objects.
971 ///
972 /// This should always be constructed using \ref compose_operator "c2_function::operator()"
973 template <typename float_type=double> class c2_composed_function_p : public c2_binary_function<float_type> {
974 public:
975 
976  /// \brief construct \a outer( \a inner (x))
977  /// \note See c2_binary_function for discussion of ownership.
978  /// \param outer the outer function
979  /// \param inner the inner function
981  c2_binary_function<float_type>(combine, outer, inner) { this->set_domain(inner.xmin(), inner.xmax()); }
982  /// \brief Create a stub just for the combiner to avoid statics.
984 
985  /// \brief execute math necessary to do composition
987  float_type x, float_type *yprime, float_type *yprime2) throw(c2_exception)
988  {
989  float_type y0, y1;
990  if(yprime || yprime2) {
991  float_type yp0, ypp0, yp1, ypp1;
992  y0=right.value_with_derivatives(x, &yp0, &ypp0);
993  y1=left.value_with_derivatives(y0, &yp1, &ypp1);
994  if(yprime) *yprime=yp1*yp0;
995  if(yprime2) *yprime2=ypp0*yp1+yp0*yp0*ypp1;
996  } else {
997  y0=right(x);
998  y1=left(y0);
999  }
1000  return y1;
1001  }
1002 };
1003 
1004 /// \brief create a c2_function which is the sum of two other c2_function objects.
1005 /// \ingroup arithmetic_functions
1006 /// This should always be constructed using c2_function::operator+()
1007 template <typename float_type=double> class c2_sum_p : public c2_binary_function<float_type> {
1008 public:
1009  /// \brief construct \a left + \a right
1010  /// \param left the left function
1011  /// \param right the right function
1013  /// \brief Create a stub just for the combiner to avoid statics.
1014  c2_sum_p() : c2_binary_function<float_type>(combine) {} ; // create a stub just for the combiner to avoid statics
1015 
1016  /// \brief execute math necessary to do addition
1018  float_type x, float_type *yprime, float_type *yprime2) throw(c2_exception)
1019  {
1020  float_type y0, y1;
1021  if(yprime || yprime2) {
1022  float_type yp0, ypp0, yp1, ypp1;
1023  y0=left.value_with_derivatives(x, &yp0, &ypp0);
1024  y1=right.value_with_derivatives(x, &yp1, &ypp1);
1025  if(yprime) *yprime=yp0+yp1;
1026  if(yprime2) *yprime2=ypp0+ypp1;
1027  } else {
1028  y0=left(x);
1029  y1=right(x);
1030  }
1031  return y0+y1;
1032  }
1033 };
1034 
1035 
1036 /// \brief create a c2_function which is the difference of two other c2_functions.
1037 /// \ingroup arithmetic_functions
1038 /// This should always be constructed using c2_function::operator-()
1039 template <typename float_type=double> class c2_diff_p : public c2_binary_function<float_type> {
1040 public:
1041  /// \brief construct \a left - \a right
1042  /// \param left the left function
1043  /// \param right the right function
1045  /// \brief Create a stub just for the combiner to avoid statics.
1046  c2_diff_p() : c2_binary_function<float_type>(combine) {} ; // create a stub just for the combiner to avoid statics
1047 
1048  /// \brief execute math necessary to do subtraction
1050  float_type x, float_type *yprime, float_type *yprime2) throw(c2_exception)
1051  {
1052  float_type y0, y1;
1053  if(yprime || yprime2) {
1054  float_type yp0, ypp0, yp1, ypp1;
1055  y0=left.value_with_derivatives(x, &yp0, &ypp0);
1056  y1=right.value_with_derivatives(x, &yp1, &ypp1);
1057  if(yprime) *yprime=yp0-yp1;
1058  if(yprime2) *yprime2=ypp0-ypp1;
1059  } else {
1060  y0=left(x);
1061  y1=right(x);
1062  }
1063  return y0-y1;
1064  }
1065 };
1066 
1067 
1068 /// \brief create a c2_function which is the product of two other c2_functions.
1069 /// \ingroup arithmetic_functions
1070 /// This should always be constructed using c2_function::operator*()
1071 template <typename float_type=double> class c2_product_p : public c2_binary_function<float_type> {
1072 public:
1073  /// \brief construct \a left * \a right
1074  /// \param left the left function
1075  /// \param right the right function
1077  /// \brief Create a stub just for the combiner to avoid statics.
1078  c2_product_p() : c2_binary_function<float_type>(combine) {} ; // create a stub just for the combiner to avoid statics
1079 
1080  /// \brief execute math necessary to do multiplication
1082  float_type x, float_type *yprime, float_type *yprime2) throw(c2_exception)
1083  {
1084  float_type y0, y1;
1085  if(yprime || yprime2) {
1086  float_type yp0, ypp0, yp1, ypp1;
1087  y0=left.value_with_derivatives(x, &yp0, &ypp0);
1088  y1=right.value_with_derivatives(x, &yp1, &ypp1);
1089  if(yprime) *yprime=y1*yp0+y0*yp1;
1090  if(yprime2) *yprime2=ypp0*y1+2.0*yp0*yp1+ypp1*y0;
1091  } else {
1092  y0=left(x);
1093  y1=right(x);
1094  }
1095  return y0*y1;
1096  }
1097 };
1098 
1099 
1100 /// \brief create a c2_function which is the ratio of two other c2_functions.
1101 /// \ingroup arithmetic_functions
1102 /// This should always be constructed using c2_function::operator/()
1103 template <typename float_type=double> class c2_ratio_p : public c2_binary_function<float_type> {
1104 public:
1105  /// \brief construct \a left / \a right
1106  /// \param left the left function
1107  /// \param right the right function
1109  /// \brief Create a stub just for the combiner to avoid statics.
1110  c2_ratio_p() : c2_binary_function<float_type>(combine) {} ; // create a stub just for the combiner to avoid statics
1111 
1112  /// \brief execute math necessary to do division
1114  float_type x, float_type *yprime, float_type *yprime2) throw(c2_exception)
1115  {
1116  float_type y0, y1;
1117  if(yprime || yprime2) {
1118  float_type yp0, ypp0, yp1, ypp1;
1119  y0=left.value_with_derivatives(x, &yp0, &ypp0);
1120  y1=right.value_with_derivatives(x, &yp1, &ypp1);
1121  if(yprime) *yprime=(yp0*y1-y0*yp1)/(y1*y1); // first deriv of ratio
1122  if(yprime2) *yprime2=(y1*y1*ypp0+y0*(2*yp1*yp1-y1*ypp1)-2*y1*yp0*yp1)/(y1*y1*y1);
1123  } else {
1124  y0=left(x);
1125  y1=right(x);
1126  }
1127  return y0/y1;
1128  }
1129 
1130 };
1131 
1132 /// \brief a c2_function which is constant
1133 /// \ingroup parametric_functions
1134 ///
1135 /// The factory function c2_factory::constant() creates *new c2_constant_p()
1136 template <typename float_type> class c2_constant_p : public c2_function<float_type> {
1137 public:
1138  c2_constant_p(float_type x) : c2_function<float_type>(), value(x) {}
1139  void reset(float_type val) { value=val; }
1140  virtual float_type value_with_derivatives(float_type, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1141  { if(yprime) *yprime=0; if(yprime2) *yprime2=0; return value; }
1142 
1143 private:
1144  float_type value;
1145 };
1146 
1147 /// \brief a transformation of a coordinate, including an inverse
1148 /// \ingroup transforms
1149 template <typename float_type> class c2_transformation {
1150 public:
1151  /// \brief initialize all our function pointers
1152  /// \param transformed true if this function is not the identity
1153  /// \param xin input X transform
1154  /// \param xinp input X transform derivative
1155  /// \param xinpp input X transform second derivative
1156  /// \param xout output X transform, which MUST be the inverse of \a xin
1157  c2_transformation(bool transformed,
1158  float_type (*xin)(float_type), float_type (*xinp)(float_type), float_type (*xinpp)(float_type), float_type (*xout)(float_type)
1159  ) :
1160  fTransformed(transformed), fHasStaticTransforms(true),
1161  pIn(xin), pInPrime(xinp), pInDPrime(xinpp), pOut(xout) { }
1162 
1163  /// \brief initialize all our function pointers so that only the (overridden) virtual functions can be called without an error
1164  /// \param transformed true if this function is nonlinear
1165  c2_transformation(bool transformed) :
1166  fTransformed(transformed), fHasStaticTransforms(false),
1168  /// \brief the destructor
1169  virtual ~c2_transformation() { }
1170  /// \brief flag to indicate if this transform is not the identity
1171  const bool fTransformed;
1172  /// \brief flag to indicate if the static function pointers can be used for efficiency
1174 
1175  /// \note the pointers to functions allow highly optimized access when static functions are available.
1176  /// They are only used inside value_with_derivatives(), which is assumed to be the most critical routine.
1177  /// \brief non-virtual pointer to input X transform
1178  float_type (* const pIn)(float_type);
1179  /// \brief non-virtual pointer to input X transform derivative
1180  float_type (* const pInPrime)(float_type);
1181  /// \brief non-virtual pointer to input X transform second derivative
1182  float_type (* const pInDPrime)(float_type);
1183  /// \brief non-virtual pointer to output X transform
1184  float_type (* const pOut)(float_type);
1185 
1186  /// \brief virtual input X transform
1187  virtual float_type fIn(float_type x) const { return pIn(x); }
1188  /// \brief virtual input X transform derivative
1189  virtual float_type fInPrime(float_type x) const { return pInPrime(x); }
1190  /// \brief virtual input X transform second derivative
1191  virtual float_type fInDPrime(float_type x) const { return pInDPrime(x); }
1192  /// \brief virtual output X transform
1193  virtual float_type fOut(float_type x) const { return pOut(x); }
1194 
1195 protected:
1196  /// \brief utility function for unimplemented conversion
1197  static float_type report_error(float_type x) { throw c2_exception("use of improperly constructed axis transform"); return x; }
1198  /// \brief utility function f(x)=x useful in axis transforms
1199  static float_type ident(float_type x) { return x; }
1200  /// \brief utility function f(x)=1 useful in axis transforms
1201  static float_type one(float_type) { return 1; }
1202  /// \brief utility function f(x)=0 useful in axis transforms
1203  static float_type zero(float_type) { return 0; }
1204  /// \brief utility function f(x)=1/x useful in axis transforms
1205  static float_type recip(float_type x) { return 1.0/x; }
1206  /// \brief utility function f(x)=-1/x**2 useful in axis transforms
1207  static float_type recip_prime(float_type x) { return -1/(x*x); }
1208  /// \brief utility function f(x)=2/x**3 useful in axis transforms
1209  static float_type recip_prime2(float_type x) { return 2/(x*x*x); }
1210 
1211 };
1212 
1213 /// \brief the identity transform
1214 /// \ingroup transforms
1215 template <typename float_type> class c2_transformation_linear : public c2_transformation<float_type> {
1216 public:
1217  /// \brief constructor
1218  c2_transformation_linear() : c2_transformation<float_type>(false, this->ident, this->one, this->zero, this->ident) { }
1219  /// \brief destructor
1221 };
1222 /// \brief log axis transform
1223 /// \ingroup transforms
1224 template <typename float_type> class c2_transformation_log : public c2_transformation<float_type> {
1225 public:
1226  /// \brief constructor
1227  c2_transformation_log() : c2_transformation<float_type>(true, std::log, this->recip, this->recip_prime, std::exp) { }
1228  /// \brief destructor
1230 };
1231 /// \brief reciprocal axis transform
1232 /// \ingroup transforms
1233 template <typename float_type> class c2_transformation_recip : public c2_transformation<float_type> {
1234 public:
1235  /// \brief constructor
1236  c2_transformation_recip() : c2_transformation<float_type>(true, this->recip, this->recip_prime, this->recip_prime2, this->recip) { }
1237  /// \brief destructor
1239 };
1240 
1241 /// \brief a transformation of a function in and out of a coordinate space, using 2 c2_transformations
1242 ///
1243 /// This class is a container for two axis transforms, but also provides the critical evaluate()
1244 /// function which converts a result in internal coordinates (with derivatives) into the external representation
1245 /// \ingroup transforms
1246 template <typename float_type>
1248 public:
1249  /// \brief construct this from two c2_transformation instances
1250  /// \param xx the X axis transform
1251  /// \param yy the Y axis transform
1254  isIdentity(!(xx.fTransformed || yy.fTransformed)), X(xx), Y(yy) { }
1255  /// \brief destructor
1256  virtual ~c2_function_transformation() { delete &X; delete &Y; }
1257  /// \brief evaluate the transformation from internal coordinates to external coordinates
1258  /// \param xraw the value of \a x in external cordinates at which the transform is taking place
1259  /// \param y the value of the function in internal coordinates
1260  /// \param yp0 the derivative in internal coordinates
1261  /// \param ypp0 the second derivative in internal coordinates
1262  /// \param [out] yprime pointer to the derivative, or NULL, in external coordinates
1263  /// \param [out] yprime2 pointer to the second derivative, or NULL, in external coordinates
1264  /// \return the value of the function in external coordinates
1265  virtual float_type evaluate(float_type xraw,
1266  float_type y, float_type yp0, float_type ypp0,
1267  float_type *yprime, float_type *yprime2) const;
1268  /// \brief flag indicating of the transform is the identity, and can be skipped for efficiency
1269  const bool isIdentity;
1270  /// \brief the X axis transform
1272  /// \brief the Y axis transform
1274 };
1275 
1276 /// \brief a transformation of a function in and out of lin-lin space
1277 ///
1278 /// \ingroup transforms
1279 template <typename float_type> class c2_lin_lin_function_transformation :
1280  public c2_function_transformation<float_type> {
1281 public:
1283  c2_function_transformation<float_type>(
1284  *new c2_transformation_linear<float_type>,
1285  *new c2_transformation_linear<float_type>
1286  ) { }
1288 };
1289 
1290 /// \brief a transformation of a function in and out of log-log space
1291 ///
1292 /// \ingroup transforms
1293 template <typename float_type> class c2_log_log_function_transformation :
1294  public c2_function_transformation<float_type> {
1295 public:
1297  c2_function_transformation<float_type>(
1298  *new c2_transformation_log<float_type>,
1299  *new c2_transformation_log<float_type>
1300  ) { }
1302 };
1303 
1304 /// \brief a transformation of a function in and out of lin-log space
1305 ///
1306 /// \ingroup transforms
1307 template <typename float_type> class c2_lin_log_function_transformation :
1308  public c2_function_transformation<float_type> {
1309 public:
1311  c2_function_transformation<float_type>(
1312  *new c2_transformation_linear<float_type>,
1313  *new c2_transformation_log<float_type>
1314  ) { }
1316 };
1317 
1318 /// \brief a transformation of a function in and out of log-lin space
1319 ///
1320 /// \ingroup transforms
1321 template <typename float_type> class c2_log_lin_function_transformation :
1322  public c2_function_transformation<float_type> {
1323 public:
1325  c2_function_transformation<float_type>(
1326  *new c2_transformation_log<float_type>,
1327  *new c2_transformation_linear<float_type>
1328  ) { }
1330 };
1331 
1332 /// \brief a transformation of a function in and out of Arrhenius (1/x vs. log(y)) space
1333 ///
1334 /// \ingroup transforms
1335 template <typename float_type> class c2_arrhenius_function_transformation :
1336  public c2_function_transformation<float_type> {
1337 public:
1339  c2_function_transformation<float_type>(
1340  *new c2_transformation_recip<float_type>,
1341  *new c2_transformation_log<float_type>
1342  ) { }
1344 };
1345 
1346 /**
1347  \brief create a cubic spline interpolation of a set of (x,y) pairs
1348  \ingroup interpolators
1349  This is one of the main reasons for c2_function objects to exist.
1350 
1351  It provides support for cubic spline interpolation of data provides from tables of \a x, \a y pairs.
1352  It supports automatic, transparent linearization of the data before storing in its tables (through
1353  subclasses such as
1354  log_lin_interpolating_function, lin_log_interpolating_function, and
1355  log_log_interpolating_function) to permit very high accuracy representations of data which have a suitable
1356  structure. It provides utility functions LinearInterpolatingGrid() and LogLogInterpolatingGrid()
1357  to create grids for mapping other functions onto a arithmetic or geometric grid.
1358 
1359  In its simplest form, an untransformed cubic spline of a data set, using natural boundary conditions
1360  (vanishing second derivative), is created as: \n
1361  \code
1362  c2_ptr<double> c2p;
1363  c2_factory<double> c2;
1364  std::vector<double> xvals(10), yvals(10);
1365  // < fill in xvals and yvals >
1366  c2p myfunc=c2.interpolating_function().load(xvals, yvals,true,0,true,0);
1367  // and it can be evaluated at a point for its value only by:
1368  double y=myfunc(x);
1369  // or it can be evaluated with its derivatives by
1370  double yprime, yprime2;
1371  double y=myfunc(x,&yprime, &yprime2);
1372  \endcode
1373 
1374  The factory function c2_factory::interpolating_function() creates *new interpolating_function_p()
1375 */
1376 
1377 template <typename float_type=double> class interpolating_function_p : public c2_function<float_type> {
1378 public:
1379  /// \brief an empty linear-linear cubic-spline interpolating_function_p
1380  ///
1381  /// lots to say here, but see Numerical Recipes for a discussion of cubic splines.
1382  ///
1384  fTransform(*new c2_lin_lin_function_transformation<float_type>) { }
1385 
1386  /// \brief an empty cubic-spline interpolating_function_p with a specific transform
1387  ///
1389  fTransform(transform) { }
1390 
1391  /// \brief do the dirty work of constructing the spline from a function.
1392  /// \param x the list of abscissas. Must be either strictly increasing or strictly decreasing.
1393  /// Strictly increasing is preferred, as less memory is used since a copy is not required for the sampling grid.
1394  /// \param f the list of function values.
1395  /// \param lowerSlopeNatural if true, set y''(first point)=0, otherwise compute it from \a lowerSope
1396  /// \param lowerSlope derivative of the function at the lower bound, used only if \a lowerSlopeNatural is false
1397  /// \param upperSlopeNatural if true, set y''(last point)=0, otherwise compute it from \a upperSope
1398  /// \param upperSlope derivative of the function at the upper bound, used only if \a upperSlopeNatural is false
1399  /// \param splined if true (default), use cubic spline, if false, use linear interpolation.
1400  /// \return the same interpolating function, filled
1401  interpolating_function_p<float_type> & load(const std::vector<float_type> &x, const std::vector<float_type> &f,
1402  bool lowerSlopeNatural, float_type lowerSlope,
1403  bool upperSlopeNatural, float_type upperSlope, bool splined=true
1404  ) throw(c2_exception);
1405 
1406  /// \brief do the dirty work of constructing the spline from a function.
1407  /// \param data std::vector of std::pairs of x,y. Will be sorted into x increasing order in place.
1408  /// \param lowerSlopeNatural if true, set y''(first point)=0, otherwise compute it from \a lowerSope
1409  /// \param lowerSlope derivative of the function at the lower bound, used only if \a lowerSlopeNatural is false
1410  /// \param upperSlopeNatural if true, set y''(last point)=0, otherwise compute it from \a upperSope
1411  /// \param upperSlope derivative of the function at the upper bound, used only if \a upperSlopeNatural is false
1412  /// \param splined if true (default), use cubic spline, if false, use linear interpolation.
1413  /// \return the same interpolating function, filled
1414  interpolating_function_p<float_type> & load_pairs(
1415  std::vector<std::pair<float_type, float_type> > &data,
1416  bool lowerSlopeNatural, float_type lowerSlope,
1417  bool upperSlopeNatural, float_type upperSlope, bool splined=true
1418  ) throw(c2_exception);
1419 
1420  /// \brief do the dirty work of constructing the spline from a function.
1421  /// \param func a function without any requirement of valid derivatives to sample into an interpolating function.
1422  /// Very probably a c2_classic_function.
1423  /// \param amin the lower bound of the region to sample
1424  /// \param amax the upper bound of the region to sample
1425  /// \param abs_tol the maximum absolute error permitted when linearly interpolating the points.
1426  /// the real error will be much smaller, since this uses cubic splines at the end.
1427  /// \param rel_tol the maximum relative error permitted when linearly interpolating the points.
1428  /// the real error will be much smaller, since this uses cubic splines at the end.
1429  /// \param lowerSlopeNatural if true, set y'(first point) from 3-point parabola, otherwise compute it from \a lowerSope
1430  /// \param lowerSlope derivative of the function at the lower bound, used only if \a lowerSlopeNatural is false
1431  /// \param upperSlopeNatural if true, set y'(last point) from 3-point parabola, otherwise compute it from \a upperSope
1432  /// \param upperSlope derivative of the function at the upper bound, used only if \a upperSlopeNatural is false
1433  /// \return the same interpolating function, filled
1434  /// \note If the interpolator being filled has a log vertical axis, put the desired relative error in
1435  /// \a abs_tol, and 0 in \a rel_tol since the absolute error on the log of a function is the relative error
1436  /// on the function itself.
1437  interpolating_function_p<float_type> & sample_function(const c2_function<float_type> &func,
1438  float_type amin, float_type amax, float_type abs_tol, float_type rel_tol,
1439  bool lowerSlopeNatural, float_type lowerSlope,
1440  bool upperSlopeNatural, float_type upperSlope
1441  ) throw(c2_exception);
1442 
1443 
1444  /// \brief initialize from a grid of points and a c2_function (un-normalized) to an
1445  /// interpolator which, when evaluated with a uniform random variate on [0,1] returns random numbers
1446  /// distributed as the input function.
1447  /// \see \ref random_subsec "Arbitrary random generation"
1448  /// inverse_integrated_density starts with a probability density std::vector, generates the integral,
1449  /// and generates an interpolating_function_p of the inverse function which, when evaluated using a uniform random on [0,1] returns values
1450  /// with a density distribution equal to the input distribution
1451  /// If the data are passed in reverse order (large X first), the integral is carried out from the big end.
1452  /// \param bincenters the positions at which to sample the function \a binheights
1453  /// \param binheights a function which describes the density of the random number distribution to be produced.
1454  /// \return an initialized interpolator, which
1455  /// if evaluated randomly with a uniform variate on [0,1] produces numbers
1456  /// distributed according to \a binheights
1458  const std::vector<float_type> &bincenters, const c2_function<float_type> &binheights)
1459  throw(c2_exception);
1460 
1461  /// \brief initialize from a grid of points and an std::vector of probability densities (un-normalized) to an
1462  /// interpolator which, when evaluated with a uniform random variate on [0,1] returns random numbers
1463  /// distributed as the input histogram.
1464  /// \see \ref random_subsec "Arbitrary random generation"
1465  /// inverse_integrated_density starts with a probability density std::vector, generates the integral,
1466  /// and generates an interpolating_function_p of the inverse function which, when evaluated using a uniform random on [0,1] returns values
1467  /// with a density distribution equal to the input distribution
1468  /// If the data are passed in reverse order (large X first), the integral is carried out from the big end.
1469  /// \param bins if \a bins .size()==\a binheights .size(), the centers of the bins. \n
1470  /// if \a bins .size()==\a binheights .size()+1, the edges of the bins
1471  /// \param binheights a vector which describes the density of the random number distribution to be produced.
1472  /// Note density... the numbers in the bins are not counts, but counts/unit bin width.
1473  /// \param splined if true (default), use cubic spline, if false, use linear interpolation.
1474  /// This can often be used to fix ringing if there are very sharp features in the generator.
1475  /// \return an initialized interpolator, which
1476  /// if evaluated randomly with a uniform variate on [0,1] produces numbers
1477  /// distributed according to \a binheights
1479  const std::vector<float_type> &bins, const std::vector<float_type> &binheights, bool splined=true)
1480  throw(c2_exception);
1481 
1482  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception);
1483 
1484  /// \brief destructor
1485  virtual ~interpolating_function_p() { delete &fTransform; }
1486 
1487  /// \brief create a new, empty interpolating function of this type (virtual constructor)
1489  { return *new interpolating_function_p<float_type>(); }
1490 
1491  /// \brief retrieve copies of the x & y tables from which this was built
1492  ///
1493  /// This is often useful in the creation of new interpolating functions with transformed data.
1494  /// The vectors will have their sizes set correctly on return.
1495  /// \param [in, out] xvals the abscissas
1496  /// \param [in, out] yvals the ordinates
1497  void get_data(std::vector<float_type> &xvals, std::vector<float_type> &yvals) const throw() ;
1498 
1499  /// \brief retrieve copies of the transformed x, y and y2 tables from which this was built
1500  ///
1501  /// The vectors will have their sizes set correctly on return.
1502  /// \param [in, out] xvals the transformed abscissas
1503  /// \param [in, out] yvals the transformed ordinates
1504  /// \param [in, out] y2vals the second derivatives
1506  std::vector<float_type> &xvals,
1507  std::vector<float_type> &yvals,
1508  std::vector<float_type> &y2vals) const
1509  { xvals=X; yvals=F; y2vals=y2; }
1510 
1511  /// \brief enable extrapolation of the function below the tabulated data.
1512  ///
1513  /// This allows the interpolator to be extrapolated outside the bounds of the data,
1514  /// using whatever derivatives it already had at the lower bound.
1515  /// \param bound the abscissa to which the function should be extended.
1516  void set_lower_extrapolation(float_type bound);
1517  /// \brief enable extrapolation of the function above the tabulated data.
1518  ///
1519  /// This allows the interpolator to be extrapolated outside the bounds of the data,
1520  /// using whatever derivatives it already had at the upper bound.
1521  /// \param bound the abscissa to which the function should be extended.
1522  void set_upper_extrapolation(float_type bound);
1523 
1524  // these functions correctly combine the interpolating function with another interpolating function
1525  // preserving the X bounds and mapping functions of the host (left hand) function.
1526 
1527  /// \brief create a new interpolating_function_p which is the \a source
1528  /// function applied to every point in the interpolating tables
1529  ///
1530  /// This carefully manages the derivative of the composed function at the two ends.
1531  /// \param source the function to apply
1532  /// \return a new interpolating_function_p with the same mappings for x and y
1534 
1535  /// \brief create a new interpolating_function_p which is the parent interpolating_function_p
1536  /// combined with \a rhs using \a combiner at every point in the interpolating tables
1537  ///
1538  /// This carefully manages the derivative of the composed function at the two ends.
1539  /// \param rhs the function to apply
1540  /// \param combining_stub a function which defines which binary operation to use.
1541  /// \return a new interpolating_function_p with the same mappings for x and y
1543  const c2_binary_function<float_type> *combining_stub
1544  ) const;
1545  /// \brief produce a newly resampled interpolating_function_p which is the specified sum.
1546  /// \param rhs the function to add, pointwise
1547  /// \return a new interpolating_function_p
1549  return binary_operator(rhs, new c2_sum_p<float_type>()); }
1550  /// \brief produce a newly resampled interpolating_function_p which is the specified difference.
1551  /// \param rhs the function to subtract, pointwise
1552  /// \return a new interpolating_function_p
1554  return binary_operator(rhs, new c2_diff_p<float_type>()); }
1555  /// \brief produce a newly resampled interpolating_function_p which is the specified product.
1556  /// \param rhs the function to multiply, pointwise
1557  /// \return a new interpolating_function_p
1559  return binary_operator(rhs, new c2_product_p<float_type>()); }
1560  /// \brief produce a newly resampled interpolating_function_p which is the specified ratio.
1561  /// \param rhs the function to divide, pointwise
1562  /// \return a new interpolating_function_p
1564  return binary_operator(rhs, new c2_ratio_p<float_type>()); }
1565  /// \brief copy data from another interpolating function. This only makes sense if the source
1566  /// function has the same transforms as the destination.
1567  /// \param rhs interpolating_function_p to copy from
1569  Xraw=rhs.Xraw; X=rhs.X; F=rhs.F; y2=rhs.y2;
1571  }
1572 
1574 
1575 protected:
1576  /// \brief create the spline coefficients
1577  void spline(
1578  bool lowerSlopeNatural, float_type lowerSlope,
1579  bool upperSlopeNatural, float_type upperSlope
1580  ) throw(c2_exception);
1581 
1582  // This is for sorting the data. It must be static if it's going to be a class member.
1583  static bool comp_pair(std::pair<float_type,float_type> const &i, std::pair<float_type,float_type> const &j) {return i.first<j.first;}
1584 
1585  std::vector<float_type> Xraw, X, F, y2;
1588  mutable size_t lastKLow;
1589 };
1590 
1591 /// \brief A spline with X transformed into log space.
1592 /// \ingroup interpolators
1593 /// Most useful for functions looking like y=log(x) or any other function with a huge X dynamic range,
1594 /// and a slowly varying Y.
1595 ///
1596 /// The factory function c2_factory::log_lin_interpolating_function() creates *new log_lin_interpolating_function_p()
1597 template <typename float_type=double> class log_lin_interpolating_function_p : public interpolating_function_p <float_type> {
1598 public:
1599  /// \brief an empty log-linear cubic-spline interpolating_function_p
1600  ///
1602  { }
1603  /// \brief create a new, empty interpolating function of this type (virtual constructor)
1606 };
1607 
1608 
1609 /// \brief A spline with Y transformed into log space.
1610 /// \ingroup interpolators
1611 /// Most useful for functions looking like y=exp(x)
1612 ///
1613 /// The factory function c2_factory::lin_log_interpolating_function() creates *new lin_log_interpolating_function_p()
1614 template <typename float_type=double> class lin_log_interpolating_function_p : public interpolating_function_p <float_type> {
1615 public:
1616  /// \brief an empty linear-log cubic-spline interpolating_function_p
1617  ///
1619  { }
1620  /// \brief create a new, empty interpolating function of this type (virtual constructor)
1623 };
1624 
1625 
1626 /// \brief A spline with X and Y transformed into log space.
1627 /// \ingroup interpolators
1628 /// Most useful for functions looking like y=x^n or any other function with a huge X and Y dynamic range.
1629 ///
1630 /// The factory function c2_factory::log_log_interpolating_function() creates *new log_log_interpolating_function_p()
1631 template <typename float_type=double> class log_log_interpolating_function_p : public interpolating_function_p <float_type> {
1632 public:
1633  /// \brief an empty log-log cubic-spline interpolating_function_p
1634  ///
1636  { }
1637  /// \brief create a new, empty interpolating function of this type (virtual constructor)
1640 };
1641 
1642 
1643 /// \brief A spline with X in reciprocal space and Y transformed in log space.
1644 /// \ingroup interpolators
1645 /// Most useful for thermodynamic types of data where Y is roughly A*exp(-B/x).
1646 /// Typical examples are reaction rate data, and thermistor calibration data.
1647 ///
1648 /// The factory function c2_factory::arrhenius_interpolating_function() creates *new arrhenius_interpolating_function_p()
1649 template <typename float_type=double> class arrhenius_interpolating_function_p : public interpolating_function_p <float_type> {
1650 public:
1651  /// \brief an empty arrhenius cubic-spline interpolating_function_p
1652  ///
1654  { }
1655  /// \brief create a new, empty interpolating function of this type (virtual constructor)
1658 };
1659 
1660 /// \brief compute sin(x) with its derivatives.
1661 /// \ingroup math_functions
1662 ///
1663 /// The factory function c2_factory::sin() creates *new c2_sin_p
1664 template <typename float_type=double> class c2_sin_p : public c2_function<float_type> {
1665 public:
1666  /// \brief constructor.
1667  c2_sin_p() : c2_function<float_type>() {}
1668 
1669  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1670  { float_type q=std::sin(x); if(yprime) *yprime=std::cos(x); if(yprime2) *yprime2=-q; return q; }
1671 
1672  /// \brief return a grid dynamically, suitable for use with trig functions with period 2*pi
1673  /// \param amin the lower bound for the grid
1674  /// \param amax upper bound for the grid
1675  /// \param [in, out] grid the sampling grid.
1676  virtual void get_sampling_grid(float_type amin, float_type amax, std::vector<float_type> &grid) const;
1677 };
1678 
1679 /// \brief compute cos(x) with its derivatives.
1680 /// \ingroup math_functions
1681 ///
1682 /// The factory function c2_factory::cos() creates *new c2_cos_p
1683 template <typename float_type=double> class c2_cos_p : public c2_sin_p<float_type> {
1684 public:
1685  /// \brief constructor.
1686  c2_cos_p() : c2_sin_p<float_type>() {}
1687 
1688  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1689  { float_type q=std::cos(x); if(yprime) *yprime=-std::sin(x); if(yprime2) *yprime2=-q; return q; }
1690 };
1691 
1692 /// \brief compute tan(x) with its derivatives.
1693 /// \ingroup math_functions
1694 ///
1695 /// The factory function c2_factory::tan() creates *new c2_tan_p
1696 template <typename float_type=double> class c2_tan_p : public c2_function<float_type> {
1697 public:
1698  /// \brief constructor.
1699  c2_tan_p() : c2_function<float_type>() {}
1700 
1701  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1702  {
1703  float_type c=std::cos(x), ss=std::sin(x);
1704  float_type t=ss/c;
1705  float_type yp=1/(c*c);
1706  if(yprime) *yprime=yp; if(yprime2) *yprime2=2*t*yp;
1707  return t;
1708  }
1709 };
1710 
1711 /// \brief compute log(x) with its derivatives.
1712 /// \ingroup math_functions
1713 ///
1714 /// The factory function c2_factory::log() creates *new c2_log_p
1715 template <typename float_type=double> class c2_log_p : public c2_function<float_type> {
1716 public:
1717  /// \brief constructor.
1718  c2_log_p() : c2_function<float_type>() {}
1719 
1720  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1721  { if(yprime) *yprime=1.0/x; if(yprime2) *yprime2=-1.0/(x*x); return std::log(x); }
1722 };
1723 
1724 /// \brief compute exp(x) with its derivatives.
1725 /// \ingroup math_functions
1726 ///
1727 /// The factory function c2_factory::exp() creates *new c2_exp_p
1728 template <typename float_type=double> class c2_exp_p : public c2_function<float_type> {
1729 public:
1730  /// \brief constructor.
1731  c2_exp_p() : c2_function<float_type>() {}
1732 
1733  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1734  { float_type q=std::exp(x); if(yprime) *yprime=q; if(yprime2) *yprime2=q; return q; }
1735 };
1736 
1737 /// \brief compute sqrt(x) with its derivatives.
1738 /// \ingroup math_functions
1739 ///
1740 /// The factory function c2_factory::sqrt() creates *new c2_sqrt_p()
1741 template <typename float_type=double> class c2_sqrt_p : public c2_function<float_type> {
1742 public:
1743  /// \brief constructor.
1744  c2_sqrt_p() : c2_function<float_type>() {}
1745 
1746  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1747  { float_type q=std::sqrt(x); if(yprime) *yprime=0.5/q; if(yprime2) *yprime2=-0.25/(x*q); return q; }
1748 };
1749 
1750 /// \brief compute scale/x with its derivatives.
1751 /// \ingroup parametric_functions
1752 ///
1753 /// The factory function c2_factory::recip() creates *new c2_recip_p
1754 template <typename float_type=double> class c2_recip_p : public c2_function<float_type> {
1755 public:
1756  /// \brief constructor.
1757  c2_recip_p(float_type scale) : c2_function<float_type>(), rscale(scale) {}
1758 
1759  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1760  {
1761  float_type q=1.0/x;
1762  float_type y=rscale*q;
1763  if(yprime) *yprime=-y*q;
1764  if(yprime2) *yprime2=2*y*q*q;
1765  return y;
1766  }
1767  /// \brief reset the scale factor
1768  /// \param scale the new numerator
1769  void reset(float_type scale) { rscale=scale; }
1770 private:
1771  float_type rscale;
1772 };
1773 
1774 /// \brief compute x with its derivatives.
1775 /// \ingroup math_functions
1776 ///
1777 /// The factory function c2_factory::identity() creates *new c2_identity_p
1778 template <typename float_type=double> class c2_identity_p : public c2_function<float_type> {
1779 public:
1780  /// \brief constructor.
1781  c2_identity_p() : c2_function<float_type>() {}
1782 
1783  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1784  { if(yprime) *yprime=1.0; if(yprime2) *yprime2=0; return x; }
1785 };
1786 
1787 /**
1788  \brief create a linear mapping of another function
1789  \ingroup parametric_functions
1790  for example, given a c2_function \a f
1791  \code
1792  c2_function<double> &F=c2_linear<double>(1.2, 2.0, 3.0)(f);
1793  \endcode
1794  produces a new c2_function F=2.0+3.0*(\a f - 1.2)
1795 
1796  The factory function c2_factory::linear() creates *new c2_linear_p
1797 */
1798 template <typename float_type=double> class c2_linear_p : public c2_function<float_type> {
1799 public:
1800  /// \brief Construct the operator f=y0 + slope * (x-x0)
1801  /// \param x0 the x offset
1802  /// \param y0 the y-intercept i.e. f(x0)
1803  /// \param slope the slope of the mapping
1804  c2_linear_p(float_type x0, float_type y0, float_type slope) :
1805  c2_function<float_type>(), xint(x0), intercept(y0), m(slope) {}
1806  /// \brief Change the slope and intercepts after construction.
1807  /// \param x0 the x offset
1808  /// \param y0 the y-intercept
1809  /// \param slope the slope of the mapping
1810  void reset(float_type x0, float_type y0, float_type slope) { xint=x0; intercept=y0; m=slope; }
1811  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1812  { if(yprime) *yprime=m; if(yprime2) *yprime2=0; return m*(x-xint)+intercept; }
1813 
1814 private:
1815  float_type xint, intercept, m;
1816 protected:
1817  c2_linear_p() {} // do not allow naked construction... it is usually an accident.
1818 };
1819 
1820 /**
1821 \brief create a quadratic mapping of another function
1822  \ingroup parametric_functions
1823  for example, given a c2_function \a f
1824  \code
1825  c2_function<double> &F=c2_quadratic<double>(1.2, 2.0, 3.0, 4.0)(f);
1826  \endcode
1827  produces a new c2_function F=2.0 + 3.0*(f-1.2) + 4.0*(f-1.2)^2
1828 
1829  note that the parameters are overdetermined, but allows the flexibility of two different representations
1830 
1831  The factory function c2_factory::quadratic() creates *new c2_quadratic_p
1832  */
1833 template <typename float_type=double> class c2_quadratic_p : public c2_function<float_type> {
1834 public:
1835  /// \brief Construct the operator
1836  /// \param x0 the center around which the powers are computed
1837  /// \param y0 the value of the function at \a x = \a x0
1838  /// \param xcoef the scale on the (\a x - \a x0) term
1839  /// \param x2coef the scale on the (\a x - \a x0)^2 term
1840  c2_quadratic_p(float_type x0, float_type y0, float_type xcoef, float_type x2coef) :
1841  c2_function<float_type>(), intercept(y0), center(x0), a(x2coef), b(xcoef) {}
1842  /// \brief Modify the coefficients after construction
1843  /// \param x0 the new center around which the powers are computed
1844  /// \param y0 the new value of the function at \a x = \a x0
1845  /// \param xcoef the new scale on the (\a x - \a x0) term
1846  /// \param x2coef the new scale on the (\a x - \a x0)^2 term
1847  void reset(float_type x0, float_type y0, float_type xcoef, float_type x2coef) { intercept=y0; center=x0; a=x2coef; b=xcoef; }
1848  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1849  { float_type dx=x-center; if(yprime) *yprime=2*a*dx+b; if(yprime2) *yprime2=2*a; return a*dx*dx+b*dx+intercept; }
1850 
1851 private:
1852  float_type intercept, center, a, b;
1853 protected:
1854  c2_quadratic_p() {} // do not allow naked construction... it is usually an accident.
1855 };
1856 
1857 /**
1858 \brief create a power law mapping of another function
1859  \ingroup parametric_functions
1860  for example, given a c2_function \a f
1861  \code
1862  c2_power_law_p<double> PLaw(1.2, 2.5);
1863  c2_composed_function_p<double> &F=PLaw(f);
1864  \endcode
1865  produces a new c2_function F=1.2 * f^2.5
1866 
1867  The factory function c2_factory::power_law() creates *new c2_power_law_p
1868  */
1869 template <typename float_type=double> class c2_power_law_p : public c2_function<float_type> {
1870 public:
1871  /// \brief Construct the operator
1872  /// \param scale the multipler
1873  /// \param power the exponent
1874  c2_power_law_p(float_type scale, float_type power) :
1875  c2_function<float_type>(), a(scale), b(power) {}
1876  /// \brief Modify the mapping after construction
1877  /// \param scale the new multipler
1878  /// \param power the new exponent
1879  void reset(float_type scale, float_type power) { a=scale; b=power; }
1880  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception)
1881  { float_type q=a*std::pow(x,b-2); if(yprime) *yprime=b*q*x; if(yprime2) *yprime2=b*(b-1)*q; return q*x*x; }
1882 
1883 private:
1884  float_type a, b;
1885 protected:
1886  c2_power_law_p() {} // do not allow naked construction... it is usually an accident.
1887 };
1888 
1889 /**
1890 \brief create the formal inverse function of another function
1891  \ingroup containers
1892  for example, given a c2_function \a f
1893  \code
1894  c2_inverse_function<double> inv(f);
1895  a=f(x);
1896  x1=inv(a);
1897  \endcode
1898  will return x1=x to machine precision. The important part of this
1899  is that the resulting function is a first-class c2_function, so it knows its
1900  derivatives, too, unlike the case of a simple root-finding inverse. This means
1901  it can be integrated (for example) quite efficiently.
1902 
1903  \see \ref combined_inversion_hinting_sampling
1904 
1905  The factory function c2_factory::inverse_function() creates *new c2_inverse_function_p
1906 */
1907 template <typename float_type=double> class c2_inverse_function_p : public c2_function<float_type> {
1908 public:
1909  /// \brief Construct the operator
1910  /// \param source the function to be inverted
1912  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw(c2_exception);
1913 
1914  /// \brief give the function a hint as to where to look for its inverse
1915  /// \param hint the likely value of the inverse, which defaults to whatever the evaluation returned.
1916  void set_start_hint(float_type hint) const { start_hint=hint; }
1917 
1918  /// \brief get the starting hint.
1919  ///
1920  /// This is virtual so if there is a better way, this can be easily overridden.
1921  /// It is used in value_with_derivatives() to guess where to start the root finder.
1922  /// \param x the abscissa for which an estimate is needed
1923  virtual float_type get_start_hint(float_type x) const
1924  { return hinting_function.valid()? hinting_function(x) : start_hint; }
1925 
1926  /// \brief set or unset the approximate function used to start the root finder
1927  /// \anchor set_hinting_function_discussion
1928  /// A hinting function is mostly useful if the evaluation of this inverse is
1929  /// going to be carried out in very non-local order, so the root finder has to start over
1930  /// for each step. If most evaluations are going to be made in fairly localized clusters (scanning
1931  /// through the function, for example), the default mechanism used (which just remembers the last point)
1932  /// is almost certainly faster.
1933  ///
1934  /// Typically, the hinting function is likely to be set up by creating the inverse function,
1935  /// and then adaptively sampling an interpolating function from it, and then using the result
1936  /// to hint it. Another way, if the parent function is already an interpolating function, is just to create
1937  /// a version of the parent with the x & y coordinates reversed.
1938  ///
1939  /// \see \ref combined_inversion_hinting_sampling
1940  ///
1941  /// \param hint_func the function that is an approximate inverse of the parent of this inverse_function
1943  { hinting_function.set_function(hint_func); }
1944  /// \brief set the hinting function from a pointer.
1945  ///
1946  /// See \ref set_hinting_function_discussion "discussion"
1947  /// \param hint_func the container holding the function
1949  { hinting_function=hint_func; }
1950 
1951 protected:
1952  c2_inverse_function_p() {} // do not allow naked construction... it is usually an accident.
1953  mutable float_type start_hint;
1956 };
1957 
1958 /**
1959  \brief
1960  An interpolating_function_p which is the cumulative integral of a histogram.
1961  \ingroup interpolators
1962  Note than binedges should be one element longer than binheights, since the lower & upper edges are specified.
1963  Note that this is a malformed spline, since the second derivatives are all zero, so it has less continuity.
1964  Also, note that the bin edges can be given in backwards order to generate the
1965  reversed accumulation (starting at the high end)
1966 */
1967 
1968 template <typename float_type=double> class accumulated_histogram : public interpolating_function_p <float_type> {
1969 public:
1970  /// \brief Construct the integrated histogram
1971  /// \param binedges the edges of the bins in \a binheights. It should have one more element than \a binheights
1972  /// \param binheights the number of counts in each bin.
1973  /// \param normalize if true, normalize integral to 1
1974  /// \param inverse_function if true, drop zero channels, and return inverse function for random generation
1975  /// \param drop_zeros eliminate null bins before integrating, so integral is strictly monotonic.
1976  accumulated_histogram(const std::vector<float_type>binedges, const std::vector<float_type> binheights,
1977  bool normalize=false, bool inverse_function=false, bool drop_zeros=true);
1978 
1979 };
1980 
1981 /// \brief create a c2_function which smoothly connects two other c2_functions.
1982 /// \ingroup parametric_functions
1983 /// This takes two points and generates a polynomial which matches two c2_function arguments
1984 /// at those two points, with two derivatives at each point, and an arbitrary value at the center of the
1985 /// region. It is useful for splicing together functions over rough spots (0/0, for example).
1986 ///
1987 /// If \a auto_center is true, the value at the midpoint is computed so that the resulting polynomial is
1988 /// of order 5. If \a auto_center is false, the value \a y1 is used at the midpoint, resulting in a
1989 /// polynomial of order 6.
1990 ///
1991 /// This is usually used in conjunction with c2_piecewise_function_p to assemble an apparently seamless
1992 /// function from a series of segments.
1993 /// \see \ref piecewise_applications_subsec "Sample Applications" and \ref c2_function::adaptively_sample() "Adaptive sampling"
1994 ///
1995 /// The factory function c2_factory::connector_function() creates *new c2_connector_function_p
1996 template <typename float_type=double> class c2_connector_function_p : public c2_function<float_type> {
1997 public:
1998  /// \brief construct the container from two functions
1999  /// \param x0 the point at which to match \a f1 and its derivatives
2000  /// \param f0 the function on the left side to be connected
2001  /// \param x2 the point at which to match \a f2 and its derivatives
2002  /// \param f2 the function on the right side to be connected
2003  /// \param auto_center if true, no midpoint value is specified. If false, match the value \a y1 at the midpoint
2004  /// \param y1 the value to match at the midpoint, if \a auto_center is false
2005  /// \return a c2_function with domain (\a x0,\a x2) which smoothly connects \a f0(x0) and \a f2(x2)
2006  c2_connector_function_p(float_type x0, const c2_function<float_type> &f0, float_type x2, const c2_function<float_type> &f2,
2007  bool auto_center, float_type y1);
2008  /// \brief construct the container from numerical values
2009  /// \param x0 the position of the left edge
2010  /// \param y0 the function derivative on the left boundary
2011  /// \param yp0 the function second derivative on the left boundary
2012  /// \param ypp0 the function value on the left boundary
2013  /// \param x2 the position of the right edge
2014  /// \param y2 the function derivative on the right boundary
2015  /// \param yp2 the function second derivative on the right boundary
2016  /// \param ypp2 the function value on the right boundary
2017  /// \param auto_center if true, no midpoint value is specified. If false, match the value \a y1 at the midpoint
2018  /// \param y1 the value to match at the midpoint, if \a auto_center is false
2019  /// \return a c2_function with domain (\a x0,\a x2) which smoothly connects the points described
2020  /// \anchor c2_connector_raw_init_docs
2022  float_type x0, float_type y0, float_type yp0, float_type ypp0,
2023  float_type x2, float_type y2, float_type yp2, float_type ypp2,
2024  bool auto_center, float_type y1);
2025  /// \brief construct the container from c2_fblock<float_type> objects
2026  /// \param fb0 the left edge
2027  /// \param fb2 the right edge
2028  /// \param auto_center if true, no midpoint value is specified. If false, match the value \a y1 at the midpoint
2029  /// \param y1 the value to match at the midpoint, if \a auto_center is false
2030  /// \return a c2_function with domain (\a fb0.x,\a fb2.x) which smoothly connects \a fb0 and \a fb2
2032  const c2_fblock<float_type> &fb0,
2033  const c2_fblock<float_type> &fb2,
2034  bool auto_center, float_type y1);
2035 
2036  /// \brief destructor
2037  virtual ~c2_connector_function_p();
2038  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw (c2_exception);
2039 protected:
2040  /// \brief fill container numerically
2041  void init(
2042  const c2_fblock<float_type> &fb0,
2043  const c2_fblock<float_type> &fb2,
2044  bool auto_center, float_type y1);
2045 
2046  float_type fhinv, fy1, fa, fb, fc, fd, fe, ff;
2047 };
2048 
2049 /// \brief create a c2_function which is a piecewise assembly of other c2_functions.
2050 /// \ingroup containers
2051 /// The functions must have increasing, non-overlapping domains. Any empty space
2052 /// between functions will be filled with a linear interpolation.
2053 ///
2054 /// \note If you want a smooth connection, instead of the default linear interpolation,
2055 /// create a c2_connector_function_p to bridge the gap. The linear interpolation is intended
2056 /// to be a barely intelligent bridge, and may never get used by anyone.
2057 ///
2058 /// \note The creation of the container results in the creation of an explicit sampling grid.
2059 /// If this is used with functions with a large domain, or which generate very dense sampling grids,
2060 /// it could eat a lot of memory. Do not abuse this by using functions which can generate gigantic grids.
2061 ///
2062 /// \see \ref piecewise_applications_subsec "Sample Applications" \n
2063 /// c2_plugin_function_p page \n
2064 /// c2_connector_function_p page \n
2065 /// \ref c2_function::adaptively_sample() "Adaptive sampling"
2066 ///
2067 /// The factory function c2_factory::piecewise_function() creates *new c2_piecewise_function_p
2068 template <typename float_type=double> class c2_piecewise_function_p : public c2_function<float_type> {
2069 public:
2070  /// \brief construct the container
2072  /// \brief destructor
2073  virtual ~c2_piecewise_function_p();
2074  virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const throw (c2_exception);
2075  /// \brief append a new function to the sequence
2076  ///
2077  /// This takes a c2_function, and appends it onto the end of the piecewise collection.
2078  /// The domain of the function (which MUST be set) specifies the place it will be used in
2079  /// the final function. If the domain exactly abuts the domain of the previous function, it
2080  /// will be directly attached. If there is a gap, the gap will be filled in by linear interpolation.
2081  /// \param func a c2_function with a defined domain to be appended to the collection
2082  void append_function(const c2_function<float_type> &func) throw (c2_exception);
2083 protected:
2084  std::vector<c2_const_ptr<float_type> > functions;
2085  mutable int lastKLow;
2086 };
2087 
2088 #include "c2_function.icc"
2089 
2090 #endif
const c2_function< float_type > * func
Definition: c2_function.hh:674
c2_recip_p(float_type scale)
constructor.
interpolating_function_p< float_type > & load_random_generator_bins(const std::vector< float_type > &bins, const std::vector< float_type > &binheights, bool splined=true)
initialize from a grid of points and an std::vector of probability densities (un-normalized) to an in...
create a non-generic container for a c2_function which handles the reference counting.
Definition: c2_function.hh:725
static float_type zero(float_type)
utility function f(x)=0 useful in axis transforms
float_type operator()(float_type x) const
evaluate the function in the classic way, ignoring derivatives.
Definition: c2_function.hh:176
const std::string cvs_file_vers() const
get versioning information for the source file
void reset(float_type scale)
set a new scale factor
Definition: c2_function.hh:910
c2_typed_ptr(c2_class< float_type > &f)
construct the container with a pre-defined function
Definition: c2_function.hh:731
bool ypbad
flag, filled in by c2_function::fill_fblock(), indicating the derivative is NaN of Inf ...
Definition: c2_function.hh:113
bool stub
if true, we don't own any functions, we are just a source of a combining function.
Definition: c2_function.hh:891
log_log_interpolating_function_p()
an empty log-log cubic-spline interpolating_function_p
c2_exp_p()
constructor.
virtual float_type fInDPrime(float_type x) const
virtual input X transform second derivative
c2_plugin_function_p(c2_function< float_type > &f)
construct the container with a pre-defined function
Definition: c2_function.hh:788
virtual float_type fOut(float_type x) const
virtual output X transform
std::vector< c2_const_ptr< float_type > > functions
c2_linear_p(float_type x0, float_type y0, float_type slope)
Construct the operator f=y0 + slope * (x-x0)
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
#define inline
Definition: internal.h:71
void set_start_hint(float_type hint) const
give the function a hint as to where to look for its inverse
interpolating_function_p< float_type > & multiply_pointwise(const c2_function< float_type > &rhs) const
produce a newly resampled interpolating_function_p which is the specified product.
c2_function< float_type > * get_ptr() const
get an unchecked pointer to our owned function
Definition: c2_function.hh:698
c2_log_p()
constructor.
create a c2_function which is the ratio of two other c2_functions.This should always be constructed u...
Definition: c2_function.hh:83
virtual void get_sampling_grid(float_type amin, float_type amax, std::vector< float_type > &grid) const
return a grid dynamically, suitable for use with trig functions with period 2*pi
void get_internal_data(std::vector< float_type > &xvals, std::vector< float_type > &yvals, std::vector< float_type > &y2vals) const
retrieve copies of the transformed x, y and y2 tables from which this was built
const bool fTransformed
flag to indicate if this transform is not the identity
float_type x
the abscissa
Definition: c2_function.hh:105
a transformation of a function in and out of Arrhenius (1/x vs. log(y)) space
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
function to manage the binary operation, used by c2_binary_function::value_with_derivatives() ...
Definition: c2_function.hh:850
Create a very lightweight method to return a scalar multiple of another function. \ \The factory func...
Definition: c2_function.hh:899
a container into which any conventional c-style function can be dropped, to create a degenerate c2_fu...
Definition: c2_function.hh:534
const c2_const_ptr< float_type > Left
Definition: c2_function.hh:889
the identity transform
compute tan(x) with its derivatives.The factory function c2_factory::tan() creates *new c2_tan_p ...
interpolating_function_p< float_type > & add_pointwise(const c2_function< float_type > &rhs) const
produce a newly resampled interpolating_function_p which is the specified sum.
float_type fXMax
Definition: c2_function.hh:449
float_type ypp
the second derivative at x
Definition: c2_function.hh:111
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
Provides function composition (nesting)This allows evaluation of f(g(x)) where f and g are c2_functio...
Definition: c2_function.hh:79
create the formal inverse function of another functionfor example, given a c2_function f ...
c2_ratio_p()
Create a stub just for the combiner to avoid statics.
c2_const_plugin_function_p(const c2_function< float_type > &f)
construct the container with a pre-defined function
Definition: c2_function.hh:829
c2_const_ptr< float_type > sampler_function
log axis transform
void claim_ownership() const
increment our reference count. Destruction is only legal if the count is zero.
Definition: c2_function.hh:406
std::vector< float_type > y2
virtual ~c2_piecewise_function_p()
destructor
size_t get_evaluations() const
this is a counter owned by the function but which can be used to monitor efficiency of algorithms...
Definition: c2_function.hh:307
virtual float_type value_with_derivatives(float_type, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
create a linear mapping of another functionfor example, given a c2_function f
virtual ~c2_connector_function_p()
destructor
c2_ratio_p< float_type > & operator/(const c2_function< float_type > &rhs) const
factory function to create a c2_ratio_p from a regular algebraic expression.
Definition: c2_function.hh:391
const c2_const_ptr< float_type > Right
Definition: c2_function.hh:889
void reset(float_type x0, float_type y0, float_type xcoef, float_type x2coef)
Modify the coefficients after construction.
A spline with X in reciprocal space and Y transformed in log space.Most useful for thermodynamic type...
void preen_sampling_grid(std::vector< float_type > *result) const
clean up endpoints on a grid of points
interpolating_function_p< float_type > & divide_pointwise(const c2_function< float_type > &rhs) const
produce a newly resampled interpolating_function_p which is the specified ratio.
float_type(*const combine)(const c2_function< float_type > &left, const c2_function< float_type > &right, float_type x, float_type *yprime, float_type *yprime2)
Definition: c2_function.hh:885
c2_function< float_type > & normalized_function(float_type amin, float_type amax, float_type norm=1.0) const
create a new c2_function from this one which is normalized on the interval
c2_diff_p< float_type > & operator-(const c2_function< float_type > &rhs) const
factory function to create a c2_diff_p from a regular algebraic expression.
Definition: c2_function.hh:655
virtual const char * what() const
Definition: c2_function.hh:73
c2_function< float_type > & square_normalized_function(float_type amin, float_type amax, float_type norm=1.0) const
create a new c2_function from this one which is square-normalized on the interval ...
c2_cos_p()
constructor.
const c2_ptr< float_type > & operator=(const c2_ptr< float_type > &f)
fill the container from another container
Definition: c2_function.hh:705
interpolating_function_p< float_type > & subtract_pointwise(const c2_function< float_type > &rhs) const
produce a newly resampled interpolating_function_p which is the specified difference.
const c2_transformation< float_type > & Y
the Y axis transform
c2_const_plugin_function_p()
construct the container with no function
Definition: c2_function.hh:827
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
static float_type recip(float_type x)
utility function f(x)=1/x useful in axis transforms
interpolating_function_p(const c2_function_transformation< float_type > &transform)
an empty cubic-spline interpolating_function_p with a specific transform
a transformation of a function in and out of log-log space
interpolating_function_p< float_type > & load(const std::vector< float_type > &x, const std::vector< float_type > &f, bool lowerSlopeNatural, float_type lowerSlope, bool upperSlopeNatural, float_type upperSlope, bool splined=true)
do the dirty work of constructing the spline from a function.
float_type get_trouble_point() const
Find out where a calculation ran into trouble, if it got a nan. If the most recent computation did no...
Definition: c2_function.hh:403
c2_diff_p< float_type > & operator-(const c2_function< float_type > &rhs) const
factory function to create a c2_diff_p from a regular algebraic expression.
Definition: c2_function.hh:381
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
virtual void get_sampling_grid(float_type amin, float_type amax, std::vector< float_type > &grid) const
return the grid of 'interesting' points along this function which lie in the region requested ...
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
Definition: c2_function.hh:799
c2_cached_function_p(const c2_function< float_type > &f)
construct the container
Definition: c2_function.hh:942
a transformation of a function in and out of lin-lin space
#define c2_isnan
Definition: c2_function.hh:51
void append_function(const c2_function< float_type > &func)
append a new function to the sequence
compute cos(x) with its derivatives.The factory function c2_factory::cos() creates *new c2_cos_p ...
c2_classic_function_p(const float_type(*c_func)(float_type))
construct the container
Definition: c2_function.hh:538
size_t evaluations
Definition: c2_function.hh:450
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
compute exp(x) with its derivatives.The factory function c2_factory::exp() creates *new c2_exp_p ...
create a c2_function which is the product of two other c2_functions.This should always be constructed...
Definition: c2_function.hh:82
virtual void get_sampling_grid(float_type amin, float_type amax, std::vector< float_type > &grid) const
return the grid of 'interesting' points along this function which lie in the region requested ...
Definition: c2_function.hh:810
const c2_function< float_type > * operator->() const
get a checked pointer to our owned function
Definition: c2_function.hh:622
static float_type ident(float_type x)
utility function f(x)=x useful in axis transforms
float_type(*const pOut)(float_type)
non-virtual pointer to output X transform
void set_function(const c2_function< float_type > *f)
fill the container with a new function, or clear it with a null pointer
Definition: c2_function.hh:832
reciprocal axis transform
c2_class< float_type > * operator->() const
get a checked pointer to our owned function
Definition: c2_function.hh:744
c2_constant_p(float_type x)
float_type partial_integrals(std::vector< float_type > xgrid, std::vector< float_type > *partials=0, float_type abs_tol=1e-12, float_type rel_tol=1e-12, int derivs=2, bool adapt=true, bool extrapolate=true) const
for points in xgrid, adaptively return Integral[f(x),{x,xgrid[i],xgrid[i+1]}] and return in vector...
const c2_transformation< float_type > & X
the X axis transform
virtual float_type evaluate(float_type xraw, float_type y, float_type yp0, float_type ypp0, float_type *yprime, float_type *yprime2) const
evaluate the transformation from internal coordinates to external coordinates
c2_sum_p()
Create a stub just for the combiner to avoid statics.
c2_piecewise_function_p()
construct the container
interpolating_function_p< float_type > & load_pairs(std::vector< std::pair< float_type, float_type > > &data, bool lowerSlopeNatural, float_type lowerSlope, bool upperSlopeNatural, float_type upperSlope, bool splined=true)
do the dirty work of constructing the spline from a function.
c2_product_p< float_type > & operator*(const c2_function< float_type > &rhs) const
factory function to create a c2_product_p from a regular algebraic expression.
Definition: c2_function.hh:660
std::vector< float_type > F
c2_plugin_function_p()
construct the container with no function
Definition: c2_function.hh:786
static float_type one(float_type)
utility function f(x)=1 useful in axis transforms
bool check_monotonicity(const std::vector< float_type > &data, const char message[]) const
check that a vector is monotonic, throw an exception if not, and return a flag if it is reversed ...
virtual ~c2_const_plugin_function_p()
destructor
Definition: c2_function.hh:835
virtual ~c2_transformation()
the destructor
c2_function(const c2_function< float_type > &src)
Definition: c2_function.hh:428
interpolating_function_p< float_type > & binary_operator(const c2_function< float_type > &rhs, const c2_binary_function< float_type > *combining_stub) const
create a new interpolating_function_p which is the parent interpolating_function_p combined with rhs ...
std::vector< float_type > * get_sampling_grid_pointer() const
get the sampling grid, which may be a null pointer
Definition: c2_function.hh:332
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
virtual interpolating_function_p< float_type > & clone() const
create a new, empty interpolating function of this type (virtual constructor)
void set_hinting_function(const c2_function< float_type > *hint_func)
set or unset the approximate function used to start the root finder
float_type(*const pInDPrime)(float_type)
non-virtual pointer to input X transform second derivative
virtual void set_sampling_grid(const std::vector< float_type > &grid)
establish a grid of 'interesting' points on the function.
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
size_t release_ownership_for_return() const
decrement our reference count. Do not destroy at zero.
Definition: c2_function.hh:409
A spline with X and Y transformed into log space.Most useful for functions looking like y=x^n or any ...
c2_function_transformation(const c2_transformation< float_type > &xx, const c2_transformation< float_type > &yy)
construct this from two c2_transformation instances
compute scale/x with its derivatives.The factory function c2_factory::recip() creates *new c2_recip_p...
void spline(bool lowerSlopeNatural, float_type lowerSlope, bool upperSlopeNatural, float_type upperSlope)
create the spline coefficients
void reset(float_type scale)
reset the scale factor
static float_type combine(const c2_function< float_type > &left, const c2_function< float_type > &right, float_type x, float_type *yprime, float_type *yprime2)
execute math necessary to do multiplication
const c2_function_transformation< float_type > & fTransform
create a cubic spline interpolation of a set of (x,y) pairsThis is one of the main reasons for c2_fun...
c2_composed_function_p()
Create a stub just for the combiner to avoid statics.
Definition: c2_function.hh:983
log_lin_interpolating_function_p()
an empty log-linear cubic-spline interpolating_function_p
bool valid() const
check if we have a valid function
Definition: c2_function.hh:625
A spline with X transformed into log space.Most useful for functions looking like y=log(x) or any oth...
void set_hinting_function(const c2_const_ptr< float_type > hint_func)
set the hinting function from a pointer.
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
Definition: c2_function.hh:915
A container into which any other c2_function can be dropped.It allows a function to be pre-evaluated ...
Definition: c2_function.hh:937
c2_quadratic_p(float_type x0, float_type y0, float_type xcoef, float_type x2coef)
Construct the operator.
structure used to hold evaluated function data at a point.
Definition: c2_function.hh:101
void reset(float_type x0, float_type y0, float_type slope)
Change the slope and intercepts after construction.
a c2_plugin_function_p which promises not to fiddle with the plugged function.The factory function c2...
Definition: c2_function.hh:824
c2_piecewise_function_p< float_type > * adaptively_sample(float_type amin, float_type amax, float_type abs_tol=1e-12, float_type rel_tol=1e-12, int derivs=2, std::vector< float_type > *xvals=0, std::vector< float_type > *yvals=0) const
create a c2_piecewise_function_p from c2_connector_function_p segments which is a representation of t...
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
c2_function< float_type > * operator->() const
get a checked pointer to our owned function
Definition: c2_function.hh:701
void set_function(c2_function< float_type > *f)
fill the container with a new function, or clear it with a null pointer and copy our domain ...
Definition: c2_function.hh:792
compute x with its derivatives.The factory function c2_factory::identity() creates *new c2_identity_p...
void set_upper_extrapolation(float_type bound)
enable extrapolation of the function above the tabulated data.
c2_const_ptr(const c2_const_ptr< float_type > &src)
copy constructor
Definition: c2_function.hh:576
a transformation of a coordinate, including an inverse
float_type(*const pIn)(float_type)
non-virtual pointer to input X transform
c2_transformation(bool transformed, float_type(*xin)(float_type), float_type(*xinp)(float_type), float_type(*xinpp)(float_type), float_type(*xout)(float_type))
initialize all our function pointers
virtual ~c2_transformation_recip()
destructor
void set_domain(float_type amin, float_type amax)
set the domain for this function.
Definition: c2_function.hh:301
virtual interpolating_function_p< float_type > & clone() const
create a new, empty interpolating function of this type (virtual constructor)
c2_binary_function(float_type(*combiner)(const c2_function< float_type > &left, const c2_function< float_type > &right, float_type x, float_type *yprime, float_type *yprime2), const c2_function< float_type > &left, const c2_function< float_type > &right)
construct the binary function
Definition: c2_function.hh:865
virtual ~c2_function()
destructor
Definition: c2_function.hh:152
virtual ~c2_transformation_log()
destructor
a transformation of a function in and out of lin-log space
c2_tan_p()
constructor.
void release_for_return()
release the function without destroying it, so it can be returned from a function ...
Definition: c2_function.hh:601
std::vector< float_type > X
const c2_const_ptr< float_type > func
Definition: c2_function.hh:962
c2_diff_p(const c2_function< float_type > &left, const c2_function< float_type > &right)
construct left - right
const c2_const_ptr< float_type > & operator=(const c2_const_ptr< float_type > &f)
fill the container from another container
Definition: c2_function.hh:589
a container into which any other c2_function can be dropped, to allow expressions with replacable com...
Definition: c2_function.hh:782
const c2_function< float_type > & get() const
get a reference to our owned function
Definition: c2_function.hh:614
float_type yp
the derivative at x
Definition: c2_function.hh:109
static float_type combine(const c2_function< float_type > &left, const c2_function< float_type > &right, float_type x, float_type *yprime, float_type *yprime2)
execute math necessary to do subtraction
create a container for a c2_function which handles the reference counting.It is useful as a smart con...
Definition: c2_function.hh:566
static float_type combine(const c2_function< float_type > &left, const c2_function< float_type > &right, float_type x, float_type *yprime, float_type *yprime2)
execute math necessary to do division
c2_transformation(bool transformed)
initialize all our function pointers so that only the (overridden) virtual functions can be called wi...
c2_transformation_log()
constructor
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
c2_transformation_linear()
constructor
std::vector< float_type > Xraw
const c2_const_ptr< float_type > func
const float_type(* func)(float_type)
pointer to our function
Definition: c2_function.hh:553
interpolating_function_p< float_type > & sample_function(const c2_function< float_type > &func, float_type amin, float_type amax, float_type abs_tol, float_type rel_tol, bool lowerSlopeNatural, float_type lowerSlope, bool upperSlopeNatural, float_type upperSlope)
do the dirty work of constructing the spline from a function.
const std::string cvs_header_vers() const
get versioning information for the header file
Definition: c2_function.hh:142
the exception class for c2_function operations.
Definition: c2_function.hh:65
float_type integral(float_type amin, float_type amax, std::vector< float_type > *partials=0, float_type abs_tol=1e-12, float_type rel_tol=1e-12, int derivs=2, bool adapt=true, bool extrapolate=true) const
a fully-automated integrator which uses the information provided by the get_sampling_grid() function ...
c2_ptr()
construct the container with no function
Definition: c2_function.hh:685
float_type find_root(float_type lower_bracket, float_type upper_bracket, float_type start, float_type value, int *error=0, float_type *final_yprime=0, float_type *final_yprime2=0) const
solve f(x)==value very efficiently, with explicit knowledge of derivatives of the function ...
const c2_function< float_type > * get_ptr() const
get an unchecked pointer to our owned function
Definition: c2_function.hh:620
create a c2_function which smoothly connects two other c2_functions.This takes two points and generat...
static float_type recip_prime(float_type x)
utility function f(x)=-1/x**2 useful in axis transforms
a transformation of a function in and out of a coordinate space, using 2 c2_transformations ...
c2_product_p()
Create a stub just for the combiner to avoid statics.
bool yppbad
flag, filled in by c2_function::fill_fblock(), indicating the second derivative is NaN of Inf ...
Definition: c2_function.hh:115
c2_transformation_recip()
constructor
the parent class for all c2_functions.c2_functions know their value, first, and second derivative at ...
Definition: c2_function.hh:138
arrhenius_interpolating_function_p()
an empty arrhenius cubic-spline interpolating_function_p
c2_sin_p()
constructor.
c2_product_p(const c2_function< float_type > &left, const c2_function< float_type > &right)
construct left * right
virtual ~c2_transformation_linear()
destructor
c2_typed_ptr(const c2_typed_ptr< float_type, c2_class > &src)
copy constructor
Definition: c2_function.hh:735
T max(const T t1, const T t2)
brief Return the largest of the two arguments
c2_ptr(c2_function< float_type > &f)
construct the container with a pre-defined function
Definition: c2_function.hh:688
c2_const_ptr(const c2_function< float_type > &f)
construct the container with a pre-defined function
Definition: c2_function.hh:572
c2_exception(const char msgcode[])
construct the exception with an error message
Definition: c2_function.hh:69
c2_ratio_p(const c2_function< float_type > &left, const c2_function< float_type > &right)
construct left / right
static float_type recip_prime2(float_type x)
utility function f(x)=2/x**3 useful in axis transforms
float_type fXMin
Definition: c2_function.hh:449
interpolating_function_p< float_type > & unary_operator(const c2_function< float_type > &source) const
create a new interpolating_function_p which is the source function applied to every point in the inte...
c2_class< float_type > * get_ptr() const
get an unchecked pointer to our owned function
Definition: c2_function.hh:747
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
Definition: c2_function.hh:542
A spline with Y transformed into log space.Most useful for functions looking like y=exp(x) ...
c2_product_p< float_type > & operator*(const c2_function< float_type > &rhs) const
factory function to create a c2_product_p from a regular algebraic expression.
Definition: c2_function.hh:386
c2_power_law_p(float_type scale, float_type power)
Construct the operator.
c2_composed_function_p(const c2_function< float_type > &outer, const c2_function< float_type > &inner)
construct outer( inner (x))
Definition: c2_function.hh:980
float_type(*const pInPrime)(float_type)
non-virtual pointer to input X transform derivative
const bool isIdentity
flag indicating of the transform is the identity, and can be skipped for efficiency ...
float_type operator()(float_type x, float_type *yprime, float_type *yprime2) const
convenience operator to make us look like a function
Definition: c2_function.hh:645
c2_binary_function(float_type(*combiner)(const c2_function< float_type > &left, const c2_function< float_type > &right, float_type x, float_type *yprime, float_type *yprime2))
construct a 'stub' c2_binary_function, which provides access to the combine() function ...
Definition: c2_function.hh:879
c2_scaled_function_p(const c2_function< float_type > &outer, float_type scale)
construct the function with its scale factor.
Definition: c2_function.hh:905
virtual interpolating_function_p< float_type > & clone() const
create a new, empty interpolating function of this type (virtual constructor)
~c2_const_ptr()
destructor
Definition: c2_function.hh:611
c2_sum_p(const c2_function< float_type > &left, const c2_function< float_type > &right)
construct left + right
c2_ptr< float_type > func
Definition: c2_function.hh:817
virtual ~c2_function_transformation()
destructor
const c2_const_ptr< float_type > func
the scaling factor for the function
Definition: c2_function.hh:926
c2_sum_p< float_type > & operator+(const c2_function< float_type > &rhs) const
factory function to create a c2_sum_p from a regular algebraic expression.
Definition: c2_function.hh:650
c2_const_ptr()
construct the container with no function
Definition: c2_function.hh:569
An interpolating_function_p which is the cumulative integral of a histogram.Note than binedges should...
virtual float_type get_start_hint(float_type x) const
get the starting hint.
size_t count_owners() const
get the reference count, mostly for debugging
Definition: c2_function.hh:425
virtual interpolating_function_p< float_type > & clone() const
create a new, empty interpolating function of this type (virtual constructor)
lin_log_interpolating_function_p()
an empty linear-log cubic-spline interpolating_function_p
float_type y
the value of the function at x
Definition: c2_function.hh:107
c2_connector_function_p(float_type x0, const c2_function< float_type > &f0, float_type x2, const c2_function< float_type > &f2, bool auto_center, float_type y1)
construct the container from two functions
const XML_Char XML_Encoding * info
virtual ~c2_exception()
Definition: c2_function.hh:70
const XML_Char int const XML_Char * value
void operator=(const c2_typed_ptr< float_type, c2_class > &f)
fill the container from another container
Definition: c2_function.hh:753
accumulated_histogram(const std::vector< float_type >binedges, const std::vector< float_type > binheights, bool normalize=false, bool inverse_function=false, bool drop_zeros=true)
Construct the integrated histogram.
static float_type combine(const c2_function< float_type > &left, const c2_function< float_type > &right, float_type x, float_type *yprime, float_type *yprime2)
execute math necessary to do addition
void reset(float_type val)
virtual ~c2_plugin_function_p()
destructor
Definition: c2_function.hh:805
void refine_sampling_grid(std::vector< float_type > &grid, size_t refinement) const
refine a grid by splitting each interval into more intervals
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
Definition: c2_function.hh:948
std::vector< float_type > * sampling_grid
Definition: c2_function.hh:446
compute sqrt(x) with its derivatives.The factory function c2_factory::sqrt() creates *new c2_sqrt_p()...
interpolating_function_p< float_type > & load_random_generator_function(const std::vector< float_type > &bincenters, const c2_function< float_type > &binheights)
initialize from a grid of points and a c2_function (un-normalized) to an interpolator which...
create a container for a c2_function which handles the reference counting.
Definition: c2_function.hh:86
virtual ~c2_binary_function()
destructor releases ownership of member functions
Definition: c2_function.hh:858
c2_ratio_p< float_type > & operator/(const c2_function< float_type > &rhs) const
factory function to create a c2_ratio_p from a regular algebraic expression.
Definition: c2_function.hh:665
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
#define c2_isfinite
Definition: c2_function.hh:52
void increment_evaluations() const
count evaluations
Definition: c2_function.hh:311
float_type bad_x_point
this point may be used to record where a calculation ran into trouble
Definition: c2_function.hh:452
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const =0
get the value and derivatives.
bool no_overwrite_grid
Definition: c2_function.hh:447
static float_type report_error(float_type x)
utility function for unimplemented conversion
void release_ownership() const
decrement our reference count. If the count reaches zero, destroy ourself.
Definition: c2_function.hh:420
create a c2_function which is the sum of two other c2_function objects.This should always be construc...
Definition: c2_function.hh:80
void clone_data(const interpolating_function_p< float_type > &rhs)
copy data from another interpolating function. This only makes sense if the source ...
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
void fill_fblock(c2_fblock< float_type > &fb) const
fill in a c2_fblock... a shortcut for the integrator & sampler
Definition: c2_function.hh:456
c2_diff_p()
Create a stub just for the combiner to avoid statics.
float_type operator()(float_type x) const
convenience operator to make us look like a function
Definition: c2_function.hh:636
virtual interpolating_function_p< float_type > & clone() const
create a new, empty interpolating function of this type (virtual constructor)
create a quadratic mapping of another functionfor example, given a c2_function f
Definition: c2_function.hh:85
float_type xmax() const
return the upper bound of the domain for this function as set by set_domain()
Definition: c2_function.hh:299
void reset_evaluations() const
reset the counter
Definition: c2_function.hh:309
void set_lower_extrapolation(float_type bound)
enable extrapolation of the function below the tabulated data.
c2_identity_p()
constructor.
void reset(float_type scale, float_type power)
Modify the mapping after construction.
virtual float_type fIn(float_type x) const
virtual input X transform
compute sin(x) with its derivatives.The factory function c2_factory::sin() creates *new c2_sin_p ...
interpolating_function_p()
an empty linear-linear cubic-spline interpolating_function_p
void unset_function()
clear our function
Definition: c2_function.hh:808
c2_sqrt_p()
constructor.
c2_const_ptr< float_type > hinting_function
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
float_type xmin() const
return the lower bound of the domain for this function as set by set_domain()
Definition: c2_function.hh:297
virtual float_type fInPrime(float_type x) const
virtual input X transform derivative
void init(const c2_fblock< float_type > &fb0, const c2_fblock< float_type > &fb2, bool auto_center, float_type y1)
fill container numerically
const XML_Char const XML_Char * data
virtual ~c2_classic_function_p()
Definition: c2_function.hh:549
static float_type combine(const c2_function< float_type > &left, const c2_function< float_type > &right, float_type x, float_type *yprime, float_type *yprime2)
execute math necessary to do composition
Definition: c2_function.hh:986
c2_typed_ptr()
construct the container with no function
Definition: c2_function.hh:728
c2_ptr(const c2_ptr< float_type > &src)
copy constructor
Definition: c2_function.hh:692
create a c2_function which is the difference of two other c2_functions.This should always be construc...
Definition: c2_function.hh:81
create a c2_function which is a piecewise assembly of other c2_functions.The functions must have incr...
Definition: c2_function.hh:84
create a power law mapping of another functionfor example, given a c2_function f
virtual float_type value_with_derivatives(float_type x, float_type *yprime, float_type *yprime2) const
get the value and derivatives.
Provides support for c2_function objects which are constructed from two other c2_function objects...
Definition: c2_function.hh:845
virtual void set_sampling_grid_pointer(std::vector< float_type > &grid)
Definition: c2_function.hh:440
void get_data(std::vector< float_type > &xvals, std::vector< float_type > &yvals) const
retrieve copies of the x & y tables from which this was built
a transformation of a function in and out of log-lin space
a c2_function which is constantThe factory function c2_factory::constant() creates *new c2_constant_p...
static bool comp_pair(std::pair< float_type, float_type > const &i, std::pair< float_type, float_type > const &j)
const bool fHasStaticTransforms
flag to indicate if the static function pointers can be used for efficiency
void unset_function(void)
clear the function
Definition: c2_function.hh:609
compute log(x) with its derivatives.The factory function c2_factory::log() creates *new c2_log_p ...
void set_function(const c2_function< float_type > *f)
fill the container with a new function, or clear it with a null pointer
Definition: c2_function.hh:580