Spline Class Reference

#include <spline.h>

Inheritance diagram for Spline:

Mathobject AnchorEnabled List of all members.

Public Types

enum  childrenToo { no, all, thoseReady, yes }

Public Member Functions

void add (int, double)
 add x to ydat[k]
void addChild (Spline *c)
 adds another Child to the Child-map
void arm (childrenToo=no, bool keepSize=false, int=0)
double average ()
double average (double b)
double average (double a, double b, double to=1e-6)
double back (int k) const
double back () const
void checkConvolutionRange (Spline *p, double *, double *, const double, const double)
void checkSpace (int=-2)
void checksum ()
double ck (double x)
double conv (const double)
 romberg-convolution integrand
void convD (const double, const double *, double *)
 usual convolution integrand
double convolution (Spline *, double, double, double=1e-4, double=1.0, double=1.0)
double convolutionResult ()
double convOneDim (const double)
void convVektor2 (const double, const double *, double *)
 new version
void createChecksum ()
void derive (Spline &s, childrenToo=no)
void disarm (childrenToo=no, bool keepSize=false)
 deletes the spline-interpolation table
void div (int k, double x)
 call mul(1/x)
void dump (ofstream &, bool=true)
void dump ()
void dump (string, childrenToo InterpolatedAlso=no)
void dump (string, bool)
void dumpConvolution (Spline *, double, double, double, double, const char *)
void error (string, string, int=1234) const
 shows this spline, mother and a lot of information
void explizit ()
double fastSplint (double xa[], double x, int *guess)
 Spline interpolation Nrecipes + guess for next Spline interpolation, modified version that uses sequential search starting from a guess, instead of interval halfing.
double fastY (const double x)
 stub to operator()
double findZero (double a, double b, double=1e-5)
double findZeroBetween (const int a, const int b, double=1e-5, const double zl=0.0)
int first () const
 index of first data-point
double fitNeeds (const pair< double, double >)
 Fits cubic splines to y and interpolates first derivs at grid points dy.
void fitToSingle (const pair< double, double >)
void flip ()
 flip the x axis, i.e. make first last and so on
double front (int k) const
double front () const
double generatedAt () const
void generateSplineXXL (Spline *e, vector< Spline * > &v, vector< Spline * > &p, const double k)
SplinegetChild (string)
 return child with name s. If there are more than a few, this is costy !!!
void getData (const Spline &, childrenToo=no)
 Copy x (if necessary) and y data from another spline, set n. If spline does not have own x-data, the copy uses mothers-xdata.
void getExtrema (list< double > *, list< double > *, int a=0, int b=0, const double tol=1e-5)
void getMaxima (list< double > *, int a=0, int b=0, const double tol=1e-5)
void getMinima (list< double > *, int a=0, int b=0, const double tol=1e-5)
int getZeroArray (int a, const int b, double *, const int, const double tol=1e-5, double zl=0.0)
list< double > * getZeroList (const double zl=0.0, const double tol=1e-2, int a=0, int b=0)
vector< double > getZeroVector (const double zl=0.0, const double tol=1e-2, int a=0, int b=0)
 convenience wrapper for getZeroList()
double inBetweenY (double)
int indexToLeft (const double x) const
 return the largest i with x(i) <= x
int indexToRight (const double x) const
 index of last data-point return the smallest i with x(i) >= x
double inte (const double)
 romberg-integration integrand
void inteD (const double, const double *, double *)
 odeint integration integrand
double integrate (double tol=1e-6)
double integrate (double a, double b, double tol=1e-6)
double inteOneDim (const double)
bool isArmed () const
bool isWithinBounds (const double x)
 true, if x value is within the boundaries of spline data
int last () const
double maximum ()
 return x of the largest maximum in this spline *including* the boundaries, if they are higher
void merge (Spline &s, const bool=true)
void merge (const double *, const int, const bool=true, const double *=0, const double keepSize=true)
double minimum ()
 return x of the smallest minimum in this spline *including* the boundaries, if they are below
void motherKilled ()
 if mother is killed
void mul (int k, double x)
 multiply ydat[k] by x
const Splineoperator *= (const double x)
 Scalar multiply all ydat[] by const double x.
double operator() (const double)
 Fast and caching access to the interpolation value. Should usesually be the one to call (not splint() or fastSplint()).
const Splineoperator+= (const double x)
 add x to all ydat[]
double operator[] (const int n)
void origSpline (double x[], double y[], int, int n, double yp1, double ypn, double y2[])
void printStatus () const
void proper ()
double * ptrY ()
double * ptrY (int)
 returns pointer to element i of ydat. Also adjusts n.
double randomY (double)
 return interpolation at point x. use this instead of operator (), if you know that the access is random
double range (int step=1)
 return (last-first)/step (useful for for-loops)
void rearm (bool keepSize=true)
void removeChild (Spline *c)
 removes a child from the childmap
void resize (int)
 allocates arrays for x and y and copies existing data to these
double safeY (const double x)
 used by derive to keep within spline boundaries
void save (const string &)
 save spline x and y data to a file
void save (ofstream &)
 save spline x and y data to file
void set (const map< double, double > &m)
void set (const map< float, float > &m)
void set (const double)
 new data point (y) for a spline that shares the xdata with mother and sisters
void set (const double, const double)
 new data point (x,y) for a spline that owns its xdata
void setChildrensN (int k=-2)
 Call setN() for all children. All precautions of setN are true here, too.
void setForce (double x, double y)
void setKillChildrenWhenDying (childrenToo c)
void setMother (Spline &)
void setN (int k=-2)
 Set n of this spline.
void setName (string s)
 set's the splines name
void setX (double x)
 sets next X
void setX (int, double)
 set xdat[index] to x
void setY (int, double)
 sets ydat[index] to y
int size () const
 the size i.e. n+1
void smoothen (int k=1, int m=2)
 smoothen the spline within k neigbouring points to the left and k to the right and 2*m*k + 1 points in this range
void Splder (double *, double *, const int)
 initialize whatever...
void splder (Spline &s)
 Spline (moSingle, const Mathobject &, Spline *, string name, Anchor *=0)
 use function moSingle to create a spline at x-data of another spline
 Spline (ifstream &i, string="unnamed", Anchor *=0)
 read data from file (i.e re-create a spline from the save() data )
 Spline (const Spline &)
 Create a Spline that is an exact copy of the Spline in the argument.
 Spline (const Spline &, string, Anchor *=0)
 Construct a Spline that has its own data, is un-armed, has no children and otherwise gets its x and y data from the Spline in the argument.
 Spline (Spline *, string="unnamed", Anchor *=0)
 Spline (const int=1000, string="unnamed", Anchor *=0)
void splini ()
double splint (double xa[], double x, int *guess)
double start (int k) const
 Return last x.
double start () const
double stepInt (int a, int b)
double stepIntegrate ()
double stepIntegrate (double a, double b)
double stop (int k) const
double stop () const
 Return first x.
double x (int k) const
 return xdat[k]
double y (int) const
 return ydat[k]
 ~Spline ()

Static Public Member Functions

static void armVector (vector< Spline * > &)
static void disarmVector (vector< Spline * > &)
static vector< Spline * > * expandVector (vector< Spline * > &, vector< Spline * > &, int total)
 Equally -spaced expansion of the spline vector.
static vector< Spline * > * expandVectorBetween (vector< Spline * > &, vector< Spline * > &, int between)
 Generate "between" splines between each existing pair of spines in v.
static void generateAndDump (vector< Spline * > &v, vector< Spline * > &p, string, const double k)
static void generateChildren (vector< Spline * > &, Spline *, const unsigned int, const string)
static void generateFamily (vector< Spline * > &, const unsigned int, const int, const string)
static void generateOne2OneChild (vector< Spline * > &, vector< Spline * > &, const string)
static void generateSisters (vector< Spline * > &, const unsigned int, const int, const string)
static void generateSpline (Spline *e, vector< Spline * > &v, vector< Spline * > &p, const double k, bool=false)
static void kill (Spline *)
static void printWatchMap ()
 output how many splines with what name are currently alive
static void vectorConvolution2 (Spline *, double, double, double=1e-4, double=1.0, double=1.0, double hnext=0)
 convolution function for family convolution type: moDerivs

Public Attributes

double a
double a3a
int acc
 Counter for operator() this is for perfomance monitoring.
bool armed
double b
 caching these for children's use
double b3b
int cC
double checkSm
double checkSm2
bool childExists
 if the spline data is calculated if children is not empty
map< Spline *, Spline * > children
 a map of children, if there are any
double convResult
 the result of a family convolution, if requested
double gen
int generated
int guess
 set by splint and fastsplint to have faster access in subconsequent calls
double h
double h26
int khi
 caching these for children's use
childrenToo killChildrenWhenDying
int klo
int lastarm
 keep track of last arm-position for rearm(), if you want to use this feature...
