CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/Alignment/SurveyAnalysis/src/SurveyPxbDicer.cc

Go to the documentation of this file.
00001 #include "Alignment/SurveyAnalysis/interface/SurveyPxbImage.h"
00002 #include "Alignment/SurveyAnalysis/interface/SurveyPxbDicer.h"
00003 
00004 //#include <stdexcept>
00005 #include <vector>
00006 #include <cmath>
00007 #include <string>
00008 #include <sstream>
00009 #include <functional>
00010 #include <algorithm>
00011 #include <fstream>
00012 #include "DataFormats/GeometryVector/interface/Point3DBase.h"
00013 #include "Math/SMatrix.h"
00014 #include "Math/SVector.h"
00015 #include <CLHEP/Random/RandGauss.h>
00016 #include <CLHEP/Random/Random.h>
00017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00018 
00019 #include <iostream>
00020 
00021 SurveyPxbDicer::SurveyPxbDicer(std::vector<edm::ParameterSet> pars, unsigned int seed)
00022 {
00023         CLHEP::HepRandom::setTheSeed( seed );
00024         mean_a0 = getParByName("a0", "mean", pars);
00025         sigma_a0 = getParByName("a0", "sigma", pars);
00026         mean_a1 = getParByName("a1", "mean", pars);
00027         sigma_a1 = getParByName("a1", "sigma", pars);
00028         mean_scale = getParByName("scale", "mean", pars);
00029         sigma_scale = getParByName("scale", "sigma", pars);
00030         mean_phi = getParByName("phi", "mean", pars);
00031         sigma_phi = getParByName("phi", "sigma", pars);
00032         mean_x = getParByName("x", "mean", pars);
00033         sigma_x = getParByName("x", "sigma", pars);
00034         mean_y = getParByName("y", "mean", pars);
00035         sigma_y = getParByName("y", "sigma", pars);
00036 }
00037 
00038 std::string SurveyPxbDicer::doDice(const fidpoint_t &fidpointvec, const idPair_t &id, const bool rotate)
00039 {
00040         // Dice the local parameters
00041         const value_t a0 = ranGauss(mean_a0, sigma_a0);
00042         const value_t a1 = ranGauss(mean_a1, sigma_a1);
00043         const value_t scale = ranGauss(mean_scale, sigma_scale);
00044         const value_t phi = ranGauss(mean_phi, sigma_phi);
00045         const value_t a2 = scale * cos(phi);
00046         const value_t a3 = scale * sin(phi);
00047         const coord_t p0 = transform(fidpointvec[0],a0,a1,a2,a3);
00048         const coord_t p1 = transform(fidpointvec[1],a0,a1,a2,a3);
00049         const coord_t p2 = transform(fidpointvec[2],a0,a1,a2,a3);
00050         const coord_t p3 = transform(fidpointvec[3],a0,a1,a2,a3);
00051         const value_t sign = rotate ? -1 : 1;
00052         std::ostringstream oss;
00053         oss << id.first << " " 
00054                  << sign*ranGauss(p0.x(),sigma_x) << " " << -sign*ranGauss(p0.y(),sigma_y) << " "
00055                  << sign*ranGauss(p1.x(),sigma_x) << " " << -sign*ranGauss(p1.y(),sigma_y) << " "
00056                  << id.second << " "
00057                  << sign*ranGauss(p2.x(),sigma_x) << " " << -sign*ranGauss(p2.y(),sigma_y) << " "
00058                  << sign*ranGauss(p3.x(),sigma_x) << " " << -sign*ranGauss(p3.y(),sigma_y) << " "
00059                  << sigma_x << " " << sigma_y << " "
00060                  << rotate
00061                  << " # MC-truth:"
00062                  << " a0-a3: " << a0 << " " << a1 << " " << a2 << " " << a3
00063                  << " S: " << scale
00064                  << " phi: " << phi
00065                  << " x0: " << fidpointvec[0].x() << " " << fidpointvec[0].y()
00066                  << " x1: " << fidpointvec[1].x() << " " << fidpointvec[1].y()
00067                  << " x2: " << fidpointvec[2].x() << " " << fidpointvec[2].y()
00068                  << " x3: " << fidpointvec[3].x() << " " << fidpointvec[3].y()
00069                  << std::endl;
00070         return oss.str();
00071 }
00072 
00073 void SurveyPxbDicer::doDice(const fidpoint_t &fidpointvec, const idPair_t &id, std::ofstream &outfile, const bool rotate)
00074 {
00075         outfile << doDice(fidpointvec, id, rotate);
00076 }
00077 
00078 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)
00079 {
00080         return coord_t(a0+a2*x.x()+a3*x.y(), a1-a3*x.x()+a2*x.y());
00081 }
00082 
00083 SurveyPxbDicer::value_t SurveyPxbDicer::getParByName(const std::string &name, const std::string &par, const std::vector<edm::ParameterSet>& pars)
00084 {
00085         std::vector<edm::ParameterSet>::const_iterator it;
00086         it = std::find_if(pars.begin(), pars.end(), std::bind1st(findParByName(),name));
00087         if (it == pars.end()) { throw std::runtime_error("Parameter not found in SurveyPxbDicer::getParByName"); }
00088         return (*it).getParameter<value_t>(par);
00089 }
00090