CMS 3D CMS Logo

TkMSParameterizationBuilder.cc
Go to the documentation of this file.
2 
6 
11 
12 
16 
19 
21 
22 
23 
28 
30 
31 
32 namespace {
33  inline
35  {
36  GlobalVector zAxis = zDir.unit();
37  GlobalVector yAxis( zAxis.y(), -zAxis.x(), 0);
38  GlobalVector xAxis = yAxis.cross( zAxis);
39  return Surface::RotationType( xAxis, yAxis, zAxis);
40  }
41 
42  struct MSData {
43  int stid;
44  int lid;
45  float zi;
46  float zo;
47  float uerr;
48  float verr;
49  };
50 }
51 inline
52 std::ostream & operator<<(std::ostream & os, MSData d) {
53  os << d.stid<<'>' <<d.lid <<'|'<<d.zi<<'/'<<d.zo<<':'<<d.uerr<<'/'<<d.verr;
54  return os;
55 }
56 
57 
58 
60 public:
64  // desc.add<std::string>("ComponentName","TkMSParameterization");
65  desc.add<std::string>("navigationSchool","SimpleNavigationSchool");
66  descriptions.add("TkMSParameterizationBuilder",desc);
67  }
68 
69  using ReturnType = std::unique_ptr<TkMSParameterization>;
70  ReturnType produce(TkMSParameterizationRecord const&);
71 
73 };
74 
76  theNavigationSchoolName(pset.getParameter<std::string>("navigationSchool")){
77  setWhatProduced(this,"");
78 }
79 
82 
83  using namespace tkMSParameterization;
84 
85  ReturnType product = std::make_unique<TkMSParameterization>();
86 
87  auto & msParam = *product;
88 
89  //
91  iRecord.getRecord<NavigationSchoolRecord>().get(theNavigationSchoolName, navSchoolH);
92  TkNavigationSchool const & navSchool = *(TkNavigationSchool const *)navSchoolH.product();
93  auto const & searchGeom = navSchool.searchTracker();
94  auto const & magfield = navSchool.field();
95 
96  //
97  edm::ESHandle<Propagator> propagatorHandle;
98  iRecord.getRecord<TrackingComponentsRecord>().get("PropagatorWithMaterial", propagatorHandle);
99  auto const & ANprop = *propagatorHandle;
100 
101 
102  // error (very very small)
103  ROOT::Math::SMatrixIdentity id;
105  C *= 1.e-16;
106 
108 
109  Chi2MeasurementEstimator estimator(30.,-3.0,0.5,2.0,0.5,1.e12); // same as defauts....
110  KFUpdator kfu;
111  LocalError he(0.001*0.001,0,0.002*0.002);
112 
113 
114  // loop over lambdas
115  bool debug = false;
116  float tl=0;
117  for (decltype(nLmBins()) ib=0; ib<nLmBins(); ++ib, tl+=lmBin()) {
118 
119  float pt = 1.0f;
120  float tanlmd = tl; // 0.1f;
121  auto sinth2 = 1.f/(1.f+tanlmd*tanlmd);
122  auto sinth = std::sqrt(sinth2);
123  auto overp = sinth/pt;
124  auto pz = pt*tanlmd;
125 
126  // debug= (tl>2.34f && tl<2.55f);
127  if(debug) std::cout << tl << ' ' << pz << ' ' << 1./overp << std::endl;
128 
129 
130  float lastzz=-18.f;
131  bool goFw=false;
132  std::string loc=" Barrel";
133  for (int iz=0;iz<2; ++iz) {
134  if (iz>0) goFw=true;
135  for (float zz=lastzz; zz<18.1f; zz+=0.2f) {
136  float z = zz;
137  GlobalPoint startingPosition(0,0,z);
138 
139  constexpr int maxLayers=6;
140 
141  std::vector<MSData> mserr[maxLayers][maxLayers];
142 
143  // define propagation from/to
144  std::function<void(int, TrajectoryStateOnSurface, DetLayer const *, float, int)>
145  propagate = [&](int from, TrajectoryStateOnSurface tsos, DetLayer const * layer, float z1, int stid) {
146 
147  for (auto il=from+1; il<maxLayers; ++il) {
148 
149  auto compLayers = navSchool.nextLayers(*layer,*tsos.freeState(),alongMomentum);
150  std::stable_sort(compLayers.begin(),compLayers.end(),[](auto a, auto b){return a->seqNum()<b->seqNum();});
151  layer = nullptr;
152  for(auto it : compLayers){
153  if (it->basicComponents().empty()) {
154  //this should never happen. but better protect for it
155  edm::LogError("TkMSParameterizationBuilder") <<"a detlayer with no components: I cannot figure out a DetId from this layer. please investigate.";
156  continue;
157  }
158  if (debug) std::cout << il << (it->isBarrel() ? " Barrel" : " Forward") << " layer " << it->seqNum() << " SubDet " << it->subDetector()<< std::endl;
159  auto const & detWithState = it->compatibleDets(tsos,ANprop,estimator);
160  if(detWithState.empty()) {
161  if(debug) std::cout << "no det on this layer" << it->seqNum() << std::endl;
162  continue;
163  }
164  layer = it;
165  auto did = detWithState.front().first->geographicalId();
166  if (debug) std::cout << "arrived at " << int(did) << std::endl;
167  tsos = detWithState.front().second;
168  if(debug) std::cout << tsos.globalPosition() << ' ' << tsos.localError().positionError() << std::endl;
169 
170 
171  // save errs
172  auto globalError = ErrorFrameTransformer::transform(tsos.localError().positionError(), tsos.surface());
173  auto gp = tsos.globalPosition();
174  float r = gp.perp();
175  float errorPhi = std::sqrt(float(globalError.phierr(gp)));
176  float errorR = std::sqrt(float(globalError.rerr(gp)));
177  float errorZ = std::sqrt(float(globalError.czz()));
178  float xerr = overp*errorPhi;
179  float zerr = overp*(layer->isBarrel() ? errorZ : errorR);
180  auto zo = layer->isBarrel() ? gp.z() : r;
181  // std::cout << tanlmd << ' ' << z1 << ' ' << it->seqNum() << ':' << xerr <<'/'<<zerr << std::endl;
182  mserr[from][il-1].emplace_back(MSData{stid,it->seqNum(),z1,zo,xerr,zerr});
183 
184  if (from==0) {
185  // spawn prop from this layer
186  // constrain it to this location (relevant for layer other than very first)
188  SiPixelRecHit hitpx(tsos.localPosition(),he,1.,*detWithState.front().first,pref);
189  auto tsosl = kfu.update(tsos, hitpx);
190  auto z1l = layer->isBarrel() ? tsos.globalPosition().z() : tsos.globalPosition().perp();
191  auto stidl = layer->seqNum();
192  propagate(il,tsosl,layer,z1l, stidl);
193  }
194  break;
195  } // loop on compLayers
196  if (!layer) break;
197  // if (!goFw) lastbz=z1;
198  if(from==0) lastzz=zz;
199 
200  } // layer loop
201 
202  if (debug && mserr[from][from].empty()) {
203  std::cout << "tl " << tanlmd << loc << ' ' <<from<< std::endl;
204  for (auto il=from; il<maxLayers; ++il) { std::cout << il << ' ';
205  for ( auto const & e : mserr[from][il]) std::cout << e<<' '; // << '-' <<e.uerr*sinth <<'/'<<e.verr*sinth <<' ';
206  std::cout << std::endl;
207  }
208  }
209  }; // propagate
210 
211  float phi = 0.0415f;
212  for (int ip=0;ip<5; ++ip) {
213  phi += (3.14159f/6.f)/5.f;
214  GlobalVector startingMomentum(pt*std::sin(phi),pt*std::cos(phi),pz);
215 
216  // make TSOS happy
217  //Define starting plane
218  PlaneBuilder pb;
219  auto rot = rotation(startingMomentum);
220  auto startingPlane = pb.plane( startingPosition, rot);
221 
222  TrajectoryStateOnSurface startingStateP( GlobalTrajectoryParameters(startingPosition,
223  startingMomentum, 1, &magfield),
224  err, *startingPlane);
225  auto tsos0 = startingStateP;
226 
227  DetLayer const * layer0 = searchGeom.pixelBarrelLayers()[0];
228  if (goFw) layer0 = searchGeom.posPixelForwardLayers()[0];
229  int stid0 = layer0->seqNum();
230 
231  {
232  auto it = layer0;
233  if(debug) std::cout << "first layer " << (it->isBarrel() ? " Barrel" : " Forward") << " layer " << it->seqNum() << " SubDet " << it->subDetector()<< std::endl;
234  }
235 
236  auto const & detWithState = layer0->compatibleDets(tsos0,ANprop,estimator);
237  if(detWithState.empty()) {
238  if(debug) std::cout << "no det on first layer" << layer0->seqNum() << std::endl;
239  continue;
240  }
241  tsos0 = detWithState.front().second;
242  if(debug) std::cout << "arrived at " << int(detWithState.front().first->geographicalId()) << ' ' << tsos0.globalPosition() << ' ' << tsos0.localError().positionError() << std::endl;
243 
244  // for barrel
245  float z1l = tsos0.globalPosition().z();
246  if (goFw) {
247  loc = " Forward";
248  z1l = tsos0.globalPosition().perp();
249  }
250 
251  propagate(0,tsos0,layer0,z1l,stid0);
252 
253  } // loop on phi
254 
255  // fill
256  for (auto from=0; from<maxLayers; ++from)
257  for (auto il=from; il<maxLayers; ++il) {
258  if (mserr[from][il].empty()) continue;
259  auto stid=mserr[from][il].front().stid;
260  auto lid=mserr[from][il].front().lid;
261  auto zi = mserr[from][il].front().zi;
262  float xerr =0;
263  float zerr =0;
264  float zo=0;
265  for (auto const & e : mserr[from][il]) {
266  if (e.stid!=stid) continue;
267  lid = std::min(lid,e.lid);
268  }
269  float nn=0;
270  for (auto const & e : mserr[from][il]) {
271  if (e.stid!=stid || lid!=e.lid) continue;
272  xerr+=e.uerr; zerr+=e.verr; zo+=e.zo; ++nn;
273  }
274  xerr/=nn; zerr/=nn; zo/=nn;
275  auto iid = packLID(stid,lid);
276  auto & md = msParam.data[iid].data[ib].data;
277  if (md.empty()
278  || std::abs(xerr-md.back().uerr)>0.2f*xerr
279  || std::abs(zerr-md.back().verr)>0.2f*zerr
280  ) md.emplace_back(Elem{zi,zo,xerr,zerr});
281  } // loop on layers
282 
283  }} // loop on z
284 
285  } // loop on tanLa
286 
287  // sort
288  for (auto & e : msParam.data)
289  for (auto & d : e.second.data) std::stable_sort(d.data.begin(),d.data.end(),[](auto const & a, auto const & b){return a.vo<b.vo;});
290 
291  return product;
292 }
293 
294 
295 
296 
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:124
unsigned short packLID(unsigned int id, unsigned int od)
static GlobalError transform(const LocalError &le, const Surface &surf)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ReturnType produce(TkMSParameterizationRecord const &)
void zerr(int)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
ReturnType plane(const PositionType &pos, const RotationType &rot) const
Definition: PlaneBuilder.h:22
T y() const
Definition: PV3DBase.h:63
const GeometricSearchTracker & searchTracker() const
GlobalPoint globalPosition() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
int seqNum() const
Definition: DetLayer.h:36
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TrackingRecHit &) const override
Definition: KFUpdator.cc:169
TkMSParameterizationBuilder(edm::ParameterSet const &)
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::unique_ptr< TkMSParameterization > ReturnType
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Vector3DBase unit() const
Definition: Vector3DBase.h:57
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
const MagneticField & field() const
#define debug
Definition: HDRShower.cc:19
std::vector< const DetLayer * > nextLayers(const DetLayer &detLayer, Args &&...args) const
NavigationDirection.
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::ostream & operator<<(std::ostream &os, MSData d)
double a
Definition: hdecay.h:121
TkRotation< float > RotationType
T x() const
Definition: PV3DBase.h:62
T const * product() const
Definition: ESHandle.h:86
#define constexpr
ib
Definition: cuy.py:662
Our base class.
Definition: SiPixelRecHit.h:23