CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
SurveyPxbDicer Class Reference

#include <SurveyPxbDicer.h>

Classes

struct  findParByName
 Function object for searching for a parameter in the VPSet. More...
 

Public Types

typedef SurveyPxbImage::coord_t coord_t
 
typedef unsigned int count_t
 
typedef std::vector< coord_tfidpoint_t
 
typedef SurveyPxbImage::id_t id_t
 
typedef std::pair< id_t, id_tidPair_t
 
typedef double value_t
 

Public Member Functions

std::string doDice (const fidpoint_t &fidpointvec, const idPair_t &id, const bool rotate=false)
 
void doDice (const fidpoint_t &fidpointvec, const idPair_t &id, std::ofstream &outfile, const bool rotate=false)
 
 SurveyPxbDicer ()
 
 SurveyPxbDicer (const std::vector< edm::ParameterSet > &pars, unsigned int seed)
 

Private Member Functions

value_t getParByName (const std::string &name, const std::string &par, const std::vector< edm::ParameterSet > &pars)
 
value_t ranGauss (value_t mean, value_t sigma)
 
coord_t transform (const coord_t &x, const value_t &a0, const value_t &a1, const value_t &a2, const value_t &a3)
 

Private Attributes

value_t mean_a0
 
value_t mean_a1
 
value_t mean_phi
 
value_t mean_scale
 
value_t mean_x
 
value_t mean_y
 
value_t sigma_a0
 
value_t sigma_a1
 
value_t sigma_phi
 
value_t sigma_scale
 
value_t sigma_x
 
value_t sigma_y
 

Detailed Description

Class to dice a picture from a given set of fiducial points This class has its use for toy MC simulations to validate the PXB survey

Definition at line 21 of file SurveyPxbDicer.h.

Member Typedef Documentation

Definition at line 23 of file SurveyPxbDicer.h.

typedef unsigned int SurveyPxbDicer::count_t

Definition at line 26 of file SurveyPxbDicer.h.

typedef std::vector<coord_t> SurveyPxbDicer::fidpoint_t

Definition at line 25 of file SurveyPxbDicer.h.

Definition at line 27 of file SurveyPxbDicer.h.

typedef std::pair<id_t, id_t> SurveyPxbDicer::idPair_t

Definition at line 28 of file SurveyPxbDicer.h.

typedef double SurveyPxbDicer::value_t

Definition at line 24 of file SurveyPxbDicer.h.

Constructor & Destructor Documentation

SurveyPxbDicer::SurveyPxbDicer ( )
inline

Definition at line 31 of file SurveyPxbDicer.h.

31 {};
SurveyPxbDicer::SurveyPxbDicer ( const std::vector< edm::ParameterSet > &  pars,
unsigned int  seed 
)

Definition at line 21 of file SurveyPxbDicer.cc.

References getParByName(), mean_a0, mean_a1, mean_phi, mean_scale, mean_x, mean_y, sigma_a0, sigma_a1, sigma_phi, sigma_scale, sigma_x, and sigma_y.

21  {
22  CLHEP::HepRandom::setTheSeed(seed);
23  mean_a0 = getParByName("a0", "mean", pars);
24  sigma_a0 = getParByName("a0", "sigma", pars);
25  mean_a1 = getParByName("a1", "mean", pars);
26  sigma_a1 = getParByName("a1", "sigma", pars);
27  mean_scale = getParByName("scale", "mean", pars);
28  sigma_scale = getParByName("scale", "sigma", pars);
29  mean_phi = getParByName("phi", "mean", pars);
30  sigma_phi = getParByName("phi", "sigma", pars);
31  mean_x = getParByName("x", "mean", pars);
32  sigma_x = getParByName("x", "sigma", pars);
33  mean_y = getParByName("y", "mean", pars);
34  sigma_y = getParByName("y", "sigma", pars);
35 }
value_t mean_scale
value_t sigma_scale
value_t getParByName(const std::string &name, const std::string &par, const std::vector< edm::ParameterSet > &pars)

Member Function Documentation

std::string SurveyPxbDicer::doDice ( const fidpoint_t fidpointvec,
const idPair_t id,
const bool  rotate = false 
)

Invoke the dicer

Parameters
fidpointvecvector with fiducial points where values need to be diced for and transformed to the photo fram
idPairpair of the id values

Definition at line 37 of file SurveyPxbDicer.cc.

References a0, isotrackTrainRegressor::a1, isotrackTrainRegressor::a2, funct::cos(), mean_a0, mean_a1, mean_phi, mean_scale, fireworks::p1, fireworks::p2, phi, ranGauss(), pileupReCalc_HLTpaths::scale, edm::second(), sigma_a0, sigma_a1, sigma_phi, sigma_scale, sigma_x, sigma_y, jetcorrextractor::sign(), funct::sin(), transform(), x, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by MillePedeAlignmentAlgorithm::addPxbSurvey(), and doDice().