double lastX
int maxDim
Splinemother
 if there is a mother, its this
int n
string name
 the name (useful for debugging)
int nvc
bool own
 if the xdat belongs to the spline, i.e. mother ==0. Own does not say, that xdat is a valid pointer, it may be 0, if for instance, mother is killed, own will be true, but the pointer will be set to zero
double * sdat
bool splgValid
 derive() needs the splg[] array, which is initialized if it is not valid
int splintcC
bool valid
 if the xdat pointer is valid, i.e. own or mother alive And valid
bool validCache
int vc
double * xdat
double * ydat

Static Public Attributes

static vector< Spline * > convolutionVektor
static SplineconvTwin
 for the specialised twin convolution, the second twin
static double * convXdat
 for vectorConvolution this is the pointer to the actual x-data. Either taken fom the first spline in convolutionVektor or its mother
static double faktor
 and a faktor for x in convolution
static unsigned int nvar
 The size of the convolution vector.
static Splinepartner
 convolutions partner spline
static double shift
 a shift in integration
static double splg [20010]
 splini and splder need it...
static int WatchCount = 0
 If spline is compiled with ENABLE_SPLINE_WATCH, counts number of splines alive.
static map< string, int > WatchMap
 If compiled with ENABLE_SPLINE_WATCH, tracks number of splines with that name.
static double zeroLevel

Detailed Description

Class for cubic spline interpolation of data. Supports shared x-data among splines that have this data in common, fast access to the right interpolation interval and caching of precalculated data.

It is fast, efficient and stable.

Constructors
There are several constructors:
  Spline(const int = 1000, string = "unnamed",Anchor* =0);
  Spline(Spline *, string = "unnamed",Anchor* = 0);
  Spline(const Spline&);
  Spline(ifstream &i, string = "unnamed",Anchor* =0); 
  Spline(moSingle, const Mathobject&, Spline*, string name, Anchor* =0);
Common to all is that you can specify a name (which you should for debugging is easier with this) and an Anchor. Splines are AnchorEnabled and hence the destruction of the Anchor will destroy the Spline, if you give an Anchor. This makes it very easy to prevent memory leaks. The first constructor creates a spline with default 1000 data points. If you fill more than 1000 points, this is no problem, as Splines will double their size each time you set more data points than allocated. On calling arm(), the overhead memory is freed. The second constructor is for constructing a child Spline from a mother given as the first argument. Children share the same x-data with their mother spline, in fact, the x-data is only allocated once for the mother and never again. Using child splines provides very effecient caching possibilities: Say you have two splines, mother m and child c and you usually access both in the same routine at the same x-value. If you then use m(x) and subsequently c(x), the c(x) call will notice that already several variables have been calculated for the mother and these will be used for calculating c(x). The fourth constructor reads data from file. Use save() to write this file beforehand. The last constructor is for constructing child splines that given a function moSingle func which you specify returns a spline with the same x-data as mother (naturally) and the y-data set to func(x-data)
Basic usage
Using Splines is deceivingly simple. After you have created a spline, you set the data points using the set(x,y) or if the spline is child spline the set(y) function. If you just don't bother, use setForce(x,y), this will decide for you. Setting is always in some order, either decending or ascending. If you set in descending order, you will have to call flip() to flip the data in ascending order, before you arm() the spline, cause spline x-data comes in increasing values. It is not allowed to arm() twice (see below, why you should want to do this). Call disarm() first. Then you may call arm() again to get the proper Spline-interpolation table. So as an example, we create a sine - spline
  Spline example(1000,"example");
  for (double x = -3; x < 5 ; x += 0.1) 
  example.set(x, sin(x));
In order to access the interpolation capabilities, you now have to call arm(). This builds up the table of second derivatives. If you try to access and the spline is not arm()ed, you will be thrown an Bad_Error().
  example.arm();
This is all. If you now would like to know the value at x=0.12345, you call
  double y = example(0.12345);
As a synonym, you can call fastY() this is exactly the same and only provided for cases where you have a pointer to a spline. In this case you have two possibilities: either use the operator () which sometimes looks akward or use fastY. We also use an anchor here:
  Anchor anchor;
  Spline *pointerSpline = new Spline(1000,"pspline",&anchor);
  for (double x = -3; x < 5 ; x += 0.1) 
  pointerSpline->set(x, sin(x));
  pointerSpline->arm();
  double y  = pointerSpline->fastY(0.12345); // 1st possib.
  double y2= (*pointerSpline)(0.12345);        // 2nd possib
Arithmetics
You can multiply or add a usual number to the whole spline y-data, say we already have a spline example with values set, then we could
  example *= 0.3;    // multiply all y-data by 0.3
  example += 5;      // and add a constant 5 to it
If you would like to act differently on different data points, you can call
  example.mul(7,0.3);
  example.mul(8,0.5);
This will multiply data point #7 by 0.3 and #8 by 0.5. If for instance you would like to multiply the spline by a function sin(x*x), you could do:
  for (int i = 0; i < example.size(); i++) 
  example.mul(i, sin(example.x(i)*example.x(i));
The size() just returns the number of data points, x(i) returns x-data points number i and hence you multiply by sin(x*x). After doing arithmetics, you will have to call arm() again. Keep in mind that if you have already called arm(), you will have to disarm() first.
Differentiation and Integration
You may differentiate a given spline, say the spline holding sine data
  Spline Sine, Cosine;  // we dont care for names here
  .... fill Sine with sine values
  Sine.arm();
  Sine.derive(&Cosine);
  Cosine.dump("cosine");  // write y-data to file "cosine.dat"
  Sine.dump("sine");        
So both splines data be written (human readably) to the files sine.dat and cosine.dat and you may take a look at them in e.g. gnuplot.
If you like to know the integral of the Sine spline above ,call
  double y  = Sine.integrate();          // whole x -range
  double y2  = Sine.integrate(-1,1);  // x = -1... 1 
  double y3 = Sine.average();     // integrate and divide by x range
Notice the average() call above, returning the integral divided by the integration interval (in this case, no interval is given, which defaults to the entire Spline).

Definition at line 185 of file spline.h.


Member Enumeration Documentation

enum Spline::childrenToo
 

Enumerator:
no 
all 
thoseReady 
yes 

Definition at line 187 of file spline.h.


Constructor & Destructor Documentation

Spline::Spline const   int = 1000,
string  nom = "unnamed",
Anchor a = 0
 

Constructor for an independent spline (owning its x-data). You may specify the initial size of the x and y data arrays as well as a name, which is useful if something goes wrong, cause error can tell you which spline actually crashed.

The initial size is not too important, as the spline will dynamically allocate more memory if needed, also, when arm()ed, the spline will free access space.

Definition at line 19 of file mathobject/spline.cc.

References name, WatchCount, WatchMap, xdat, and ydat.

Referenced by expandVectorBetween(), generateAndDump(), generateChildren(), generateFamily(), generateOne2OneChild(), generateSisters(), and getExtrema().

Spline::Spline Spline ,
string  = "unnamed",
Anchor = 0
[explicit]
 

Definition at line 35 of file mathobject/spline.cc.

References addChild(), error(), maxDim, mother, name, own, WatchCount, WatchMap, and ydat.

Spline::Spline const Spline s,
string  nom,
Anchor a = 0
[explicit]
 

Construct a Spline that has its own data, is un-armed, has no children and otherwise gets its x and y data from the Spline in the argument.

This constructor is in some ways superior to the Spline(const Spline&), because it only is a copy w.r.t to the data.

Definition at line 68 of file mathobject/spline.cc.

References getData(), maxDim, name, WatchCount, WatchMap, xdat, and ydat.

Spline::Spline const Spline  )  [explicit]
 

Create a Spline that is an exact copy of the Spline in the argument.

Definition at line 52 of file mathobject/spline.cc.

References all, armed, error(), getData(), maxDim, name, own, sdat, WatchCount, WatchMap, xdat, and ydat.

Spline::Spline ifstream &  i,
string  = "unnamed",
Anchor = 0
 

read data from file (i.e re-create a spline from the save() data )

Definition at line 80 of file mathobject/spline.cc.

References maxDim, set(), x(), xdat, y(), and ydat.

Spline::Spline moSingle  func,
const Mathobject mobj,
Spline s,
string  nom,
Anchor a = 0
 

use function moSingle to create a spline at x-data of another spline

Create a spline as child of mother s using (naturally) mother x-data and y-data coming from calls of moSingle func(x-data) This is an ultra - convenient way of generating mappings of splines using arbitrary functions :-)

Definition at line 107 of file mathobject/spline.cc.

References addChild(), error(), maxDim, mother, n, name, own, WatchCount, WatchMap, x(), and ydat.

Spline::~Spline  ) 
 

