gambit is hosted by Hepforge, IPPP Durham
GAMBIT  v1.5.0-2191-ga4742ac
a Global And Modular Bsm Inference Tool

#include <daFunk.hpp>

Inheritance diagram for daFunk::FunkIntegrate_gsl1d:
Collaboration diagram for daFunk::FunkIntegrate_gsl1d:

Public Member Functions

 FunkIntegrate_gsl1d (Funk f0, std::string arg, Funk f1, Funk f2)
 
 FunkIntegrate_gsl1d (Funk f0, std::string arg, double x, Funk f)
 
 FunkIntegrate_gsl1d (Funk f0, std::string arg, double x, double y)
 
 FunkIntegrate_gsl1d (Funk f0, std::string arg, Funk f, double x)
 
 FunkIntegrate_gsl1d (Funk f0, std::string arg, std::string x, Funk f)
 
 FunkIntegrate_gsl1d (Funk f0, std::string arg, std::string x, std::string y)
 
 FunkIntegrate_gsl1d (Funk f0, std::string arg, Funk f, std::string x)
 
 FunkIntegrate_gsl1d (Funk f0, std::string arg, std::string x, double y)
 
 FunkIntegrate_gsl1d (Funk f0, std::string arg, double y, std::string x)
 
void resolve (std::map< std::string, size_t > datamap, size_t &datalen, size_t bindID, std::map< std::string, size_t > &argmap)
 
 ~FunkIntegrate_gsl1d ()
 
shared_ptr< FunkIntegrate_gsl1dset_epsrel (double epsrel)
 
shared_ptr< FunkIntegrate_gsl1dset_epsabs (double epsabs)
 
shared_ptr< FunkIntegrate_gsl1dset_limit (size_t limit)
 
shared_ptr< FunkIntegrate_gsl1dset_singularity_factor (double f)
 
shared_ptr< FunkIntegrate_gsl1dset_use_log_fallback (bool flag)
 
double value (const std::vector< double > &data, size_t bindID)
 
- Public Member Functions inherited from daFunk::FunkBase
 FunkBase ()
 
virtual ~FunkBase ()
 
template<typename... Args>
Funk set (std::string arg, Funk g, Args... args)
 
template<typename... Args>
Funk set (std::string arg, double x, Args... args)
 
template<typename... Args>
Funk set (std::string arg, std::string arg1, Args... args)
 
template<typename... Args>
Funk set ()
 
template<typename... Args>
shared_ptr< FunkBoundbind (Args... args)
 
const std::vector< std::string > & getArgs ()
 
std::size_t getNArgs ()
 
bool hasArg (std::string)
 
bool hasArgs ()
 
Funk help ()
 
template<typename... Args>
bool assert_args (Args... args)
 
Singularities getSingl ()
 
Funk set_singularity (std::string arg, double pos, double width)
 
Funk print (std::string arg)
 
template<typename... Args>
shared_ptr< FunkIntegrate_gsl1dgsl_integration (Args... args)
 
PlainPtrs1 plain (std::string)
 
PlainPtrs2 plain (std::string, std::string)
 
PlainPtrs3 plain (std::string, std::string, std::string)
 
PlainPtrs4 plain (std::string, std::string, std::string, std::string)
 
template<typename T >
PlainPtr1 plain (std::string)
 
template<typename T >
PlainPtr2 plain (std::string, std::string)
 
template<typename T >
PlainPtr3 plain (std::string, std::string, std::string)
 
template<typename T >
PlainPtr4 plain (std::string, std::string, std::string, std::string)
 
template<>
Funk set ()
 

Private Member Functions

void setup (Funk f0, std::string arg, Funk f1, Funk f2)
 

Static Private Member Functions

static double invoke (double x, void *params)
 

Private Attributes

std::vector< doublelocal_data
 
size_t local_bindID
 
std::vector< std::pair< Funk, Funk > > my_singularities
 
std::string arg
 
gsl_integration_workspace * gsl_workspace
 
size_t limit
 
std::vector< size_t > index
 
double epsrel
 
double epsabs
 
bool use_log_fallback
 
double singl_factor
 

Additional Inherited Members

- Protected Attributes inherited from daFunk::FunkBase
std::vector< Funkfunctions
 
ArgsType arguments
 
std::vector< std::vector< size_t > > indices
 
size_t datalen
 
Singularities singularities
 

Detailed Description

Definition at line 1475 of file daFunk.hpp.

Constructor & Destructor Documentation

◆ FunkIntegrate_gsl1d() [1/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
Funk  f1,
Funk  f2 
)
inline