37  {
38  // Dice the local parameters
43  const value_t a2 = scale * cos(phi);
44  const value_t a3 = scale * sin(phi);
45  const coord_t p0 = transform(fidpointvec[0], a0, a1, a2, a3);
46  const coord_t p1 = transform(fidpointvec[1], a0, a1, a2, a3);
47  const coord_t p2 = transform(fidpointvec[2], a0, a1, a2, a3);
48  const coord_t p3 = transform(fidpointvec[3], a0, a1, a2, a3);
49  const value_t sign = rotate ? -1 : 1;
50  std::ostringstream oss;
51  oss << id.first << " " << sign * ranGauss(p0.x(), sigma_x) << " " << -sign * ranGauss(p0.y(), sigma_y) << " "
52  << sign * ranGauss(p1.x(), sigma_x) << " " << -sign * ranGauss(p1.y(), sigma_y) << " " << id.second << " "
53  << sign * ranGauss(p2.x(), sigma_x) << " " << -sign * ranGauss(p2.y(), sigma_y) << " "
54  << sign * ranGauss(p3.x(), sigma_x) << " " << -sign * ranGauss(p3.y(), sigma_y) << " " << sigma_x << " "
55  << sigma_y << " " << rotate << " # MC-truth:"
56  << " a0-a3: " << a0 << " " << a1 << " " << a2 << " " << a3 << " S: " << scale << " phi: " << phi
57  << " x0: " << fidpointvec[0].x() << " " << fidpointvec[0].y() << " x1: " << fidpointvec[1].x() << " "
58  << fidpointvec[1].y() << " x2: " << fidpointvec[2].x() << " " << fidpointvec[2].y()
59  << " x3: " << fidpointvec[3].x() << " " << fidpointvec[3].y() << std::endl;
60  return oss.str();
61 }
value_t mean_scale
const TString p2
Definition: fwPaths.cc:13
double sign(double x)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
value_t sigma_scale
coord_t transform(const coord_t &x, const value_t &a0, const value_t &a1, const value_t &a2, const value_t &a3)
U second(std::pair< T, U > const &p)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
value_t ranGauss(value_t mean, value_t sigma)
const TString p1
Definition: fwPaths.cc:12
static constexpr float a0
SurveyPxbImage::coord_t coord_t
def rotate
Definition: svgfig.py:705
void SurveyPxbDicer::doDice ( const fidpoint_t fidpointvec,
const idPair_t id,
std::ofstream &  outfile,
const bool  rotate = false 
)

Definition at line 63 of file SurveyPxbDicer.cc.

References doDice().

66  {
67  outfile << doDice(fidpointvec, id, rotate);
68 }
std::string doDice(const fidpoint_t &fidpointvec, const idPair_t &id, const bool rotate=false)
def rotate
Definition: svgfig.py:705
SurveyPxbDicer::value_t SurveyPxbDicer::getParByName ( const std::string &  name,
const std::string &  par,
const std::vector< edm::ParameterSet > &  pars 
)
private

Gets parameter by name from the VPSet

Parameters
namename of the parameter to be searched for in field 'name' of the VPSet
parselects the value, i.e. mean or sigma
parsreference to VPSet

Definition at line 75 of file SurveyPxbDicer.cc.

References c, and mergeVDriftHistosByStation::name.

Referenced by SurveyPxbDicer().

77  {
78  std::vector<edm::ParameterSet>::const_iterator it;
79  it = std::find_if(pars.begin(), pars.end(), [&name](auto const &c) { return findParByName()(name, c); });
80  if (it == pars.end()) {
81  throw std::runtime_error("Parameter not found in SurveyPxbDicer::getParByName");
82  }
83  return (*it).getParameter<value_t>(par);
84 }
const edm::EventSetup & c
value_t SurveyPxbDicer::ranGauss ( value_t  mean,
value_t  sigma 
)
inlineprivate

invoke the RNG to geat a gaussian smeared value

Parameters
meanmean value
sigma

Definition at line 45 of file SurveyPxbDicer.h.

Referenced by doDice().

45 { return CLHEP::RandGauss::shoot(mean, sigma); };
SurveyPxbDicer::coord_t SurveyPxbDicer::transform ( const coord_t x,
const value_t a0,
const value_t a1,
const value_t a2,
const value_t a3 
)
private

Transform the diced values to the frame of a toy photograph

Parameters
xcoordinate to be transformed from local frame to photo frame
a0Transformation parameter, shift in x
a1Transformation parameter, shift in y
a2Transformation parameter, scale*cos(phi)
a3Transformation parameter, scale*sin(phi)

Definition at line 70 of file SurveyPxbDicer.cc.

References PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by doDice().

71  {
72  return coord_t(a0 + a2 * x.x() + a3 * x.y(), a1 - a3 * x.x() + a2 * x.y());
73 }
static constexpr float a0
SurveyPxbImage::coord_t coord_t

Member Data Documentation

value_t SurveyPxbDicer::mean_a0
private

Definition at line 45 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::mean_a1
private

Definition at line 48 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::mean_phi
private

Definition at line 50 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::mean_scale
private

Definition at line 49 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::mean_x
private

Definition at line 51 of file SurveyPxbDicer.h.

Referenced by SurveyPxbDicer().

value_t SurveyPxbDicer::mean_y
private

Definition at line 52 of file SurveyPxbDicer.h.

Referenced by SurveyPxbDicer().

value_t SurveyPxbDicer::sigma_a0
private

Definition at line 45 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::sigma_a1
private

Definition at line 48 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::sigma_phi
private

Definition at line 50 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::sigma_scale
private

Definition at line 49 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::sigma_x
private

Definition at line 51 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().

value_t SurveyPxbDicer::sigma_y
private

Definition at line 52 of file SurveyPxbDicer.h.

Referenced by doDice(), and SurveyPxbDicer().