Definition at line 127 of file mathobject/spline.cc.

References all, children, and killChildrenWhenDying.


Member Function Documentation

void Spline::add int  ,
double 
 

add x to ydat[k]

Definition at line 1625 of file mathobject/spline.cc.

References maxDim, n, printStatus(), and ydat.

Referenced by AllSkyLensing::oneAngle().

void Spline::addChild Spline c  )  [inline]
 

adds another Child to the Child-map

Definition at line 245 of file spline.h.

References wmap_tt_beam_ptsrc_chisq::c, childExists, and children.

Referenced by setMother(), and Spline().

void Spline::arm childrenToo  = no,
bool  keepSize = false,
int  = 0
 

Definition at line 194 of file mathobject/spline.cc.

References all, armed, childExists, children, error(), guess, maxDim, mother, n, origSpline(), own, proper(), resize(), sdat, thoseReady, valid, xdat, and ydat.

Referenced by AnalyzeThis::ACBARChi2WithCalibration(), armVector(), AnalyzeThis::BOOMERANGInit(), AnalyzeThis::CBIChi2WithCalibration(), AnalyzeThis::CBIMosaicChi2WithCalibration(), Cosmos::createLensingPowerLimber(), Cosmos::createPower(), CmbMainWindow::drawLikeli_1d(), CmbMainWindow::endSpool(), SplineWeb::fitToSingleAlongX(), SplineWeb::fitToSingleAlongY(), FlatSkyLensing::FlatSkyLensing(), generateAndDump(), getExtrema(), Cosmos::growthFactor(), distCosmos::history(), Cosmos::history(), AnalyzeThis::isSane(), AnalyzeThis::lymanAlphaPatMcDonaldChi2(), SplineWeb::maximum_internal(), Model::outputScalarCl(), CrossoverField::prepare(), Arthur::prepare(), rearm(), PlotWidget::setSpline(), CmbMainWindow::sigma8(), AnalyzeThis::sigma8(), AnalyzeThis::Sn1aCore(), AnalyzeThis::TwoDF_bestBias(), AnalyzeThis::VSAChi2WithCalibration(), AnalyzeThis::WMAP5computeLikelihood(), AnalyzeThis::WMAPLikelihood(), and AnalyzeThis::WMAPNormalize().

void Spline::armVector vector< Spline * > &   )  [static]
 

Definition at line 1669 of file mathobject/spline.cc.

References arm().

Referenced by SplineWeb::arm(), SplineWeb::createAlongX(), SplineWeb::createAlongY(), SplineWeb::dumpAlongX(), SplineWeb::dumpAlongY(), SplineWeb::fitToSingleAlongX(), and SplineWeb::fitToSingleAlongY().

double Spline::average  )  [inline]
 

Definition at line 419 of file spline.h.

References stop().

Referenced by average().

double Spline::average double  b  )  [inline]
 

Definition at line 418 of file spline.h.

References average(), and start().

double Spline::average double  a,
double  b,
double  to = 1e-6
 

Definition at line 1064 of file mathobject/spline.cc.

References integrate().

Referenced by QuintCosmos::printStatus(), and QuintCosmos::weff().

double Spline::back int  k  )  const [inline]
 

Definition at line 369 of file spline.h.

References n, and y().

double Spline::back  )  const [inline]
 

Definition at line 367 of file spline.h.

References n, and y().

Referenced by safeY().

void Spline::checkConvolutionRange Spline p,
double *  ,
double *  ,
const   double,
const   double
 

Definition at line 984 of file mathobject/spline.cc.

References cC, faktor, max(), min(), n, partner, shift, start(), stop(), WARN, and x().

Referenced by convolution(), dumpConvolution(), and vectorConvolution2().

void Spline::checkSpace int  = -2  ) 
 

Definition at line 513 of file mathobject/spline.cc.

References error(), maxDim, n, resize(), and valid.

Referenced by set(), setX(), and setY().

void Spline::checksum  ) 
 

After you have creatChecksum()ed, you can call checksum(). It will reevaluate and if the checksums do not coincide, throw a Bad_Error

Definition at line 367 of file mathobject/spline.cc.

References ck(), n, own, and xdat.

double Spline::ck double  x  )  [inline]
 

Definition at line 252 of file spline.h.

Referenced by checksum(), and createChecksum().

double Spline::conv const   double  ) 
 

romberg-convolution integrand

void Spline::convD const   double,
const double *  ,
double * 
 

usual convolution integrand

Definition at line 709 of file mathobject/spline.cc.

References cC, faktor, and shift.

Referenced by convolution().

double Spline::convolution Spline p,
double  start,
double  end,
double  tol = 1e-4,
double  fak = 1.0,
double  shft = 1.0
 

Convolute this Spline with a partner Spline p. The Convolution starts at start and ends at end, and has the tolerance tol.

In principle this looks like: ^end this(x) times partner(fak * x + shft) dx

So fak(tor) and sh(i)ft give additional freedom to convolute shifted and streched partners.

Definition at line 699 of file mathobject/spline.cc.

References checkConvolutionRange(), convD(), Miscmath::rungeInt(), and y().

double Spline::convolutionResult  )  [inline]
 

Definition at line 463 of file spline.h.

References convResult.

double Spline::convOneDim const   double  ) 
 

Definition at line 972 of file mathobject/spline.cc.

References cC, error(), faktor, n, partner, shift, and x().

void Spline::convVektor2 const   double,
const double *  ,
double * 
 

new version

Definition at line 822 of file mathobject/spline.cc.

References a, a3a, b, b3b, convolutionVektor, faktor, fastSplint(), guess, h26, khi, klo, nvar, own, partner, shift, validCache, and xdat.

void Spline::createChecksum  ) 
 

Kept for debugging purposes. If You think something very strange is going on, you can ask the spline to build checksums of its y and s data. Hence you will be able to notice memory corruptions

Definition at line 348 of file mathobject/spline.cc.

References checkSm, ck(), n, own, sdat, xdat, and ydat.

void Spline::derive Spline s,
childrenToo  ct = no
 

Fill spline s with data of derivatives of this spline. In contrast to splder(), no equally spaced points are assumed and the spline function is evaluated, nothing else....

The x-data does not have to coincide (spline s should however lie in the definition range of this spline)

If spline s has no xdata, the xdata of this spline is taken.

Definition at line 643 of file mathobject/spline.cc.

References Miscmath::dfridr(), maxDim, n, own, safeY(), valid, x(), xdat, and ydat.

Referenced by getExtrema(), distCosmos::history(), Arthur::prepare(), and AnalyzeThis::Sn1aCore().

void Spline::disarm childrenToo  = no,
bool  keepSize = false
 

deletes the spline-interpolation table

Definition at line 237 of file mathobject/spline.cc.

References all, armed, childExists, children, own, sdat, and validCache.

Referenced by AnalyzeThis::ACBARChi2WithCalibration(), AnalyzeThis::CBIChi2WithCalibration(), AnalyzeThis::CBIMosaicChi2WithCalibration(), FlatSkyLensing::cllens_(), disarmVector(), smoothen(), AnalyzeThis::TwoDF_bestBias(), AnalyzeThis::VSAChi2WithCalibration(), and AnalyzeThis::WMAPNormalize().

void Spline::disarmVector vector< Spline * > &   )  [static]
 

Definition at line 1670 of file mathobject/spline.cc.

References disarm().

Referenced by SplineWeb::disarm(), SplineWeb::fitToSingleAlongX(), and SplineWeb::fitToSingleAlongY().

void Spline::div int  k,
double  x
[inline]
 

call mul(1/x)

Definition at line 339 of file spline.h.

References mul().

Referenced by FlatSkyLensing::cllens_().

void Spline::dump ofstream &  o,
bool  only = true
 

Write spline data to ofstream o, If bool only is false, 10 interpolated points per spline - data point are written (to give a smoother picture

Definition at line 312 of file mathobject/spline.cc.

References n, x(), and ydat.

void Spline::dump  ) 
 

Definition at line 298 of file mathobject/spline.cc.

References error(), and name.

Referenced by dump(), and proper().

void Spline::dump string  nam,
childrenToo  InterpolatedAlso = no
 

Wrapper for dump(string,bool). Reason: I can never remember, if bool has to be true or false to get in addition to name.dat a file name.plt with interpolated spline values in addition.

So I now use the plain-text attribute "yes" and "no"

Definition at line 273 of file mathobject/spline.cc.

References dump(), and no.

void Spline::dump string  nam,
bool  only
 

Write spline data to file called "nam.dat". If bool only is false, it will also create a file "nam.plt" containing 10 times as many points interpolated from this spline

Definition at line 283 of file mathobject/spline.cc.

References dump().