Definition at line 1478 of file daFunk.hpp.

References Gambit::SpecBit::setup().

1479  {
1480  setup(f0, arg, f1, f2);
1481  }
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Here is the call graph for this function:

◆ FunkIntegrate_gsl1d() [2/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
double  x,
Funk  f 
)
inline

Definition at line 1482 of file daFunk.hpp.

References daFunk::cnst(), and Gambit::SpecBit::setup().

1483  {
1484  setup(f0, arg, cnst(x), f);
1485  }
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Funk cnst(double x, Args... argss)
Definition: daFunk.hpp:606
Here is the call graph for this function:

◆ FunkIntegrate_gsl1d() [3/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
double  x,
double  y 
)
inline

Definition at line 1486 of file daFunk.hpp.

References daFunk::cnst(), and Gambit::SpecBit::setup().

1487  {
1488  setup(f0, arg, cnst(x), cnst(y));
1489  }
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Funk cnst(double x, Args... argss)
Definition: daFunk.hpp:606
Here is the call graph for this function:

◆ FunkIntegrate_gsl1d() [4/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
Funk  f,
double  x 
)
inline

Definition at line 1490 of file daFunk.hpp.

References daFunk::cnst(), and Gambit::SpecBit::setup().

1491  {
1492  setup(f0, arg, f, cnst(x));
1493  }
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Funk cnst(double x, Args... argss)
Definition: daFunk.hpp:606
Here is the call graph for this function:

◆ FunkIntegrate_gsl1d() [5/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
std::string  x,
Funk  f 
)
inline

Definition at line 1494 of file daFunk.hpp.

References Gambit::SpecBit::setup(), and daFunk::var().

1495  {
1496  setup(f0, arg, var(x), f);
1497  }
Funk var(std::string arg)
Definition: daFunk.hpp:921
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Here is the call graph for this function:

◆ FunkIntegrate_gsl1d() [6/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
std::string  x,
std::string  y 
)
inline

Definition at line 1498 of file daFunk.hpp.

References Gambit::SpecBit::setup(), and daFunk::var().

1499  {
1500  setup(f0, arg, var(x), var(y));
1501  }
Funk var(std::string arg)
Definition: daFunk.hpp:921
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Here is the call graph for this function:

◆ FunkIntegrate_gsl1d() [7/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
Funk  f,
std::string  x 
)
inline

Definition at line 1502 of file daFunk.hpp.

References Gambit::SpecBit::setup(), and daFunk::var().

1503  {
1504  setup(f0, arg, f, var(x));
1505  }
Funk var(std::string arg)
Definition: daFunk.hpp:921
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Here is the call graph for this function:

◆ FunkIntegrate_gsl1d() [8/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
std::string  x,
double  y 
)
inline

Definition at line 1506 of file daFunk.hpp.

References daFunk::cnst(), Gambit::SpecBit::setup(), and daFunk::var().

1507  {
1508  setup(f0, arg, var(x), cnst(y));
1509  }
Funk var(std::string arg)
Definition: daFunk.hpp:921
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Funk cnst(double x, Args... argss)
Definition: daFunk.hpp:606
Here is the call graph for this function:

◆ FunkIntegrate_gsl1d() [9/9]

daFunk::FunkIntegrate_gsl1d::FunkIntegrate_gsl1d ( Funk  f0,
std::string  arg,
double  y,
std::string  x 
)
inline

Definition at line 1510 of file daFunk.hpp.

References daFunk::cnst(), Gambit::SpecBit::setup(), and daFunk::var().

