CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
EcalTBMCInfoProducer.cc
Go to the documentation of this file.
1 /*
2  * \file EcalTBMCInfoProducer.cc
3  *
4  *
5  */
6 
9 
20 
24 
25 #include "Math/GenVector/Rotation3D.h"
26 #include <CLHEP/Random/RandFlat.h>
27 
28 #include <fstream>
29 #include <iostream>
30 #include <vector>
31 
33 public:
35  explicit EcalTBMCInfoProducer(const edm::ParameterSet &ps);
36 
38  ~EcalTBMCInfoProducer() override;
39 
41  void produce(edm::Event &event, const edm::EventSetup &eventSetup) override;
42 
43 private:
44  double beamEta;
45  double beamPhi;
46  double beamTheta;
47 
49 
50  double beamXoff;
51  double beamYoff;
52 
53  double partXhodo;
54  double partYhodo;
55 
57 
58  ROOT::Math::Rotation3D *fromCMStoTB;
59 
61 };
62 
64  produces<PEcalTBInfo>();
65 
67  GenVtxToken = consumes<edm::HepMCProduct>(edm::InputTag("moduleLabelVtx", "source"));
68  double fMinEta = ps.getParameter<double>("MinEta");
69  double fMaxEta = ps.getParameter<double>("MaxEta");
70  double fMinPhi = ps.getParameter<double>("MinPhi");
71  double fMaxPhi = ps.getParameter<double>("MaxPhi");
72  beamEta = (fMaxEta + fMinEta) / 2.;
73  beamPhi = (fMaxPhi + fMinPhi) / 2.;
74  beamTheta = 2.0 * atan(exp(-beamEta));
75  beamXoff = ps.getParameter<double>("BeamMeanX");
76  beamYoff = ps.getParameter<double>("BeamMeanX");
77 
78  std::string fullMapName = CrystalMapFile.fullPath();
79  theTestMap = new EcalTBCrystalMap(fullMapName);
80  crysNumber = 0;
81 
82  double deltaEta = 999.;
83  double deltaPhi = 999.;
84  for (int cryIndex = 1; cryIndex <= EcalTBCrystalMap::NCRYSTAL; ++cryIndex) {
85  double eta = 0;
86  double phi = 0.;
87  theTestMap->findCrystalAngles(cryIndex, eta, phi);
88  if (fabs(beamEta - eta) < deltaEta && fabs(beamPhi - phi) < deltaPhi) {
89  deltaEta = fabs(beamEta - eta);
90  deltaPhi = fabs(beamPhi - phi);
91  crysNumber = cryIndex;
92  } else if (fabs(beamEta - eta) < deltaEta && fabs(beamPhi - phi) > deltaPhi) {
93  if (fabs(beamPhi - phi) < 0.017) {
94  deltaEta = fabs(beamEta - eta);
95  deltaPhi = fabs(beamPhi - phi);
96  crysNumber = cryIndex;
97  }
98  } else if (fabs(beamEta - eta) > deltaEta && fabs(beamPhi - phi) < deltaPhi) {
99  if (fabs(beamEta - eta) < 0.017) {
100  deltaEta = fabs(beamEta - eta);
101  deltaPhi = fabs(beamPhi - phi);
102  crysNumber = cryIndex;
103  }
104  }
105  }
106 
107  edm::LogVerbatim("EcalTBInfo") << "Initialize TB MC ECAL info producer with parameters: \n"
108  << "Crystal map file: " << CrystalMapFile << "\n"
109  << "Beam average eta = " << beamEta << "\n"
110  << "Beam average phi = " << beamPhi << "\n"
111  << "Corresponding to crystal number = " << crysNumber << "\n"
112  << "Beam X offset = " << beamXoff << "\n"
113  << "Beam Y offset = " << beamYoff;
114 
115  // rotation matrix to move from the CMS reference frame to the test beam one
116 
117  double xx = -cos(beamTheta) * cos(beamPhi);
118  double xy = -cos(beamTheta) * sin(beamPhi);
119  double xz = sin(beamTheta);
120 
121  double yx = sin(beamPhi);
122  double yy = -cos(beamPhi);
123  double yz = 0.;
124 
125  double zx = sin(beamTheta) * cos(beamPhi);
126  double zy = sin(beamTheta) * sin(beamPhi);
127  double zz = cos(beamTheta);
128 
129  fromCMStoTB = new ROOT::Math::Rotation3D(xx, xy, xz, yx, yy, yz, zx, zy, zz);
130 
131  // random number
133  if (!rng.isAvailable()) {
134  throw cms::Exception("Configuration") << "EcalTBMCInfoProducer requires the RandomNumberGeneratorService\n"
135  "which is not present in the configuration file. You must add the "
136  "service\n"
137  "in the configuration file or remove the modules that require it.";
138  }
139 }
140 
142 
145  CLHEP::HepRandomEngine *engine = &rng->getEngine(event.streamID());
146 
147  std::unique_ptr<PEcalTBInfo> product(new PEcalTBInfo());
148 
149  // Fill the run information
150 
151  product->setCrystal(crysNumber);
152 
153  product->setBeamDirection(beamEta, beamPhi);
154  product->setBeamOffset(beamXoff, beamYoff);
155 
156  // Compute the event x,y vertex coordinates in the beam reference system
157  // e.g. in the place orthogonal to the beam average direction
158 
159  partXhodo = partYhodo = 0.;
160 
162  event.getByToken(GenVtxToken, GenEvt);
163 
164  const HepMC::GenEvent *Evt = GenEvt->GetEvent();
165  HepMC::GenEvent::vertex_const_iterator Vtx = Evt->vertices_begin();
166 
167  math::XYZPoint eventCMSVertex((*Vtx)->position().x(), (*Vtx)->position().y(), (*Vtx)->position().z());
168 
169  LogDebug("EcalTBInfo") << "Generated vertex position = " << eventCMSVertex.x() << " " << eventCMSVertex.y() << " "
170  << eventCMSVertex.z();
171 
172  math::XYZPoint eventTBVertex = (*fromCMStoTB) * eventCMSVertex;
173 
174  LogDebug("EcalTBInfo") << "Rotated vertex position = " << eventTBVertex.x() << " " << eventTBVertex.y() << " "
175  << eventTBVertex.z();
176 
177  partXhodo = eventTBVertex.x();
178  partYhodo = eventTBVertex.y();
179 
180  product->setBeamPosition(partXhodo, partYhodo);
181 
182  // Asynchronous phase shift
183  double thisPhaseShift = CLHEP::RandFlat::shoot(engine);
184 
185  product->setPhaseShift(thisPhaseShift);
186  LogDebug("EcalTBInfo") << "Asynchronous Phaseshift = " << thisPhaseShift;
187 
188  // store the object in the framework event
189 
190  event.put(std::move(product));
191 }
192 
Log< level::Info, true > LogVerbatim
edm::EDGetTokenT< edm::HepMCProduct > GenVtxToken
static const int NCRYSTAL
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
~EcalTBMCInfoProducer() override
Destructor.
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
static const double deltaEta
Definition: CaloConstants.h:8
void produce(edm::Event &event, const edm::EventSetup &eventSetup) override
Produce digis out of raw data.
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
def move
Definition: eostools.py:511
bool isAvailable() const
Definition: Service.h:40
EcalTBCrystalMap * theTestMap
Basic2DVector< T > xy() const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void findCrystalAngles(const int thisCrysIndex, double &thisEta, double &thisPhi)
StreamID streamID() const
Definition: Event.h:98
ROOT::Math::Rotation3D * fromCMStoTB
#define LogDebug(id)
EcalTBMCInfoProducer(const edm::ParameterSet &ps)
Constructor.