Referenced by AnalyzeThis::ACBARChi2WithCalibration(), AnalyzeThis::CBIMosaicChi2WithCalibration(), Cosmos::dumpPower(), SplineWeb::fitToSingleAlongX(), generateAndDump(), distCosmos::history(), and Cosmos::history().

void Spline::dumpConvolution Spline p,
double  start,
double  end,
double  fak,
double  shft,
const char *  datei
 

Print the convolution integrand of this and partner to file datei

Definition at line 907 of file mathobject/spline.cc.

References checkConvolutionRange(), par, and x().

void Spline::error string  ,
string  ,
int  = 1234
const
 

shows this spline, mother and a lot of information

Definition at line 1559 of file mathobject/spline.cc.

References mother, own, and printStatus().

Referenced by arm(), checkSpace(), convOneDim(), dump(), fastSplint(), fitNeeds(), getChild(), getData(), getZeroArray(), indexToLeft(), indexToRight(), merge(), motherKilled(), operator()(), proper(), ptrY(), randomY(), resize(), set(), setN(), setX(), Spline(), splint(), vectorConvolution2(), x(), and y().

vector< Spline * > * Spline::expandVector vector< Spline * > &  v,
vector< Spline * > &  p,
int  total
[static]
 

Equally -spaced expansion of the spline vector.

Like expandVectorBetween, however it does not generate a fixed number of new splines in between each pair of existing splines, but a equally spaced total number total between to left and rightmost existing spline

Definition at line 1783 of file mathobject/spline.cc.

References first(), generateFamily(), and generateSpline().

vector< Spline * > * Spline::expandVectorBetween vector< Spline * > &  v,
vector< Spline * > &  p,
int  between
[static]
 

Generate "between" splines between each existing pair of spines in v.

Now, it gets tricky.

In some applications (as in cmbfast) it is useful to calculate some function S(k,tau) on a grid of k and tau values.

In the later calculation, we may need this function on a much denser grid.

Now, assume, as we don't have 2-dim Splines, but only 1-dim, that there are 2 vectors<Spline*> each holding splines T and K respectively that do a T_k(tau) and K_tau(k) interpolation of S(k,tau).

The subscript _tau and _k resembles the fact that we deal with vectors T[k](tau).

On calling expandVector( vector of T splines, vector of K splines, int between), a new vector is calculated that in between of any k value, the T splines exist (i.e. T_k), has "between" new Splines that are also T splines but with k values in between two neighbouring T splines.

Now how do we know which k values to take and who will give us the right Datapoints ???

To answer this, we assume that the Splines K_tau(k) have tau values EXACTLY at the tau data-points of T_k(tau), i.e. for all values that exist in T[k]->x(i), there is a spline K[tau].

Put in other words: K[1] is a spline along the k-direction that has fixed value tau = T[anything]->x(1) K[2] = T[anything]->x(2) etc. etc.

If this sounds pretty strange then please observe that in most situation exactly these properties arise:

Example:

for (k=... ) for (tau = ...)

S( k, tau ) = some function

T[k]->set(tau, S(k,tau)) K[tau]->set(k, S(k,tau))

That's it: You now have exactly these splines and on calling e.g expandVector(T, K, 3) you get about two to three times the number of T[k] splines back that let you interpolate much better for differen k values

Technical note:

(1) expandVector() does new a vector<Spline*>, please delelte it after use (2) as you will want to know the k value these new splines belong to, call generatedAt() (3) vector[1] is mother, all others are children. This is true for the result as well as for the argument (4) all xdata is assumed to INCREASE with index i (5) the new k - data is EQUALLY SPACED, i.e. now logarithmic spacing in between each interval. if however the intervalls have already logarithmic spacing and are not to far apart, then this should not matter too much

Definition at line 1731 of file mathobject/spline.cc.

References wmap_tt_beam_ptsrc_chisq::c, first(), name, and Spline().

void Spline::explizit  ) 
 

Definition at line 1376 of file mathobject/spline.cc.

References mother, n, own, x(), xdat, and ydat.

double Spline::fastSplint double  xa[],
double  x,
int *  guess
[inline]
 

Spline interpolation Nrecipes + guess for next Spline interpolation, modified version that uses sequential search starting from a guess, instead of interval halfing.

The same as the original splint, except that the original version does a half step strategy in order to find the right interval and this one exepts a guess, and searches sequentially from that guess on for a matching interval. This is of course much faster for a nearby solution to the problem.

However, if after a few sequential searches, it didn't succeed, it calles splint() for its half-step strategy

In additon, it caches some quantities for faster access by related splines.

The guess is set to the new guess.

if too many steps, do half step

if too many steps, do half step

Definition at line 512 of file spline.h.

References a, a3a, b, b3b, childExists, error(), h, h26, khi, klo, lastX, n, own, sdat, splint(), validCache, and ydat.

Referenced by convVektor2(), generateSplineXXL(), and operator()().

double Spline::fastY const double  x  )  [inline]
 

stub to operator()

Definition at line 377 of file spline.h.

Referenced by distCosmos::angulardiameterDistance(), Cosmos::angulardiameterDistance(), distCosmos::angulardiameterDistanceDeriv(), Cosmos::angulardiameterDistanceDeriv(), CrossoverField::betaFunction(), AllSkyLensing::C_gl2(), Cosmos::createLensingPowerLimber(), Cosmos::createPower(), Recombination::cs2(), Crossover::d2wda2(), Crossover::dwda(), findZero(), fitNeeds(), CrossoverField::initialQ(), Arbitrary::initialQDot(), CrossoverField::kineticTerm(), distCosmos::luminosityDistance(), Cosmos::luminosityDistance(), distCosmos::luminosityDistanceDeriv(), Cosmos::luminosityDistanceDeriv(), maximum(), minimum(), Cosmos::Pk2Delta2(), distCosmos::propermotionDistance(), Cosmos::propermotionDistance(), distCosmos::propermotionDistanceDeriv(), Cosmos::propermotionDistanceDeriv(), Arbitrary::q(), Arbitrary::qDot(), Arbitrary::rho(), safeY(), CmbMainWindow::sig8Integrator(), AnalyzeThis::sig8Integrator(), AllSkyLensing::sigma2(), Perturbation::soundSpeed(), FlatSkyLensing::t2_(), Recombination::Tb(), Recombination::TbDeriv(), Arbitrary::Vprime(), Crossover::w(), AnalyzeThis::WMAP5computeLikelihood(), Recombination::xe(), distCosmos::Z2iH(), and Cosmos::Z2iH().

double Spline::findZero double  a,
double  b,
double  = 1e-5
 

Definition at line 1083 of file mathobject/spline.cc.

References fastY(), SPLINE_DBG_MSG, and Miscmath::zbrent().

double Spline::findZeroBetween const int  a,
const int  b,
double  = 1e-5,
const double  zl = 0.0
 

Definition at line 1088 of file mathobject/spline.cc.

References inBetweenY(), khi, klo, SPLINE_DBG_MSG, x(), Miscmath::zbrent(), and zeroLevel.

Referenced by getZeroArray(), and getZeroList().

int Spline::first  )  const [inline]
 

index of first data-point

Definition at line 302 of file spline.h.

Referenced by expandVector(), expandVectorBetween(), NewdatClChi2::init(), CBI2::init(), AllSkyLensing::LensingDifference::LensingDifference(), and AllSkyLensing::oneAngle().

double Spline::fitNeeds const pair< double, double >   ) 
 

Fits cubic splines to y and interpolates first derivs at grid points dy.

Definition at line 1907 of file mathobject/spline.cc.

References error(), fastY(), and y().

Referenced by fitToSingle(), SplineWeb::fitToSingleAlongX(), and SplineWeb::fitToSingleAlongY().

void Spline::fitToSingle const pair< double, double >   ) 
 

Definition at line 1913 of file mathobject/spline.cc.

References fitNeeds().

Referenced by SplineWeb::fitToSingleAlongX().

void Spline::flip  ) 
 

flip the x axis, i.e. make first last and so on

Definition at line 1066 of file mathobject/spline.cc.

References n, own, x(), and xdat.

Referenced by Cosmos::history().

double Spline::front int  k  )  const [inline]
 

Definition at line 368 of file spline.h.

References y().

double Spline::front  )  const [inline]
 

Definition at line 366 of file spline.h.

References y().

Referenced by Recombination::highZFraction(), maximum(), minimum(), and safeY().

void Spline::generateAndDump vector< Spline * > &  v,
vector< Spline * > &  p,
string  file,
const double  k
[static]
 

Easy to use convenience function. For example to generate Transferfunctions at some specific time and output it to a file

Input: vector<Spline*> &v must be a vector that stores function values of the same k-position as the yet to be generated spline.

vector<Spline*>&p has the same function values but knows the corresponding tau values.

The variable k in the function argument is then the TIME (not k :-) at which you would like to have a slice through this 2-D plot, so to speak. The direction is of course the k-direction of the plot, i.e. the x-values of the v splines.

