CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
FP420RecoMain.cc
Go to the documentation of this file.
1 // File: FP420RecoMain.cc
3 // Date: 12.2006
4 // Description: FP420RecoMain for FP420
5 // Modifications:
7 #include <vector>
8 #include <iostream>
11 
13 
19 
20 // #include "CLHEP/Vector/LorentzVector.h"
21 // #include "CLHEP/Random/RandFlat.h"
22 #include <cmath>
23 
24 using namespace std;
25 
26 FP420RecoMain::FP420RecoMain(const edm::ParameterSet &conf) : conf_(conf) {
27  verbosity = conf_.getUntrackedParameter<int>("VerbosityLevel");
28  m_rpp420_f = conf_.getParameter<double>("RP420f"); //mm
29  m_rpp420_b = conf_.getParameter<double>("RP420b"); //mm
30  m_zreff = conf_.getParameter<double>("zreff"); //mm
31  m_zrefb = conf_.getParameter<double>("zrefb"); //mm
32  dn0 = conf_.getParameter<int>("NumberFP420Detectors");
33 
34  if (verbosity > 0) {
35  std::cout << "FP420RecoMain constructor::" << std::endl;
36  std::cout << "m_rpp420_f=" << m_rpp420_f << " m_rpp420_b=" << m_rpp420_b << std::endl;
37  std::cout << "m_zreff=" << m_zreff << " m_zrefb=" << m_zrefb << std::endl;
38  }
39 
40  double eee1 = 11.;
41  double eee2 = 12.;
42  // zinibeg_ = (eee1-eee2)/2.;
43  zinibeg_ = 0.;
44  //
45  if (verbosity > 1) {
46  std::cout << "FP420RecoMain constructor::" << std::endl;
47  std::cout << " eee1=" << eee1 << " eee2=" << eee2 << " zinibeg =" << zinibeg_ << std::endl;
48  }
51 }
52 
54  if (finderParameters_ != nullptr) {
55  delete finderParameters_;
56  }
57 }
58 
60  edm::Handle<TrackCollectionFP420> &input, RecoCollectionFP420 *toutput, double VtxX, double VtxY, double VtxZ) {
61  // initialization
62  bool first = true;
63  // finderParameters_->clear();
64  // finderParameters_->setIP( 0., 0., 0. );
65  std::vector<TrackFP420> rhits;
66  int restracks = 10; // max # tracks
67  rhits.reserve(restracks);
68  rhits.clear();
69 
70  // loop over detunits:
71  for (int number_detunits = 1; number_detunits < dn0; number_detunits++) {
72  unsigned int StID = number_detunits;
73  std::vector<RecoFP420> rcollector;
74  int restracks = 10; // max # tracks
75  rcollector.reserve(restracks);
76  rcollector.clear();
77 
79  std::vector<TrackFP420> collector;
80  collector.clear();
81  TrackCollectionFP420::Range outputRange;
82  unsigned int StIDTrack = 1111;
83  double z420 = m_rpp420_f;
84  double zref1 = m_zreff;
85  double zinibeg = zinibeg_;
86  double VtxXcur = VtxX;
87  double VtxYcur = VtxY;
88  double VtxZcur = VtxZ;
89  if (StID == 2) {
90  StIDTrack = 2222;
91  z420 = -m_rpp420_b;
92  zref1 = -m_zrefb;
93  zinibeg = -zinibeg_;
94  // VtxXcur = -VtxX;
95  // VtxYcur = -VtxY;
96  // VtxZcur = -VtxZ;
97  }
98  double z1 = z420 + zinibeg - VtxZcur;
99  double z2 = z420 + zinibeg + zref1 - VtxZcur;
100  if (verbosity > 1) {
101  std::cout << "FP420RecoMain: StIDTrack=" << StIDTrack << std::endl;
102  }
103  outputRange = input->get(StIDTrack);
104  //
105  // fill output in collector vector (for may be sorting? or other checks)
106  //
107  TrackCollectionFP420::ContainerIterator sort_begin = outputRange.first;
108  TrackCollectionFP420::ContainerIterator sort_end = outputRange.second;
109  //
110  for (; sort_begin != sort_end; ++sort_begin) {
111  collector.push_back(*sort_begin);
112  } // for sort_begin
113  if (verbosity > 1) {
114  std::cout << "FP420RecoMain: track collector.size=" << collector.size() << std::endl;
115  }
116  std::vector<TrackFP420>::const_iterator simHitIter = collector.begin();
117  std::vector<TrackFP420>::const_iterator simHitIterEnd = collector.end();
118  for (; simHitIter != simHitIterEnd; ++simHitIter) {
119  const TrackFP420 itrack = *simHitIter;
120  double x1 = (itrack.bx() * z1 + (itrack.ax() - VtxXcur)) * 1000.; //um
121  double y1 = (itrack.by() * z1 + (itrack.ay() - VtxYcur)) * 1000.; //um
122  double x2 = (itrack.bx() * z2 + (itrack.ax() - VtxXcur)) * 1000.; //um
123  double y2 = (itrack.by() * z2 + (itrack.ay() - VtxYcur)) * 1000.; //um
125  if (verbosity == -49) {
126  std::cout << "==================================================================== " << std::endl;
127  std::cout << "FP420RecoMain: StID= " << StID << std::endl;
128  std::cout << "input coord. in mm: z1= " << z1 << std::endl;
129  std::cout << "input coord. in mm: z2= " << z2 << std::endl;
130  std::cout << "input: itrack.bx()= " << itrack.bx() << std::endl;
131  std::cout << "input: itrack.ax()= " << itrack.ax() << std::endl;
132  std::cout << "input: itrack.by()= " << itrack.by() << std::endl;
133  std::cout << "input: itrack.ay()= " << itrack.ay() << std::endl;
134 
135  std::cout << "input: in um X1noVtx= " << (itrack.bx() * (z420 + zinibeg) + itrack.ax()) * 1000. << std::endl;
136  std::cout << "input: in um Y1noVtx= " << (itrack.by() * (z420 + zinibeg) + itrack.ay()) * 1000. << std::endl;
137  std::cout << "input: in um X2noVtx= " << (itrack.bx() * (z420 + zinibeg + zref1) + itrack.ax()) * 1000.
138  << std::endl;
139  std::cout << "input: in um Y2noVtx= " << (itrack.by() * (z420 + zinibeg + zref1) + itrack.ay()) * 1000.
140  << std::endl;
141 
142  std::cout << "input: in mm VtxXcur= " << VtxXcur << std::endl;
143  std::cout << "input: in mm VtxYcur= " << VtxYcur << std::endl;
144  std::cout << "input: in mm VtxZcur= " << VtxZcur << std::endl;
145  std::cout << "input coord. in um: x1= " << x1 << std::endl;
146  std::cout << "input coord. in um: y1= " << y1 << std::endl;
147  std::cout << "input coord. in um: x2= " << x2 << std::endl;
148  std::cout << "input coord. in um: y2= " << y2 << std::endl;
149  }
150  double zz1 = fabs(z1);
151  double zz2 = fabs(z2);
152  rcollector = finderParameters_->reconstruct(
153  StID, x1, y1, x2, y2, zz1, zz2); // x1,y1,x2,y2 input coord. in um; z1, z2 in mm
154  }
156 
157  if (verbosity > 1) {
158  std::cout << "FP420RecoMain: track rcollector.size=" << rcollector.size() << std::endl;
159  }
160  if (!rcollector.empty()) {
161  RecoCollectionFP420::Range rinputRange;
162  rinputRange.first = rcollector.begin();
163  rinputRange.second = rcollector.end();
164 
165  if (first) {
166  // use it only if RecoCollectionFP420 is the RecoCollection of one event, otherwise, do not use (loose 1st cl. of 1st event only)
167  first = false;
168  unsigned int StID0 = 0;
169  toutput->put(rinputRange, StID0); // !!! put into adress 0 for detID which will not be used never
170  if (verbosity > 1) {
171  std::cout << "FP420RecoMain: put(rinputRange,StID0)" << std::endl;
172  }
173  } //if ( first )
174 
175  // !!! put !!! put
176  toutput->put(rinputRange, StID);
177  if (verbosity > 1) {
178  std::cout << "FP420RecoMain: put(rinputRange,StID)" << std::endl;
179  }
180 
181  } // if rcollector.size
182 
183  } // for loop over detunits
184 
185  if (verbosity > 1) {
186  // check of access to the zcollector:
187  for (int number_detunits = 1; number_detunits < dn0; number_detunits++) {
188  int StID = number_detunits;
189  if (verbosity > 1) {
190  std::cout << " ===" << std::endl;
191  std::cout << " ===" << std::endl;
192  std::cout << "FP420RecoMain: re-new StID= " << StID << std::endl;
193  }
194  std::vector<RecoFP420> zcollector;
195  zcollector.clear();
196  RecoCollectionFP420::Range zoutputRange;
197  zoutputRange = toutput->get(StID);
198  // fill output in zcollector vector (for may be sorting? or other checks)
199  RecoCollectionFP420::ContainerIterator sort_begin = zoutputRange.first;
200  RecoCollectionFP420::ContainerIterator sort_end = zoutputRange.second;
201  for (; sort_begin != sort_end; ++sort_begin) {
202  zcollector.push_back(*sort_begin);
203  } // for
204  std::cout << "=======FP420RecoMain:check of re-new zcollector size = " << zcollector.size() << std::endl;
205  std::cout << " ===" << std::endl;
206  std::cout << " ===" << std::endl;
207  std::vector<RecoFP420>::const_iterator simHitIter = zcollector.begin();
208  std::vector<RecoFP420>::const_iterator simHitIterEnd = zcollector.end();
209  // loop in recoess
210  for (; simHitIter != simHitIterEnd; ++simHitIter) {
211  const RecoFP420 itrack = *simHitIter;
212  //double e0, double x0, double y0, double tx0, double ty0, double q20, int direction
213  std::cout << "FP420RecoMain:check: direction = " << itrack.direction() << std::endl;
214  std::cout << " e0 = " << itrack.e0() << " q20 = " << itrack.q20() << std::endl;
215  std::cout << " tx0 = " << itrack.tx0() << " ty0 = " << itrack.ty0() << std::endl;
216  std::cout << " x0= " << itrack.x0() << " y0= " << itrack.y0() << std::endl;
217  std::cout << " ===" << std::endl;
218  std::cout << " ===" << std::endl;
219  std::cout << " =======================" << std::endl;
220  }
221  }
222 
223  //==================================
224 
225  // end of check of access to the strip collection
226  std::cout << "======= FP420RecoMain: end of check " << std::endl;
227 
228  } //if (verbos
229 }
T getUntrackedParameter(std::string const &, T const &) const
std::pair< ContainerIterator, ContainerIterator > Range
double ty0() const
Definition: RecoFP420.h:18
FP420RecoMain(const edm::ParameterSet &conf)
edm::ParameterSet conf_
Definition: FP420RecoMain.h:32
double ax() const
Definition: TrackFP420.h:14
double by() const
Definition: TrackFP420.h:19
std::pair< ContainerIterator, ContainerIterator > Range
double m_rpp420_f
Definition: FP420RecoMain.h:37
static std::string const input
Definition: EdmProvDump.cc:47
std::vector< RecoFP420 >::const_iterator ContainerIterator
double q20() const
Definition: RecoFP420.h:19
double ay() const
Definition: TrackFP420.h:18
std::vector< TrackFP420 >::const_iterator ContainerIterator
RecoProducerFP420 * finderParameters_
Definition: FP420RecoMain.h:33
double bx() const
Definition: TrackFP420.h:15
const Range get(unsigned int stationID) const
void run(edm::Handle< TrackCollectionFP420 > &input, RecoCollectionFP420 *toutput, double VtxX, double VtxY, double VtxZ)
Runs the algorithm.
int direction() const
Definition: RecoFP420.h:20
void put(Range input, unsigned int stationID)
double m_rpp420_b
Definition: FP420RecoMain.h:38
std::vector< RecoFP420 > reconstruct(int, double, double, double, double, double, double)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double x0() const
Definition: RecoFP420.h:15
double tx0() const
Definition: RecoFP420.h:17
tuple cout
Definition: gather_cfg.py:144
double y0() const
Definition: RecoFP420.h:16
double e0() const
Definition: RecoFP420.h:14