1511  {
1512  setup(f0, arg, cnst(y), var(x));
1513  }
Funk var(std::string arg)
Definition: daFunk.hpp:921
void setup(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1658
Funk cnst(double x, Args... argss)
Definition: daFunk.hpp:606
Here is the call graph for this function:

◆ ~FunkIntegrate_gsl1d()

daFunk::FunkIntegrate_gsl1d::~FunkIntegrate_gsl1d ( )
inline

Definition at line 1560 of file daFunk.hpp.

1561  {
1562  gsl_integration_workspace_free(gsl_workspace);
1563  }
gsl_integration_workspace * gsl_workspace
Definition: daFunk.hpp:1696

Member Function Documentation

◆ invoke()

static double daFunk::FunkIntegrate_gsl1d::invoke ( double  x,
void params 
)
inlinestaticprivate

Definition at line 1681 of file daFunk.hpp.

References daFunk::FunkBase::functions, index, local_bindID, and local_data.

Referenced by value().

1681  {
1682  FunkIntegrate_gsl1d * ptr = static_cast<FunkIntegrate_gsl1d*>(params);
1683  ptr->local_data[ptr->index[ptr->local_bindID]] = x;
1684  return ptr->functions[0]->value(ptr->local_data, ptr->local_bindID);
1685  }
FunkIntegrate_gsl1d(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1478
Here is the caller graph for this function:

◆ resolve()

void daFunk::FunkIntegrate_gsl1d::resolve ( std::map< std::string, size_t >  datamap,
size_t &  datalen,
size_t  bindID,
std::map< std::string, size_t > &  argmap 
)
inlinevirtual

Reimplemented from daFunk::FunkBase.

Definition at line 1515 of file daFunk.hpp.

References combine_hdf5::index.

1516  {
1517  functions[1]->resolve(datamap, datalen, bindID, argmap); // Resolve boundary 0
1518  functions[2]->resolve(datamap, datalen, bindID, argmap); // Resolve boundary 1
1519 
1520  // Set indices
1521  if(index.size() <= bindID)
1522  {
1523  index.resize(bindID+1);
1524  // TODO: Throw error if problems are encountered
1525  //std::cout << "FATAL ERROR: bind() attempts to resolve dependencies of FunkIntegrate_gsl1d" << std::endl;
1526  //std::cout << "object inconsistently." << std::endl;
1527  //std::cout << "Encountered while resolving " << arg << " and:" << std::endl;
1528  //for (auto it = datamap.begin(); it != datamap.end(); it++)
1529  //{
1530  // std::cout << " " << it->first << std::endl;
1531  //}
1532  //exit(1);
1533  }
1534  if ( argmap.find(arg) == argmap.end() )
1535  {
1536  if(datamap.find(arg) == datamap.end())
1537  {
1538  index[bindID] = datalen;
1539  argmap[arg] = datalen;
1540  ++datalen;
1541  }
1542  else
1543  {
1544  index[bindID] = datamap[arg];
1545  }
1546  }
1547  else
1548  {
1549  index[bindID] = argmap[arg];
1550  }
1551  datamap[arg] = index[bindID];
1552  functions[0]->resolve(datamap, datalen, bindID, argmap);
1553  for ( auto it = my_singularities.begin(); it != my_singularities.end(); ++it )
1554  {
1555  it->first->resolve(datamap, datalen, bindID, argmap);
1556  it->second->resolve(datamap, datalen, bindID, argmap);
1557  }
1558  }
std::vector< Funk > functions
Definition: daFunk.hpp:368
std::vector< std::pair< Funk, Funk > > my_singularities
Definition: daFunk.hpp:1690
std::vector< size_t > index
Definition: daFunk.hpp:1698
size_t datalen
Definition: daFunk.hpp:371

◆ set_epsabs()

shared_ptr<FunkIntegrate_gsl1d> daFunk::FunkIntegrate_gsl1d::set_epsabs ( double  epsabs)
inline

Definition at line 1567 of file daFunk.hpp.

1568  { this->epsabs = epsabs; return static_pointer_cast<FunkIntegrate_gsl1d>(this->shared_from_this()); }
FunkIntegrate_gsl1d(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1478

◆ set_epsrel()

shared_ptr<FunkIntegrate_gsl1d> daFunk::FunkIntegrate_gsl1d::set_epsrel ( double  epsrel)
inline

Definition at line 1565 of file daFunk.hpp.

1566  { this->epsrel = epsrel; return static_pointer_cast<FunkIntegrate_gsl1d>(this->FunkIntegrate_gsl1d::shared_from_this()); }
FunkIntegrate_gsl1d(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1478

◆ set_limit()

shared_ptr<FunkIntegrate_gsl1d> daFunk::FunkIntegrate_gsl1d::set_limit ( size_t  limit)
inline

Definition at line 1569 of file daFunk.hpp.

1570  { this->limit = limit; return static_pointer_cast<FunkIntegrate_gsl1d>(this->shared_from_this()); }
FunkIntegrate_gsl1d(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1478

◆ set_singularity_factor()

shared_ptr<FunkIntegrate_gsl1d> daFunk::FunkIntegrate_gsl1d::set_singularity_factor ( double  f)
inline

Definition at line 1571 of file daFunk.hpp.

References f.

1572  { this->singl_factor = f; return static_pointer_cast<FunkIntegrate_gsl1d>(this->shared_from_this()); }
FunkIntegrate_gsl1d(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1478

◆ set_use_log_fallback()

shared_ptr<FunkIntegrate_gsl1d> daFunk::FunkIntegrate_gsl1d::set_use_log_fallback ( bool  flag)
inline

Definition at line 1573 of file daFunk.hpp.

1574  { this->use_log_fallback = flag; return static_pointer_cast<FunkIntegrate_gsl1d>(this->shared_from_this()); }
FunkIntegrate_gsl1d(Funk f0, std::string arg, Funk f1, Funk f2)
Definition: daFunk.hpp:1478

◆ setup()

void daFunk::FunkIntegrate_gsl1d::setup ( Funk  f0,
std::string  arg,
Funk  f1,
Funk  f2 
)
inlineprivate

Definition at line 1658 of file daFunk.hpp.

References daFunk::eraseArg(), daFunk::joinArgs(), daFunk::joinSingl(), and daFunk::vec().

1659  {
1660  this->functions = vec(f0, f1, f2);
1661 
1662  singularities = joinSingl(f1->getSingl(), f2->getSingl());
1663  if ( f0->getSingl().find(arg) != f0->getSingl().end() )
1664  my_singularities = f0->getSingl()[arg];
1665  Singularities tmp_singl = f0->getSingl();
1666  tmp_singl.erase(arg);
1667  singularities = joinSingl(singularities, tmp_singl);
1668 
1669  arguments = joinArgs(eraseArg(f0->getArgs(), arg), joinArgs(f1->getArgs(), f2->getArgs()));
1670  gsl_workspace = gsl_integration_workspace_alloc (100000);
1671 
1672  this->arg = arg;
1673  limit = 100;
1674  epsrel = 1e-2;
1675  epsabs = 1e-2;
1676  use_log_fallback = false;
1677  singl_factor = 4;
1678  }
ArgsType arguments
Definition: daFunk.hpp:369
ArgsType eraseArg(ArgsType args, std::string arg)
Definition: daFunk.hpp:236
gsl_integration_workspace * gsl_workspace
Definition: daFunk.hpp:1696
std::map< std::string, std::vector< std::pair< Funk, Funk > > > Singularities
Definition: daFunk.hpp:118
std::vector< Funk > functions
Definition: daFunk.hpp:368
ArgsType joinArgs(ArgsType args1, ArgsType args2)
Definition: daFunk.hpp:228
std::vector< std::pair< Funk, Funk > > my_singularities
Definition: daFunk.hpp:1690
Singularities singularities
Definition: daFunk.hpp:372
std::vector< T > vec(std::vector< T > vector)
Definition: daFunk.hpp:142
Singularities joinSingl(Singularities s1, Singularities s2)
Definition: daFunk.hpp:243
Here is the call graph for this function:

◆ value()

double daFunk::FunkIntegrate_gsl1d::value ( const std::vector< double > &  data,
size_t  bindID 
)
inlinevirtual

Implements daFunk::FunkBase.

Definition at line 1576 of file daFunk.hpp.

References Gambit::GreAT::data, combine_hdf5::index, invoke(), daFunk::logspace(), generate_raster_scan_settings::N, and Gambit::DecayBit::SM_Z::sigma.

1577  {
1578  double result;
1579  #pragma omp critical(FunkIntegrate_gsl1d_integration)
1580  {
1581  local_data = data;
1582  local_bindID = bindID;
1583  double error;
1584  function=&FunkIntegrate_gsl1d::invoke;
1585  params=this;
1586  double x0 = functions[1]->value(data, bindID);
1587  double x1 = functions[2]->value(data, bindID);
1588  gsl_set_error_handler_off();
1589  int status = 0;
1590  if ( my_singularities.size() == 0 )
1591  {
1592  status = gsl_integration_qags(this, x0, x1, epsabs, epsrel, limit, gsl_workspace, &result, &error);
1593  }
1594  else
1595  {
1596  double s = 0;
1597  std::vector<double> ranges;
1598  ranges.push_back(x0);
1599  ranges.push_back(x1);
1600  for ( auto it = my_singularities.begin(); it != my_singularities.end(); ++it )
1601  {
1602  double mean = it->first->value(data, bindID);
1603  double sigma = it->second->value(data, bindID);
1604  double z0 = mean - singl_factor*sigma;
1605  double z1 = mean + singl_factor*sigma;
1606  if ( z0 == z1 )
1607  std::cout << "daFunk::FunkBase WARNING: Singularity width is beyond machine precision." << std::endl;
1608  if ( z0 > x0 and z0 < x1 ) ranges.push_back(z0);
1609  if ( z1 > x0 and z1 < x1 ) ranges.push_back(z1);
1610  }
1611  std::sort(ranges.begin(), ranges.end());
1612  for ( auto it = ranges.begin(); it != ranges.end()-1; ++it )
1613  {
1614  status = gsl_integration_qags(this, *it, *(it+1), epsabs, epsrel, limit, gsl_workspace, &result, &error);
1615  s += result;
1616  if (status) break;
1617  }
1618  result = s;
1619  }
1620  if (status and this->use_log_fallback)
1621  {
1622  // The last resort: A cheap integration on log grid, linear interpolation
1623  const double N = 300;
1624  std::vector<double> Xgrid =
1625  logspace(std::log10(x0), std::log10(x1), N);
1626  double sum = 0, y0, y1, dx;
1627  local_data[index[bindID]] = Xgrid[0];
1628  y0 = functions[0]->value(local_data, bindID);
1629  for (size_t i = 0; i<N-1; i++)
1630  {
1631  local_data[index[bindID]] = Xgrid[i+1];
1632  y1 = functions[0]->value(local_data, bindID);
1633  dx = Xgrid[i+1]-Xgrid[i];
1634  sum += dx*(y0+y1)/2;
1635  y0 = y1;
1636  }
1637  result = sum;
1638  }
1639  // TODO: Implement flags to optionally throw an error
1640  if (status and not this->use_log_fallback)
1641  {
1642  std::cerr << "daFunk::FunkIntegrate_gsl1d WARNING: " << gsl_strerror(status) << std::endl;
1643  std::cerr << "Attempt to integrate from " << x0 << " to " << x1 << std::endl;
1644  std::cerr << "Attempt to integrate from " << x0 << " to " << x1 << std::endl;
1645  std::cerr << "Details about the integrand:" << std::endl;
1646  functions[0]->help();
1647 // std::cout << "Dumping integrand:" << std::endl;
1648 // for ( double x = x0; x <= x1; x = (x0>0) ? x*1.01 : x+(x1-x0)/1000)
1649 // std::cerr << " " << x << " " << invoke(x, this) << std::endl;
1650  std::cerr << "Returning zero." << std::endl;
1651  result = 0.;
1652  }
1653  }
1654  return result;
1655  }
greatScanData data
Definition: great.cpp:38
std::vector< double > local_data
Definition: daFunk.hpp:1688
gsl_integration_workspace * gsl_workspace
Definition: daFunk.hpp:1696
const double sigma
Definition: SM_Z.hpp:43
std::vector< Funk > functions
Definition: daFunk.hpp:368
std::vector< double > logspace(double x0, double x1, unsigned int n)
Definition: daFunk.hpp:186
std::vector< std::pair< Funk, Funk > > my_singularities
Definition: daFunk.hpp:1690
std::vector< size_t > index
Definition: daFunk.hpp:1698
static double invoke(double x, void *params)
Definition: daFunk.hpp:1681
Here is the call graph for this function:

Member Data Documentation

◆ arg

std::string daFunk::FunkIntegrate_gsl1d::arg
private

Definition at line 1693 of file daFunk.hpp.

◆ epsabs

double daFunk::FunkIntegrate_gsl1d::epsabs
private

Definition at line 1700 of file daFunk.hpp.

◆ epsrel

double daFunk::FunkIntegrate_gsl1d::epsrel
private

Definition at line 1699 of file daFunk.hpp.

◆ gsl_workspace

gsl_integration_workspace* daFunk::FunkIntegrate_gsl1d::gsl_workspace
private

Definition at line 1696 of file daFunk.hpp.

◆ index

std::vector<size_t> daFunk::FunkIntegrate_gsl1d::index
private

Definition at line 1698 of file daFunk.hpp.

Referenced by invoke().

◆ limit

size_t daFunk::FunkIntegrate_gsl1d::limit
private

Definition at line 1697 of file daFunk.hpp.

◆ local_bindID

size_t daFunk::FunkIntegrate_gsl1d::local_bindID
private

Definition at line 1689 of file daFunk.hpp.

Referenced by invoke().

◆ local_data

std::vector<double> daFunk::FunkIntegrate_gsl1d::local_data
private

Definition at line 1688 of file daFunk.hpp.

Referenced by invoke().

◆ my_singularities

std::vector<std::pair<Funk, Funk> > daFunk::FunkIntegrate_gsl1d::my_singularities
private

Definition at line 1690 of file daFunk.hpp.

◆ singl_factor

double daFunk::FunkIntegrate_gsl1d::singl_factor
private

Definition at line 1703 of file daFunk.hpp.

◆ use_log_fallback

bool daFunk::FunkIntegrate_gsl1d::use_log_fallback
private

Definition at line 1701 of file daFunk.hpp.


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