If you need examples, look at cosmos.cc

Definition at line 1869 of file mathobject/spline.cc.

References arm(), dump(), generateSpline(), and Spline().

Referenced by SplineWeb::dumpAlongX(), and SplineWeb::dumpAlongY().

void Spline::generateChildren vector< Spline * > &  ,
Spline ,
const unsigned  int,
const   string
[static]
 

Definition at line 1659 of file mathobject/spline.cc.

References Spline().

Referenced by AllSkyLensing::lensedCls(), CmbCalc::prepareScalarCl(), and CmbCalc::prepareTensorCl().

double Spline::generatedAt  )  const
 

If this Spline has been automatically generated via expandVector() then you may very well be interested which value of "kk" of expandVector() belongs to this specific spline. In other words and in the language of expandVector():

We generated splines which are interpolating along the tau direction by taking "kk" values along the k direction in between a given spline-vector interpolating along the k direction at the tau values of the splines T.

We store thes kk values so you know to which value of k this spline belongs.

Definition at line 1890 of file mathobject/spline.cc.

References gen.

Referenced by printStatus().

void Spline::generateFamily vector< Spline * > &  ,
const unsigned  int,
const   int,
const   string
[static]
 

Definition at line 1645 of file mathobject/spline.cc.

References Spline().

Referenced by expandVector(), AllSkyLensing::lensedCls(), CmbCalc::prepareScalarCl(), and SplineWeb::SplineWeb().

void Spline::generateOne2OneChild vector< Spline * > &  ,
vector< Spline * > &  ,
const   string
[static]
 

Definition at line 1664 of file mathobject/spline.cc.

References Spline().

void Spline::generateSisters vector< Spline * > &  ,
const unsigned  int,
const   int,
const   string
[static]
 

Definition at line 1654 of file mathobject/spline.cc.

References Spline().

Referenced by FlatSkyLensing::FlatSkyLensing().

void Spline::generateSpline Spline e,
vector< Spline * > &  v,
vector< Spline * > &  p,
const double  k,
bool  linear = false
[static]
 

Give a vector<Spline*> v with splines T interpolating along some tau-axis and vector<Spline*> p of splines K interpolating along some k-axis,

fill an empty spline e with data points at tau values that are at the x-data of T containing values interpolated using the K splines.

It is assumed that the first k-spline coressponds to the first tau x-data value of the T-Splines.

See comments at expandVectorBetween()

Definition at line 1817 of file mathobject/spline.cc.

References gen, and setForce().

Referenced by expandVector(), SplineWeb::fitToSingleAlongX(), SplineWeb::fitToSingleAlongY(), and generateAndDump().

void Spline::generateSplineXXL Spline e,
vector< Spline * > &  v,
vector< Spline * > &  p,
const double  k
 

Same as generateSpline() however, it assumes that p[0] is this spline and that this spline is the mother of all other p[j] Splines.

Hence, it can efficiently cache the access and speed things up

Definition at line 1834 of file mathobject/spline.cc.

References a, a3a, b, b3b, fastSplint(), gen, guess, h26, khi, klo, sdat, setForce(), validCache, x(), xdat, y(), and ydat.

Spline * Spline::getChild string   ) 
 

return child with name s. If there are more than a few, this is costy !!!

Definition at line 182 of file mathobject/spline.cc.

References children, and error().

void Spline::getData const Spline ,
childrenToo  = no
 

Copy x (if necessary) and y data from another spline, set n. If spline does not have own x-data, the copy uses mothers-xdata.

Definition at line 459 of file mathobject/spline.cc.

References error(), maxDim, own, size(), x(), xdat, y(), and ydat.

Referenced by PlotWidget::setSpline(), and Spline().

void Spline::getExtrema list< double > *  ,
list< double > *  ,
int  a = 0,
int  b = 0,
const double  tol = 1e-5
 

Definition at line 1176 of file mathobject/spline.cc.

References arm(), armed, derive(), getZeroList(), and Spline().

Referenced by getMaxima(), and getMinima().

void Spline::getMaxima list< double > *  ,
int  a = 0,
int  b = 0,
const double  tol = 1e-5
 

Definition at line 1202 of file mathobject/spline.cc.

References getExtrema().

Referenced by maximum().

void Spline::getMinima list< double > *  ,
int  a = 0,
int  b = 0,
const double  tol = 1e-5
 

Definition at line 1211 of file mathobject/spline.cc.

References getExtrema().

Referenced by minimum().

int Spline::getZeroArray int  a,
const int  b,
double *  ,
const   int,
const double  tol = 1e-5,
double  zl = 0.0
 

Definition at line 1144 of file mathobject/spline.cc.

References error(), findZeroBetween(), and y().

list< double > * Spline::getZeroList const double  zl = 0.0,
const double  tol = 1e-2,
int  a = 0,
int  b = 0
 

get a list of zeros or actually crossings of zl (zerolevel) of this spline, starting from data point a to data point b. If the y-data at a coincides with zl, then this point is skipped. This happens with all subsequent points. Thus, if for instance you would like to know the zeros of your Spline and say from start() = -3.5 to x-data = -1.7 the Spline is zero at the x-data points of this inverval, then the search for zeros will only start from the next x-data points higher than -1.7 on.

Naturally, in between two data points, there can only be one zero and hence this procedure will never cost you any zeros *except* trailing identical zeros.

Definition at line 1109 of file mathobject/spline.cc.

References findZeroBetween(), n, SPLINE_DBG_MSG, and y().

Referenced by getExtrema(), getZeroVector(), CrossoverField::initialQ(), and Arthur::initialQ().

vector< double > Spline::getZeroVector const double  zl = 0.0,
const double  tol = 1e-2,
int  a = 0,
int  b = 0
 

convenience wrapper for getZeroList()

Definition at line 1135 of file mathobject/spline.cc.

References getZeroList().

double Spline::inBetweenY double  x  ) 
 

moSingle for Miscmath::zbrent() for finding the "zeros" of a spline. Zero is relative, cause the static double zeroLevel will be subtracted from the y(x). So for instance: if y(x) = 3 and zeroLevel = 2 then one is returnde. Usually for finding really zero zero's, zeroLevel =0

Definition at line 547 of file mathobject/spline.cc.

References a, b, h, khi, klo, mother, own, sdat, xdat, ydat, and zeroLevel.

Referenced by findZeroBetween().

int Spline::indexToLeft const double  x  )  const
 

return the largest i with x(i) <= x

Definition at line 1443 of file mathobject/spline.cc.

References error(), n, and xdat.

int Spline::indexToRight const double  x  )  const
 

index of last data-point return the smallest i with x(i) >= x

Definition at line 1459 of file mathobject/spline.cc.

References error(), n, and xdat.

double Spline::inte const   double  ) 
 

romberg-integration integrand

Definition at line 1369 of file mathobject/spline.cc.

void Spline::inteD const   double,
const double *  ,
double * 
 

odeint integration integrand

Definition at line 1030 of file mathobject/spline.cc.

Referenced by integrate().

double Spline::integrate double  tol = 1e-6  )  [inline]
 

Definition at line 415 of file spline.h.

References integrate(), max(), min(), n, and x().

double Spline::integrate double  a,
double  b,
double  tol = 1e-6
 

Integrate from double a to double b by performing a runge - kutta integration, i.e. converting the problem to the solution of a ordinary diff-eqn.

Definition at line 1018 of file mathobject/spline.cc.

References inteD(), inteOneDim(), Miscmath::oneDimOdeint(), and Miscmath::rungeInt().

Referenced by Cosmos::angulardiameterDistance(), average(), Cosmos::history(), integrate(), Cosmos::luminosityDistance(), QuintCosmos::omesf(), Cosmos::propermotionDistance(), and QuintCosmos::tau2AvOmega_q().

double Spline::inteOneDim const   double  ) 
 

Definition at line 1028 of file mathobject/spline.cc.

Referenced by integrate().

bool Spline::isArmed  )  const [inline]
 

Definition at line 295 of file spline.h.

References armed.

Referenced by CmbMainWindow::sigma8(), AnalyzeThis::sigma8(), and vectorConvolution2().

bool Spline::isWithinBounds const double  x  ) 
 

true, if x value is within the boundaries of spline data

Definition at line 683 of file mathobject/spline.cc.

References start(), and stop().

Referenced by Cosmos::angulardiameterDistance(), Cosmos::luminosityDistance(), and Cosmos::propermotionDistance().

void Spline::kill Spline  )  [static]
 

Definition at line 1636 of file mathobject/spline.cc.

References children.

int Spline::last  )  const [inline]
 

Definition at line 303 of file spline.h.

References n.

Referenced by CL::createTotalXXSpline(), NewdatClChi2::init(), CBI2::init(), AllSkyLensing::LensingDifference::LensingDifference(), and AllSkyLensing::oneAngle().

