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