CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SurveyPxbDicer.cc
Go to the documentation of this file.
3 
4 //#include <stdexcept>
5 #include <vector>
6 #include <cmath>
7 #include <string>
8 #include <sstream>
9 #include <functional>
10 #include <algorithm>
11 #include <fstream>
13 #include "Math/SMatrix.h"
14 #include "Math/SVector.h"
15 #include "CLHEP/Random/RandGauss.h"
16 #include "CLHEP/Random/Random.h"
18 
19 #include <iostream>
20 
21 SurveyPxbDicer::SurveyPxbDicer(const std::vector<edm::ParameterSet>& pars, unsigned int seed)
22 {
23  CLHEP::HepRandom::setTheSeed( seed );
24  mean_a0 = getParByName("a0", "mean", pars);
25  sigma_a0 = getParByName("a0", "sigma", pars);
26  mean_a1 = getParByName("a1", "mean", pars);
27  sigma_a1 = getParByName("a1", "sigma", pars);
28  mean_scale = getParByName("scale", "mean", pars);
29  sigma_scale = getParByName("scale", "sigma", pars);
30  mean_phi = getParByName("phi", "mean", pars);
31  sigma_phi = getParByName("phi", "sigma", pars);
32  mean_x = getParByName("x", "mean", pars);
33  sigma_x = getParByName("x", "sigma", pars);
34  mean_y = getParByName("y", "mean", pars);
35  sigma_y = getParByName("y", "sigma", pars);
36 }
37 
38 std::string SurveyPxbDicer::doDice(const fidpoint_t &fidpointvec, const idPair_t &id, const bool rotate)
39 {
40  // Dice the local parameters
41  const value_t a0 = ranGauss(mean_a0, sigma_a0);
42  const value_t a1 = ranGauss(mean_a1, sigma_a1);
45  const value_t a2 = scale * cos(phi);
46  const value_t a3 = scale * sin(phi);
47  const coord_t p0 = transform(fidpointvec[0],a0,a1,a2,a3);
48  const coord_t p1 = transform(fidpointvec[1],a0,a1,a2,a3);
49  const coord_t p2 = transform(fidpointvec[2],a0,a1,a2,a3);
50  const coord_t p3 = transform(fidpointvec[3],a0,a1,a2,a3);
51  const value_t sign = rotate ? -1 : 1;
52  std::ostringstream oss;
53  oss << id.first << " "
54  << sign*ranGauss(p0.x(),sigma_x) << " " << -sign*ranGauss(p0.y(),sigma_y) << " "
55  << sign*ranGauss(p1.x(),sigma_x) << " " << -sign*ranGauss(p1.y(),sigma_y) << " "
56  << id.second << " "
57  << sign*ranGauss(p2.x(),sigma_x) << " " << -sign*ranGauss(p2.y(),sigma_y) << " "
58  << sign*ranGauss(p3.x(),sigma_x) << " " << -sign*ranGauss(p3.y(),sigma_y) << " "
59  << sigma_x << " " << sigma_y << " "
60  << rotate
61  << " # MC-truth:"
62  << " a0-a3: " << a0 << " " << a1 << " " << a2 << " " << a3
63  << " S: " << scale
64  << " phi: " << phi
65  << " x0: " << fidpointvec[0].x() << " " << fidpointvec[0].y()
66  << " x1: " << fidpointvec[1].x() << " " << fidpointvec[1].y()
67  << " x2: " << fidpointvec[2].x() << " " << fidpointvec[2].y()
68  << " x3: " << fidpointvec[3].x() << " " << fidpointvec[3].y()
69  << std::endl;
70  return oss.str();
71 }
72 
73 void SurveyPxbDicer::doDice(const fidpoint_t &fidpointvec, const idPair_t &id, std::ofstream &outfile, const bool rotate)
74 {
75  outfile << doDice(fidpointvec, id, rotate);
76 }
77 
78 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)
79 {
80  return coord_t(a0+a2*x.x()+a3*x.y(), a1-a3*x.x()+a2*x.y());
81 }
82 
83 SurveyPxbDicer::value_t SurveyPxbDicer::getParByName(const std::string &name, const std::string &par, const std::vector<edm::ParameterSet>& pars)
84 {
85  std::vector<edm::ParameterSet>::const_iterator it;
86  it = std::find_if(pars.begin(), pars.end(), std::bind1st(findParByName(),name));
87  if (it == pars.end()) { throw std::runtime_error("Parameter not found in SurveyPxbDicer::getParByName"); }
88  return (*it).getParameter<value_t>(par);
89 }
90 
value_t mean_scale
Function object for searching for a parameter in the VPSet.
double sign(double x)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T y() const
Definition: PV3DBase.h:63
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)
T x() const
Cartesian x coordinate.
std::string doDice(const fidpoint_t &fidpointvec, const idPair_t &id, const bool rotate=false)
list outfile
Definition: EdgesToViz.py:91
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
value_t ranGauss(value_t mean, value_t sigma)
double p2[4]
Definition: TauolaWrapper.h:90
std::vector< coord_t > fidpoint_t
SurveyPxbImage::coord_t coord_t
def rotate
Definition: svgfig.py:704
value_t getParByName(const std::string &name, const std::string &par, const std::vector< edm::ParameterSet > &pars)
std::pair< id_t, id_t > idPair_t
double p1[4]
Definition: TauolaWrapper.h:89
T x() const
Definition: PV3DBase.h:62
double p3[4]
Definition: TauolaWrapper.h:91