double Spline::maximum  ) 
 

return x of the largest maximum in this spline *including* the boundaries, if they are higher

Definition at line 1220 of file mathobject/spline.cc.

References fastY(), front(), getMaxima(), max(), and start().

Referenced by AnalyzeThis::ACBARChi2WithCalibration(), AnalyzeThis::CBIChi2WithCalibration(), AnalyzeThis::CBIMosaicChi2WithCalibration(), CmbMainWindow::drawLikeli_1d(), and AnalyzeThis::VSAChi2WithCalibration().

void Spline::merge Spline s,
const   bool = true
 

Definition at line 1245 of file mathobject/spline.cc.

References error(), merge(), mother, n, own, valid, xdat, and ydat.

void Spline::merge const double *  X,
const   int,
const   bool = true,
const double *  Y = 0,
const double  keepSize = true
 

if keepSize, then the new xdata und ydata arrays do either have the old size maxDim, or if the combined size would be too large, n + max +1.

if keepSize is false, then the new arrays will always have the smallest possible size, however, this may cause doubling soon after. If however, you do not want to add anything after merging then keepSize = false is the better choice

Definition at line 1262 of file mathobject/spline.cc.

References error(), max(), maxDim, n, own, sdat, and x().

Referenced by merge().

double Spline::minimum  ) 
 

return x of the smallest minimum in this spline *including* the boundaries, if they are below

Definition at line 1232 of file mathobject/spline.cc.

References fastY(), front(), getMinima(), and start().

Referenced by AnalyzeThis::WMAPNormalize().

void Spline::motherKilled  ) 
 

if mother is killed

Definition at line 174 of file mathobject/spline.cc.

References armed, error(), mother, own, valid, and xdat.

void Spline::mul int  k,
double  x
 

multiply ydat[k] by x

Definition at line 1631 of file mathobject/spline.cc.

References n, and ydat.

Referenced by div().

const Spline& Spline::operator *= const double  x  )  [inline]
 

Scalar multiply all ydat[] by const double x.

Definition at line 387 of file spline.h.

References n, and ydat.

double Spline::operator() const   double  )  [inline]
 

Fast and caching access to the interpolation value. Should usesually be the one to call (not splint() or fastSplint()).

Definition at line 1578 of file mathobject/spline.cc.

References a, a3a, acc, armed, b, b3b, error(), fastSplint(), guess, h26, khi, klo, lastX, mother, nvc, own, sdat, valid, validCache, vc, xdat, and ydat.

const Spline& Spline::operator+= const double  x  )  [inline]
 

add x to all ydat[]

Definition at line 393 of file spline.h.

References n, and ydat.

double Spline::operator[] const int  n  )  [inline]
 

Definition at line 375 of file spline.h.

References y().

void Spline::origSpline double  x[],
double  y[],
int  ,
int  n,
double  yp1,
double  ypn,
double  y2[]
 

Definition at line 1526 of file mathobject/spline.cc.

Referenced by arm().

void Spline::printStatus  )  const
 

Definition at line 557 of file mathobject/spline.cc.

References acc, armed, checkSm, checkSm2, children, generatedAt(), guess, maxDim, mother, n, name, nvc, own, sdat, splintcC, valid, validCache, vc, xdat, and ydat.

Referenced by add(), error(), and set().

void Spline::printWatchMap  )  [static]
 

output how many splines with what name are currently alive

Definition at line 1917 of file mathobject/spline.cc.

References WatchMap.

void Spline::proper  ) 
 

Definition at line 248 of file mathobject/spline.cc.

References dump(), error(), n, name, and x().

Referenced by arm().

double* Spline::ptrY  )  [inline]
 

Definition at line 344 of file spline.h.

References n.

double * Spline::ptrY int   )  [inline]
 

returns pointer to element i of ydat. Also adjusts n.

Definition at line 486 of file spline.h.

References error(), maxDim, n, and ydat.

double Spline::randomY double   ) 
 

return interpolation at point x. use this instead of operator (), if you know that the access is random

Definition at line 526 of file mathobject/spline.cc.

References armed, error(), mother, own, splint(), valid, and xdat.

double Spline::range int  step = 1  )  [inline]
 

return (last-first)/step (useful for for-loops)

Definition at line 371 of file spline.h.

References n, and x().

void Spline::rearm bool  keepSize = true  ) 
 

Definition at line 189 of file mathobject/spline.cc.

References arm(), lastarm, n, and no.

void Spline::removeChild Spline c  )  [inline]
 

removes a child from the childmap

Definition at line 246 of file spline.h.

References wmap_tt_beam_ptsrc_chisq::c, childExists, and children.

void Spline::resize int  k  ) 
 

allocates arrays for x and y and copies existing data to these

resize xdat and ydat to fit index k, if k > maxDim. If this spline is a mother, it calls resize() for all children

Definition at line 478 of file mathobject/spline.cc.

References childExists, children, error(), maxDim, n, own, valid, and xdat.

Referenced by arm(), and checkSpace().

double Spline::safeY const double  x  ) 
 

used by derive to keep within spline boundaries

Definition at line 676 of file mathobject/spline.cc.

References back(), fastY(), front(), start(), and stop().

Referenced by derive().

void Spline::save const string &   ) 
 

save spline x and y data to a file

Definition at line 325 of file mathobject/spline.cc.

References save().

void Spline::save ofstream &   ) 
 

save spline x and y data to file

Definition at line 330 of file mathobject/spline.cc.

References n, x(), and y().

Referenced by save().

void Spline::set const map< double, double > &  m  ) 
 

Definition at line 408 of file mathobject/spline.cc.

References setForce().

void Spline::set const map< float, float > &  m  ) 
 

Definition at line 404 of file mathobject/spline.cc.

References setForce().

void Spline::set const   double  ) 
 

new data point (y) for a spline that shares the xdata with mother and sisters

Definition at line 395 of file mathobject/spline.cc.

References checkSpace(), n, own, printStatus(), and ydat.

void Spline::set const   double,
const   double
 

new data point (x,y) for a spline that owns its xdata

Definition at line 388 of file mathobject/spline.cc.

References checkSpace(), error(), n, own, xdat, and ydat.

Referenced by AnalyzeThis::ACBARChi2WithCalibration(), AnalyzeThis::CBIChi2WithCalibration(), AnalyzeThis::CBIMosaicChi2WithCalibration(), Cosmos::createLensingPowerLimber(), Cosmos::createPower(), Recombination::delAphaSpline(), CmbMainWindow::drawLikeli_1d(), FlatSkyLensing::epsilongen_(), Cosmos::fillHistorySplines(), CrossoverField::getAlpha(), Cosmos::growthFactor(), distCosmos::history(), Cosmos::history(), AnalyzeThis::lymanAlphaPatMcDonaldChi2(), main(), Model::Model(), Model::power_cdm(), CrossoverField::prepare(), Crossover::prepare(), Arthur::prepare(), Arbitrary::prepare(), Arbitrary::reconstructPhi(), setForce(), Spline(), Model::spline(), Cosmos::thermo(), AnalyzeThis::VSAChi2WithCalibration(), AnalyzeThis::WMAPLikelihood(), and AnalyzeThis::WMAPNormalize().

void Spline::setChildrensN int  k = -2  ) 
 

Call setN() for all children. All precautions of setN are true here, too.

Definition at line 454 of file mathobject/spline.cc.

References children.

void Spline::setForce double  x,
double  y
[inline]
 

Wrapper for set() that depending on wether the spline owns the xdata or not, calls the two possible set() functions

Definition at line 323 of file spline.h.

References own, and set().

Referenced by generateSpline(), generateSplineXXL(), and set().

void Spline::setKillChildrenWhenDying childrenToo  c  )  [inline]
 

Definition at line 331 of file spline.h.

References killChildrenWhenDying.

void Spline::setMother Spline  ) 
 

Definition at line 604 of file mathobject/spline.cc.

References addChild(), armed, mother, own, valid, and xdat.

void Spline::setN int  k = -2  ) 
 

Set n of this spline.

This routine lets you set the counter of the data points of this spline EXPLICITLY. No need to mention that this is dangerous, not nice and should be avoided. However, there are two closely related exceptions:

(1) If the data points starting from some x value on are all zero, this can be useful, as by default the new y-data is always zero'd when allocated.

(2) If the spline is not really initialized, but the daugther of some other spline and for the sake of simple programming it has to be accessed nevertheless (to avoid a lot of if thens), the access delivering zero, of course, then a call to SetN() WITHOUT an argument will set the counter to the n value of the mother which is assumed to contain the amount of data you wish to have controle of here

setN() will return error, if n is out of maxDim range. This is a safety precaution. use resize() explicitly, if you really need to do such strange things. In general:

Do not use this function easily, try set, setX, setY, ptrY etc, which all have checks, increase n if necessary and so on and so on and so on

Definition at line 447 of file mathobject/spline.cc.

References error(), maxDim, mother, n, and own.

void Spline::setName string  s  )  [inline]
 

set's the splines name

Definition at line 301 of file spline.h.

References name.

void Spline::setX double  x  )  [inline]
 

sets next X

Definition at line 341 of file spline.h.

References n, and setX().

void Spline::setX int  ,
double 
 

set xdat[index] to x

Definition at line 412 of file mathobject/spline.cc.

References checkSpace(), error(), n, own, and xdat.

Referenced by setX().

void Spline::setY int  ,
double 
 

sets ydat[index] to y

Definition at line 418 of file mathobject/spline.cc.

References checkSpace(), n, and ydat.

Referenced by AllSkyLensing::DlSplines::dl(), NewdatClChi2::init(), CBI2::init(), and AllSkyLensing::LensingDifference::LensingDifference().

int Spline::size  )  const [inline]
 

the size i.e. n+1

Definition at line 308 of file spline.h.

References n.

Referenced by FlatSkyLensing::cllens_(), FlatSkyLensing::epsilongen_(), CrossoverField::getAlpha(), getData(), Cosmos::growthFactor(), AnalyzeThis::lymanAlphaPatMcDonaldChi2(), Arbitrary::reconstructPhi(), PlotWidget::saveToFile(), and PlotWidget::setSpline().

void Spline::smoothen int  k = 1,
int  m = 2
 

smoothen the spline within k neigbouring points to the left and k to the right and 2*m*k + 1 points in this range

Smoothen uses a gaussian window along the spline to smoothen out wiggles. It starts k points to the left and k to the right and in total uses m-fold more points to interpolate.

In other words: k sets the range over which we average and m sets the number of points that are considered in this range.

It disarms but does not arm, so it is up to you to arm() again

Definition at line 1338 of file mathobject/spline.cc.

References disarm(), max(), min(), n, and x().

void Spline::Splder double *  ,
double *  ,
const   int
 

initialize whatever...

Definition at line 1478 of file mathobject/spline.cc.

References splg, splgValid, and splini().

Referenced by splder().

void Spline::splder Spline s  ) 
 

Definition at line 615 of file mathobject/spline.cc.

References maxDim, n, own, Splder(), valid, xdat, and ydat.

void Spline::splini  ) 
 

Definition at line 1517 of file mathobject/spline.cc.

References splg, and splgValid.

Referenced by Splder().

double Spline::splint double  xa[],
double  x,
int *  guess
 

This ist the true (well, a bit modified) NR spline interpolation, sorry for the misnomer of CMBFASTS splint which is an integration routine

If "x" is out of the range xa[1] xa[n], then the boundary values will be given. Meaning:: You will not get any errormessages, if you are out of the boundary, just the closest available value

Definition at line 1395 of file mathobject/spline.cc.

References a, a3a, b, b3b, error(), h, h26, khi, klo, lastX, mother, n, own, sdat, splintcC, validCache, and ydat.

Referenced by fastSplint(), randomY(), and stepIntegrate().

double Spline::start int  k  )  const [inline]
 

Return last x.

Definition at line 363 of file spline.h.

References x().

double Spline::start  )  const [inline]
 

Definition at line 361 of file spline.h.

References x().

Referenced by Cosmos::a_min(), average(), AllSkyLensing::C_gl2(), checkConvolutionRange(), Cosmos::createPower(), PlotWidget::drawSpline(), isWithinBounds(), maximum(), Recombination::maxRedshift(), minimum(), AllSkyLensing::oneAngle(), safeY(), AllSkyLensing::sigma2(), CmbMainWindow::sigma8(), AnalyzeThis::sigma8(), and Cosmos::t_min().

double Spline::stepInt int  a,
int  b
 

Definition at line 1050 of file mathobject/spline.cc.

References x(), and ydat.

Referenced by stepIntegrate().

double Spline::stepIntegrate  )  [inline]
 

Definition at line 423 of file spline.h.

References max(), min(), n, and x().

double Spline::stepIntegrate double  a,
double  b
 

Definition at line 1035 of file mathobject/spline.cc.

References mother, own, splint(), stepInt(), and xdat.

double Spline::stop int  k  )  const [inline]
 

Definition at line 364 of file spline.h.

References n, and x().

double Spline::stop  )  const [inline]
 

Return first x.

Definition at line 362 of file spline.h.

References n, and x().

Referenced by Cosmos::a_max(), AnalyzeThis::ACBARChi2(), average(), AllSkyLensing::C_gl2(), AnalyzeThis::CBIChi2(), AnalyzeThis::CBIDeepChi2(), AnalyzeThis::CBIMosaicChi2(), checkConvolutionRange(), Cosmos::createPower(), PlotWidget::drawSpline(), isWithinBounds(), AllSkyLensing::oneAngle(), Model::outputScalarCl(), safeY(), AnalyzeThis::SDSS_Pk(), AllSkyLensing::sigma2(), CmbMainWindow::sigma8(), AnalyzeThis::sigma8(), Cosmos::t_max(), and AnalyzeThis::VSAChi2().

void Spline::vectorConvolution2 Spline p,
double  start,
double  end,
double  tol = 1e-4,
double  fak = 1.0,
double  shft = 1.0,
double  hnext = 0
[static]
 

convolution function for family convolution type: moDerivs

Optimized convolution of several splines with S one partner spline P.

Effecient in the case that several splines have the same x-data (this need not be owned by themselves). The splines S are store in the static vector convolutionVektor[] which you have to resize and fill with the Spline S.

vektorConvolution() will then temporariliy make convolutionVektor[0] the mother of all other splines S. Hence, access to the splines is greatly sped up by caching of interpolation variables etc.

Definition at line 775 of file mathobject/spline.cc.

References checkConvolutionRange(), childExists, convolutionVektor, error(), guess, isArmed(), mother, nvar, own, xdat, and y().

Referenced by TensorIntegrator::integrate(), and ScalarIntegrator::integrate().

double Spline::x int  k  )  const [inline]
 

return xdat[k]

Definition at line 353 of file spline.h.

References error(), mother, n, own, and xdat.

Referenced by checkConvolutionRange(), FlatSkyLensing::cllens_(), convOneDim(), derive(), dump(), dumpConvolution(), FlatSkyLensing::epsilongen_(), explizit(), findZeroBetween(), flip(), generateSplineXXL(), CrossoverField::getAlpha(), getData(), NewdatClChi2::init(), CBI2::init(), integrate(), AnalyzeThis::lymanAlphaPatMcDonaldChi2(), merge(), AllSkyLensing::oneAngle(), proper(), range(), Arbitrary::reconstructPhi(), save(), PlotWidget::saveToFile(), smoothen(), Spline(), start(), stepInt(), stepIntegrate(), and stop().

double Spline::y int   )  const [inline]
 

return ydat[k]

Definition at line 481 of file spline.h.

References error(), n, and ydat.

Referenced by back(), AllSkyLensing::C_gl2(), convolution(), CmbCalc::dumpTransfer(), FlatSkyLensing::epsilongen_(), fitNeeds(), front(), generateSplineXXL(), getData(), getZeroArray(), getZeroList(), NewdatClChi2::init(), CBI2::init(), AnalyzeThis::lymanAlphaPatMcDonaldChi2(), AllSkyLensing::oneAngle(), operator[](), AnalyzeThis::quickWMAPNormalize(), Arbitrary::reconstructPhi(), save(), AllSkyLensing::sigma2(), Spline(), vectorConvolution2(), AllSkyLensing::xiLensed(), AllSkyLensing::xiLensed_minus(), AllSkyLensing::xiLensed_plus(), and AllSkyLensing::xiLensed_X().


Member Data Documentation

double Spline::a
 

Definition at line 218 of file spline.h.

Referenced by convVektor2(), fastSplint(), generateSplineXXL(), inBetweenY(), operator()(), and splint().

double Spline::a3a
 

Definition at line 219 of file spline.h.

Referenced by convVektor2(), fastSplint(), generateSplineXXL(), operator()(), and splint().

int Spline::acc
 

Counter for operator() this is for perfomance monitoring.

Definition at line 201 of file spline.h.

Referenced by operator()(), and printStatus().

bool Spline::armed
 

Definition at line 198 of file spline.h.

Referenced by arm(), disarm(), getExtrema(), isArmed(), motherKilled(), operator()(), printStatus(), randomY(), setMother(), and Spline().

double Spline::b
 

caching these for children's use

Definition at line 218 of file spline.h.

Referenced by convVektor2(), fastSplint(), generateSplineXXL(), inBetweenY(), operator()(), and splint().

double Spline::b3b
 

Definition at line 219 of file spline.h.

Referenced by convVektor2(), fastSplint(), generateSplineXXL(), operator()(), and splint().

int Spline::cC
 

Definition at line 224 of file spline.h.

Referenced by checkConvolutionRange(), convD(), and convOneDim().

double Spline::checkSm
 

Definition at line 242 of file spline.h.

Referenced by createChecksum(), and printStatus().

double Spline::checkSm2
 

Definition at line 242 of file spline.h.

Referenced by printStatus().

bool Spline::childExists
 

if the spline data is calculated if children is not empty

Definition at line 199 of file spline.h.

Referenced by addChild(), arm(), disarm(), fastSplint(), removeChild(), resize(), and vectorConvolution2().

map<Spline*, Spline*> Spline::children
 

a map of children, if there are any

Definition at line 212 of file spline.h.

Referenced by addChild(), arm(), disarm(), getChild(), kill(), printStatus(), removeChild(), resize(), setChildrensN(), and ~Spline().

vector< Spline * > Spline::convolutionVektor [static]
 

Definition at line 277 of file spline.h.

Referenced by convVektor2(), TensorIntegrator::integrate(), ScalarIntegrator::integrate(), and vectorConvolution2().

double Spline::convResult
 

the result of a family convolution, if requested

Definition at line 229 of file spline.h.

Referenced by convolutionResult().

Spline * Spline::convTwin [static]
 

for the specialised twin convolution, the second twin

Definition at line 236 of file spline.h.

double * Spline::convXdat [static]
 

for vectorConvolution this is the pointer to the actual x-data. Either taken fom the first spline in convolutionVektor or its mother

Definition at line 235 of file spline.h.

double Spline::faktor [static]
 

and a faktor for x in convolution

Definition at line 234 of file spline.h.

Referenced by checkConvolutionRange(), convD(), convOneDim(), and convVektor2().

double Spline::gen
 

Definition at line 228 of file spline.h.

Referenced by generatedAt(), generateSpline(), and generateSplineXXL().

int Spline::generated
 

Definition at line 226 of file spline.h.

int Spline::guess
 

set by splint and fastsplint to have faster access in subconsequent calls

Definition at line 193 of file spline.h.

Referenced by arm(), convVektor2(), generateSplineXXL(), operator()(), printStatus(), and vectorConvolution2().

double Spline::h
 

Definition at line 218 of file spline.h.

Referenced by fastSplint(), inBetweenY(), and splint().

double Spline::h26
 

Definition at line 218 of file spline.h.

Referenced by convVektor2(), fastSplint(), generateSplineXXL(), operator()(), and splint().

int Spline::khi
 

caching these for children's use

Definition at line 220 of file spline.h.

Referenced by convVektor2(), fastSplint(), findZeroBetween(), generateSplineXXL(), inBetweenY(), operator()(), and splint().

childrenToo Spline::killChildrenWhenDying
 

Definition at line 270 of file spline.h.

Referenced by setKillChildrenWhenDying(), and ~Spline().

int Spline::klo
 

Definition at line 220 of file spline.h.

Referenced by convVektor2(), fastSplint(), findZeroBetween(), generateSplineXXL(), inBetweenY(), operator()(), and splint().

int Spline::lastarm
 

keep track of last arm-position for rearm(), if you want to use this feature...

Definition at line 206 of file spline.h.

Referenced by rearm().

double Spline::lastX
 

Definition at line 218 of file spline.h.

Referenced by fastSplint(), operator()(), and splint().

int Spline::maxDim
 

Definition at line 192 of file spline.h.

Referenced by add(), arm(), checkSpace(), derive(), getData(), merge(), printStatus(), ptrY(), resize(), setN(), splder(), and Spline().

Spline* Spline::mother
 

if there is a mother, its this

Definition at line 194 of file spline.h.

Referenced by arm(), error(), explizit(), inBetweenY(), merge(), motherKilled(), operator()(), printStatus(), randomY(), setMother(), setN(), Spline(), splint(), stepIntegrate(), vectorConvolution2(), and x().

int Spline::n
 

Definition at line 189 of file spline.h.

Referenced by add(), arm(), back(), checkConvolutionRange(), checkSpace(), checksum(), convOneDim(), createChecksum(), derive(), dump(), explizit(), fastSplint(), flip(), getZeroList(), indexToLeft(), indexToRight(), integrate(), last(), merge(), mul(), operator *=(), operator+=(), printStatus(), proper(), ptrY(), range(), rearm(), resize(), save(), set(), setN(), setX(), setY(), size(), smoothen(), splder(), Spline(), splint(), stepIntegrate(), stop(), x(), and y().

string Spline::name
 

the name (useful for debugging)

Definition at line 190 of file spline.h.

Referenced by dump(), expandVectorBetween(), printStatus(), proper(), setName(), and Spline().

unsigned int Spline::nvar [static]
 

The size of the convolution vector.

Definition at line 208 of file spline.h.

Referenced by convVektor2(), and vectorConvolution2().

int Spline::nvc
 

Definition at line 200 of file spline.h.

Referenced by operator()(), and printStatus().

bool Spline::own
 

if the xdat belongs to the spline, i.e. mother ==0. Own does not say, that xdat is a valid pointer, it may be 0, if for instance, mother is killed, own will be true, but the pointer will be set to zero

Definition at line 196 of file spline.h.

Referenced by arm(), checksum(), convVektor2(), createChecksum(), derive(), disarm(), error(), explizit(), fastSplint(), flip(), getData(), inBetweenY(), merge(), motherKilled(), operator()(), printStatus(), randomY(), resize(), set(), setForce(), setMother(), setN(), setX(), splder(), Spline(), splint(), stepIntegrate(), vectorConvolution2(), and x().

Spline * Spline::partner [static]
 

convolutions partner spline

Definition at line 232 of file spline.h.

Referenced by checkConvolutionRange(), convOneDim(), and convVektor2().

double * Spline::sdat
 

Definition at line 195 of file spline.h.

Referenced by arm(), createChecksum(), disarm(), fastSplint(), generateSplineXXL(), inBetweenY(), merge(), operator()(), printStatus(), Spline(), and splint().

double Spline::shift [static]
 

a shift in integration

Definition at line 233 of file spline.h.

Referenced by checkConvolutionRange(), convD(), convOneDim(), and convVektor2().

double Spline::splg [static]
 

splini and splder need it...

Definition at line 203 of file spline.h.

Referenced by Splder(), and splini().

bool Spline::splgValid
 

derive() needs the splg[] array, which is initialized if it is not valid

Definition at line 204 of file spline.h.

Referenced by Splder(), and splini().

int Spline::splintcC
 

Definition at line 202 of file spline.h.

Referenced by printStatus(), and splint().

bool Spline::valid
 

if the xdat pointer is valid, i.e. own or mother alive And valid

Definition at line 197 of file spline.h.

Referenced by arm(), checkSpace(), derive(), merge(), motherKilled(), operator()(), printStatus(), randomY(), resize(), setMother(), and splder().

bool Spline::validCache
 

Definition at line 191 of file spline.h.

Referenced by convVektor2(), disarm(), fastSplint(), generateSplineXXL(), operator()(), printStatus(), and splint().

int Spline::vc
 

Definition at line 200 of file spline.h.

Referenced by operator()(), and printStatus().

int Spline::WatchCount = 0 [static]
 

If spline is compiled with ENABLE_SPLINE_WATCH, counts number of splines alive.

Definition at line 274 of file spline.h.

Referenced by Spline().

map< string, int > Spline::WatchMap [static]
 

If compiled with ENABLE_SPLINE_WATCH, tracks number of splines with that name.

Definition at line 275 of file spline.h.

Referenced by printWatchMap(), and Spline().

double* Spline::xdat
 

Definition at line 195 of file spline.h.

Referenced by arm(), checksum(), convVektor2(), createChecksum(), derive(), explizit(), flip(), generateSplineXXL(), getData(), inBetweenY(), indexToLeft(), indexToRight(), merge(), motherKilled(), operator()(), printStatus(), randomY(), resize(), set(), setMother(), setX(), splder(), Spline(), stepIntegrate(), vectorConvolution2(), and x().

double * Spline::ydat
 

Definition at line 195 of file spline.h.

Referenced by add(), arm(), createChecksum(), derive(), dump(), explizit(), fastSplint(), generateSplineXXL(), getData(), inBetweenY(), merge(), mul(), operator *=(), operator()(), operator+=(), printStatus(), ptrY(), set(), setY(), splder(), Spline(), splint(), stepInt(), and y().

double Spline::zeroLevel [static]
 

Definition at line 209 of file spline.h.

Referenced by findZeroBetween(), and inBetweenY().


The documentation for this class was generated from the following files:
Generated on Mon May 5 17:48:36 2008 for cmbeasy by  doxygen 1.4.6