ATLAS Offline Software
Loading...
Searching...
No Matches
TRT_TrackSegmentsMaker_BarrelCosmics.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6// TRT_TrackSegmentsMaker_BarrelCosmics.h, (c) ATLAS Detector software
8
10
17#include <cmath>
18
20// standard methods: constructor, initialize, finalize
22
24(const std::string& t,const std::string& n,const IInterface* p)
25 : AthAlgTool(t,n,p)
26{
27 declareInterface<ITRT_TrackSegmentsMaker>(this);
28}
29
31
32 ATH_MSG_DEBUG("InDet::TRT_TrackSegmentsMaker_BarrelCosmics::initialize(), March 2012"
33 << ", magnetic field: " << (m_magneticField?"ON":"OFF")
34 << " search bins: " << m_nBinsInX << ", " << m_nBinsInPhi);
35
36 StatusCode sc = StatusCode::SUCCESS;
37
38 // Initialize ReadHandle
39 ATH_CHECK(m_driftCirclesName.initialize());
40 // TRT
41 if (detStore()->retrieve(m_trtid, "TRT_ID").isFailure()) {
42 ATH_MSG_FATAL("Could not get TRT ID helper");
43 return StatusCode::FAILURE;
44 }
45
46 if (m_minHitsForSeed<=0) m_minHitsForSeed = (int) ( 0.601 * m_minHitsForSegment );
48
49 ATH_MSG_DEBUG("m_minHitsForSegment = " << m_minHitsForSegment);
50 ATH_MSG_DEBUG("m_minHitsForSeed = " << m_minHitsForSeed);
51 ATH_MSG_DEBUG("m_minHitsAboveTOT = " << m_minHitsAboveTOT);
52
54 ATH_MSG_WARNING("initialize(): are you sure about the MinimalTOTForSeedSearch setting? (set at " << m_minSeedTOT << ")");
55
56 ATH_MSG_DEBUG("InDet::TRT_TrackSegmentsMaker_BarrelCosmics::initialize(), jobProperties: "
57 << "MinimalNumberOfTRTHits " << m_minHitsForSegment << ", MinimalTOTForSeedSearch: " << m_minSeedTOT
58 << ", m_minHitsForSeed: " << m_minHitsForSeed << ", m_minHitsAboveTOT: " << m_minHitsAboveTOT );
59
60 return sc;
61}
62
64 return StatusCode::SUCCESS;
65}
66
68// other methods inherited from ITRT_TrackSegmentsMaker: newEvent, newRegion, endEvent
70
71std::unique_ptr<InDet::ITRT_TrackSegmentsMaker::IEventData> InDet::TRT_TrackSegmentsMaker_BarrelCosmics::newEvent(const EventContext& ctx) const {
72
73 ATH_MSG_DEBUG( "InDet::TRT_TrackSegmentsMaker_BarrelCosmics::newEvent()");
74
75
76 SG::ReadHandle<InDet::TRT_DriftCircleContainer> trtcontainer(m_driftCirclesName, ctx); // get TRT_DriftCircle list from StoreGate containers
77
78 if (not trtcontainer.isValid()) {
79 msg(MSG::ERROR) << "Could not find TRT_DriftCircles collection!" << endmsg;
80 std::stringstream msg;
81 msg << name() <<" no TRT drift circles : " << m_driftCirclesName.key();
82 throw std::runtime_error(msg.str());
83 }
84
85 std::unique_ptr<TRT_TrackSegmentsMaker_BarrelCosmics::EventData>
86 event_data = std::make_unique<TRT_TrackSegmentsMaker_BarrelCosmics::EventData>(trtcontainer.cptr());
87
88 for(const auto *it : *trtcontainer) {
89
90 const InDet::TRT_DriftCircleCollection *colNext=&(*it);
91 if (!colNext) { msg(MSG::WARNING) << "newEvent(): !colNext " << endmsg; continue; }
92
93 for (const auto *circleit : *colNext){
94
95 if ( m_trtid->is_barrel((*circleit).identify()) ) { // TRT barrel
96
97 event_data->m_listHits.push_back( circleit );
98
99 double onMyTime = circleit->timeOverThreshold();
100 if ( circleit->firstBinHigh() && (circleit->trailingEdge() != 24) && (circleit->trailingEdge() != 0) ) onMyTime = ( (double) circleit->trailingEdge() + 1 ) * 3.125;
101 if (onMyTime < m_minSeedTOT) continue;
102
103 const Amg::Vector3D &center = circleit->detectorElement()->surface(circleit->identify()).center();
104 event_data->m_listHitCenter.push_back( center );
105 }
106 }
107 } // end trtcontainer loop
108
109 if ( m_maxTotalHits && ((int)event_data->m_listHits.size()) > m_maxTotalHits ) {
110 ATH_MSG_DEBUG( "skipping high occupancy event of " << event_data->m_listHits.size() << " barrel hits, limit at "
111 << m_maxTotalHits );
112 event_data->clear();
113 }
114
115 ATH_MSG_DEBUG( "newEvent(): Number of TRT barrel hits: " << event_data->m_listHits.size()
116 << " Number of hits with TOT > " << m_minSeedTOT << ": " << event_data->m_listHitCenter.size() );
117
118 return std::unique_ptr<InDet::ITRT_TrackSegmentsMaker::IEventData>(event_data.release());
119}
120
121std::unique_ptr<InDet::ITRT_TrackSegmentsMaker::IEventData>
122InDet::TRT_TrackSegmentsMaker_BarrelCosmics::newRegion(const EventContext& ctx, const std::vector<IdentifierHash> &vTRT) const {
123
124 ATH_MSG_DEBUG("InDet::TRT_TrackSegmentsMaker_BarrelCosmics::newRegion()" );
125
127 if (not trtcontainer.isValid()) {
128 msg(MSG::ERROR) << "m_trtcontainer is empty!!!" << endmsg;
129 std::stringstream msg;
130 msg << name() <<" no TRT drift circles : " << m_driftCirclesName.key();
131 throw std::runtime_error(msg.str());
132 }
133
134 std::unique_ptr<TRT_TrackSegmentsMaker_BarrelCosmics::EventData>
135 event_data = std::make_unique<TRT_TrackSegmentsMaker_BarrelCosmics::EventData>(trtcontainer.cptr());
136
137
138 for(auto d : vTRT) {
139 for ( InDet::TRT_DriftCircleContainer::const_iterator w = trtcontainer->indexFind(d); w!=trtcontainer->end(); ++w ) {
140 for(const auto *circleit : **w) {
141
142 if(std::abs(m_trtid->barrel_ec( circleit->identify() ))!=1) continue;
143
144 event_data->m_listHits.push_back(circleit);
145
146 double onMyTime = circleit->timeOverThreshold();
147 if ( circleit->firstBinHigh() && (circleit->trailingEdge() != 24) && (circleit->trailingEdge() != 0) ) onMyTime = ( (double) circleit->trailingEdge() + 1 ) * 3.125;
148 if (onMyTime < m_minSeedTOT) continue;
149 const Amg::Vector3D &center = circleit->detectorElement()->surface(circleit->identify()).center();
150 event_data->m_listHitCenter.push_back( center );
151 }
152 }
153 }
154
155 if ( m_maxTotalHits && ((int)event_data->m_listHits.size()) > m_maxTotalHits ) {
156 ATH_MSG_DEBUG("skipping high occupancy event of " << event_data->m_listHits.size() << " barrel hits, limit at "
157 << m_maxTotalHits);
158 event_data->clear();
159 }
160
161 ATH_MSG_DEBUG( "newRegion(): Number of TRT barrel hits: " << event_data->m_listHits.size()
162 << " Number of hits with TOT > " << m_minSeedTOT << ": " << event_data->m_listHitCenter.size() );
163
164 return std::unique_ptr<InDet::ITRT_TrackSegmentsMaker::IEventData>(event_data.release());
165}
166
170
171 ATH_MSG_DEBUG("InDet::TRT_TrackSegmentsMaker_BarrelCosmics::endEvent()" );
172
173 // elements of m_segments created by new have not been passed on
174 if ( event_data.m_segmentDriftCirclesCount < event_data.m_segments.size() ) {
175 msg(MSG::WARNING) << "endEvent() you called the function to create the segments but didn't retrieve them later??" << endmsg;
176 msg(MSG::WARNING) << "endEvent() deleting remaining elements of m_segments" << endmsg;
177 for (unsigned int i=event_data.m_segmentDriftCirclesCount; i<event_data.m_segments.size(); i++) delete event_data.m_segments[i];
178 }
179
180 }
181
183// segment finding algorithm function: find, calls findSeed
185
191
192 if (!m_magneticField) { findOld(event_data); return; }
193
194 ATH_MSG_DEBUG("InDet::TRT_TrackSegmentsMaker_BarrelCosmics::find()" );
195
196
197 if ((int)event_data.m_listHitCenter.size()<m_minHitsAboveTOT) return;
198
199 if (!event_data.m_segmentDriftCircles.empty()) { // debug only
200 msg(MSG::WARNING) << "TRT_TrackSegmentsMaker_BarrelCosmics::find() probably called twice per event? or newEvent / newRegion have not been called. check your program" << endmsg;
201 event_data.clear(); return;
202 }
203
204 std::vector<double> x0, phi, pivotX, pivotY, Xparabola, cotanParabola, InverseR;
205 std::vector<int> nHitsPosY, nHitsNegY;
206
207 // first find seeds
208 int nIterations(15), countCalls(0);
209 while (countCalls<150) {
210
211 double xRange(1000.), phiRange(M_PI_4);
212 double par[] = {0., M_PI_2, 0., 0., 0., 0., 0., 0.};
213 int foundSeed(0);
214
215 for (int i=0; i<nIterations; i++) {
216
217 countCalls++;
218 foundSeed = findSeed( par[0]-xRange, par[0]+xRange, par[1]-phiRange, par[1]+phiRange, par, event_data);
219 if ( !foundSeed ) break;
220
221 xRange /= 3.;
222 phiRange /= 2.;
223 if ( xRange < 0.01 || phiRange < 0.00001) break;
224 }
225 if (!foundSeed) break;
226
227 if (m_magneticField) findSeedInverseR(par, event_data);
228
229// remove the hits associated with this region, save this region, search again
230 int countAssociatedHits[] = {0, 0};
231 double cosphi = std::cos( par[1] );
232 double sinphi = std::sqrt( 1. - cosphi*cosphi );
233 double range = (m_magneticField) ? 10. : 2.; // BE CAREFULL ABOUT THIS RANGE, IT USED TO BE 2 BUT PERHAPS IT SHOULD BE INCREASED TO 10 in the case of magnetic field!
234
235
236 double measx[200], measy[200];
237 int countMeas(0);
238
239 for (std::vector< Amg::Vector3D>::iterator it = event_data.m_listHitCenter.begin(); it!=event_data.m_listHitCenter.end(); ) {
240
241 double a = (par[3]-it->x())*sinphi+(it->y()-par[4])*cosphi;
242 double b = (m_magneticField) ? 0.5 * (std::pow(it->x()-par[3], 2.) + std::pow(it->y()-par[4], 2.) - std::pow(a, 2)) : 0.;
243 if ( std::abs(a+par[2]*b) > range ) { ++it; continue; }
244
245 if (m_magneticField && countMeas<200) { measx[countMeas] = it->x(); measy[countMeas] = it->y(); countMeas++; }
246
247 countAssociatedHits[(it->y()>0)?0:1]++;
248 it = event_data.m_listHitCenter.erase( it );
249 }
250
251 if ( countAssociatedHits[0]+countAssociatedHits[1] < m_minHitsAboveTOT ) continue;
252
253 if (m_magneticField) segFit(measx, measy, countMeas, nullptr, par+3);
254
255 ATH_MSG_DEBUG("countAssociatedHits " << countAssociatedHits[0] << " "
256 << countAssociatedHits[1] << " m_minHitsAboveTOT " << m_minHitsAboveTOT );
257 x0.push_back( par[0] ); phi.push_back( par[1] ); nHitsPosY.push_back( countAssociatedHits[0] ); nHitsNegY.push_back( countAssociatedHits[1] );
258
259 pivotX.push_back( par[3] ); pivotY.push_back( par[4] ); Xparabola.push_back( par[5] ); cotanParabola.push_back( par[6] ); InverseR.push_back( par[7] );
260 } // end for (int i=0; i<nIterations; i++) loop
261// end finding seeds
262
263// save all TRT hits for found segments
264
265 int nFoundSegments = x0.size(); if (nFoundSegments>10) nFoundSegments = 10; // number of found segments
266 {
267 for (int i=0; i<nFoundSegments; i++) {
268 event_data.m_segmentDriftCircles.emplace_back();
269 }
270 }
271
272 double window = (m_magneticField) ? 4. : 2.; // change for magnetic field
273 double residual;
274
275 for (unsigned int i=0; i<event_data.m_listHits.size(); i++) {
276
277 const Amg::Vector3D &center = event_data.m_listHits[i]->detectorElement()->surface(event_data.m_listHits[i]->identify()).center();
278 for (int j=0; j<nFoundSegments; j++) {
279
280 if (nHitsPosY[j]<5 && center.y()>0.) continue;
281 if (nHitsNegY[j]<5 && center.y()<0.) continue;
282
283 if (m_magneticField) {
284 double sinphi = std::sqrt(1./(1+cotanParabola[j]*cotanParabola[j]));
285 double cosphi = std::sqrt(1.-sinphi*sinphi); if (cotanParabola[j]<0.) cosphi *= -1.;
286 double a = (pivotX[j]+Xparabola[j]-center.x())*sinphi+(center.y()-pivotY[j])*cosphi;
287 double b = 0.5 * (std::pow(center.x()-pivotX[j], 2.) + std::pow(center.y()-pivotY[j], 2.) - std::pow(a, 2)) ;
288 residual = a + InverseR[j] * b;
289
290 } else {
291 double cosphi = std::cos(phi[j]);
292 double sinphi = std::sqrt(1.-cosphi*cosphi);
293 residual = (x0[j]-center.x())*sinphi+center.y()*cosphi;
294 }
295
296 if (std::abs(residual)<window) {
297 event_data.m_segmentDriftCircles[j].push_back(event_data.m_listHits[i]);
298 break;
299 }
300 }
301 }
302 // end saving TRT hits
303
304
305
306 // merge segments: can be simple: if 2 of them similar values of impact par -> copy hits from the second vector into the first one and clear the second one
307 // XXXXXXXX to do: determine the optimal merging parameters
308
309 if (m_mergeSegments) { // merge segments, not yet tested properly
310
311 if (x0.size()>1) {
312 int mergeI = 0;
313 int mergeJ = 0;
314 double mergeX0(100.), mergePhi(0.1);
315 for (int i=0; i<nFoundSegments; i++) {
316 for (int j=i+1; j<nFoundSegments; j++) {
317 if ( std::abs(x0[i]-x0[j])<mergeX0 && std::abs(phi[i]-phi[j])<mergePhi ) {
318 mergeI = i;
319 mergeJ = j;
320 mergeX0 = std::abs(x0[i]-x0[j]);
321 mergePhi = std::abs(phi[i]-phi[j]);
322 }
323 }
324 }
325 if (mergeI != mergeJ) {
326 ATH_MSG_DEBUG("Merge segments " << mergeI << " and " << mergeJ << " of size "
327 << event_data.m_segmentDriftCircles[mergeI].size() << ", " << event_data.m_segmentDriftCircles[mergeJ].size()
328 << "; difference in the impact par x: " << mergeX0 << " phi: " << mergePhi );
329 for (unsigned int i=0; i<event_data.m_segmentDriftCircles[mergeJ].size(); i++) {
330 event_data.m_segmentDriftCircles[mergeI].push_back(event_data.m_segmentDriftCircles[mergeJ][i]);
331 }
332 event_data.m_segmentDriftCircles[mergeJ].clear();
333 }
334 }
335 } // end merge segments
336
337
338 if (msgLvl(MSG::DEBUG)) { // debug: check how many hits per found segments
339 msg(MSG::DEBUG) << "find() debug (" << nFoundSegments << ")" ;
340 for (unsigned int i=0; i<event_data.m_segmentDriftCircles.size(); i++) {
341 msg(MSG::DEBUG) << " " << i << " " << event_data.m_segmentDriftCircles[i].size();
342 }
343 msg(MSG::DEBUG) << endmsg;
344 }
345
346 for (unsigned int i=0; i<event_data.m_segmentDriftCircles.size(); i++) { // convert to segments
347 if ((int)event_data.m_segmentDriftCircles[i].size()<m_minHitsForSegment) continue;
348 double trackpar[] = {x0[i], phi[i]};
349 convert(event_data.m_segmentDriftCircles[i], trackpar, event_data);
350 }
351 ATH_MSG_DEBUG( "find(), number of converted segments: " << event_data.m_segments.size() );
352}
353
355// the last method inherited from ITRT_TrackSegmentsMaker: next()
356// returns ptr to i-th found TrackSegment, until it returns 0 (no more segments)
357// assumes any provided segments will be deleted by those who request them!!!
359
361{
364
365 // next 6 lines: for debugging purposes only
366 ATH_MSG_DEBUG( "InDet::TRT_TrackSegmentsMaker_BarrelCosmics::next(): return "
367 << event_data.m_segmentDriftCirclesCount << " out of " << event_data.m_segments.size() );
368
369 if (event_data.m_segmentDriftCirclesCount > event_data.m_segments.size())
370 msg(MSG::ERROR) << "m_segmentDriftCirclesCount = " << event_data.m_segmentDriftCirclesCount << ", m_segments.size() = "
371 << event_data.m_segments.size() << endmsg;
372
373 return (event_data.m_segmentDriftCirclesCount<event_data.m_segments.size())
374 ?(event_data.m_segments[event_data.m_segmentDriftCirclesCount++])
375 :nullptr;
376}
377
378
379
381 double xmax,
382 double phimin,
383 double phimax,
384 double *bestParameters,
386{
387 if ((int)event_data.m_listHitCenter.size()<m_minHitsAboveTOT) return 0;
388
389 double xbin = (xmax-xmin) / (1.*m_nBinsInX);
390 double searchWindow = (xbin>2.) ? xbin : 2.; //searchWindow *= searchWindow;
391 double phibin = (phimax-phimin) / (1.*m_nBinsInPhi);
392 int maxHits(0), index, indexmin, indexmax;
393
394 int binInX[100]; for (int & i : binInX) i = 0;
395
396 for (int j=0; j<m_nBinsInPhi; j++) {
397
398 double phi = phimin+(0.5+j)*phibin;
399 double cosphi( std::cos(phi) ), sinphi( std::sin(phi) );
400 for (int i=0; i<m_nBinsInX; i++) binInX[i] = 0;
401 double transformXMin(xmin*sinphi); xbin = (xmax-xmin)*sinphi / (1.*m_nBinsInX);
402
403 for (auto & k : event_data.m_listHitCenter) {
404 double transformX = k.x() * sinphi - k.y() * cosphi;
405 indexmin = (int) ((transformX-transformXMin-searchWindow)/xbin);
406 indexmax = (int) ((transformX-transformXMin+searchWindow)/xbin);
407
408 if (indexmin<0) indexmin=0;
409 if (indexmax>99) indexmax=99;
410 for (index = indexmin; index<=indexmax; index++) binInX[index]++;
411 }
412 index = -1;
413 for (int i=0; i<m_nBinsInX; i++) if (binInX[i]>maxHits) { index = i; maxHits = binInX[i]; }
414 if (index>=0) {
415 bestParameters[0] = ( xbin * index + transformXMin) / sinphi;
416 bestParameters[1] = phi;
417 }
418 }
419
420 return (maxHits>=m_minHitsForSeed) ? 1 : 0;
421}
422
425
426 double cosphi = std::cos( par[1] );
427 double sinphi = std::sqrt( 1. - cosphi*cosphi );
428
429 double parTransformX = par[0] * sinphi;
430
431 double meanTransformY(0.);
432 int countMeanTransformY(0);
433 for (auto & k : event_data.m_listHitCenter) { // first determine the mean
434
435 double transformX = k.x() * sinphi - k.y() * cosphi;
436 if ( std::abs(transformX-parTransformX) > 2. ) continue;
437 meanTransformY += k.x() * cosphi + k.y() * sinphi;
438 countMeanTransformY++;
439 }
440 if (!countMeanTransformY) {
441 ATH_MSG_WARNING("InDet::TRT_TrackSegmentsMaker_BarrelCosmics::findSeedInverseR(), no hits in the seed region???" );
442 return;
443 }
444 meanTransformY /= (double) countMeanTransformY;
445
446
447// hard-coded parameters in the second part of this function:
448// search range: [-0.00202, 0.00202] divided into 101 bins of width 0.00004 (such that range [-0.00002, 0.00002] is bin 50)
449// 1/R = 0.002/mm corresponds to R = 0.5 m
450// circle search window: 200 mm (for R=1m, circle-to-line distance 0.5 m away from a common point is 130mm)
451
452// -> pivot at: x_pivot = meanTransformY * cosphi + parTransformX * sinphi; y_pivot = meanTransformY * sinphi - parTransformX * cosphi;
453 double mean[] = { meanTransformY * cosphi + parTransformX * sinphi, meanTransformY * sinphi - parTransformX * cosphi };
454 int scanInverseR[101]; for (int & i : scanInverseR) i = 0;
455
456 double window = 10.; // parabola search window BE CAREFULL OF THIS RANGE!
457
458 for (auto & it : event_data.m_listHitCenter) {
459
460 double transformX = it.x() * sinphi - it.y() * cosphi;
461 if ( std::abs(transformX-parTransformX) > 200. ) continue; // search circles in a broad band
462
463 double a = (mean[0]-it.x())*sinphi+(it.y()-mean[1])*cosphi; // TEST THAT THIS CONDITION ( if ( fabs(a) > 200. ) continue; ) IS THE SAME AS ABOVE!
464 double b = 0.5 * (std::pow(it.x()-mean[0], 2.) + std::pow(it.y()-mean[1], 2.) - std::pow(a, 2));
465
466 // estimated allowed range in radius: x_1,2 = ( +- window - a ) / b -> get integer as ((int) (ceil(x/bin-0.5))) + 50
467 b *= 0.00004; // multiply by the bin width
468
469 int i1(0), i2(0);
470 if (std::abs(b)>0.0000001) {
471 double i1_tmp = ( window - a ) / b + 0.5;
472 double i2_tmp = ( -window - a ) / b + 0.5;
473 if (std::abs(i1_tmp)<1000.) { i1 = (int) (std::ceil( i1_tmp )); i1 += 50; }
474 if (std::abs(i2_tmp)<1000.) { i2 = (int) (std::ceil( i2_tmp )); i2 += 50; }
475 }
476 if (i1>100) { i1 = 100; } else { if (i1<0) i1 = 0; }
477 if (i2>100) { i2 = 100; } else { if (i2<0) i2 = 0; }
478 if (i1+i2==0 || i1+i2==200) continue; // out of search range
479 if (i1<=i2) { for (int i=i1; i<=i2; i++) scanInverseR[i]++; }
480 else { for (int i=i2; i<=i1; i++) scanInverseR[i]++; }
481 }
482
483 int iMax(-1), hitsMax(0);
484 for (int i=0; i<101; i++) if (scanInverseR[i]>hitsMax) { iMax = i; hitsMax = scanInverseR[i]; } // find max bin
485
486 par[2] = 0.00004 * (iMax-50); // save estimated inverse radius
487 for (int i=0; i<2; i++) par[3+i] = mean[i]; }
488
490// convert list of TRT hits to "TrackSegment": convert()
491//
492// the only TrackSegment constructor one can use for the purpose:
493//
494// TrackSegment( const LocalParameters* locpars, // LocalParameters(radius, z0, phi0, theta, QoverP) at Surface point "sf"
495// const ErrorMatrix* locerr, // error matrix for above params
496// const Surface* sf, //
497// DataVector<const MeasurementBase>* crots, // TRT_DriftCircleOnTrack : RIO_OnTrack : MeasurementBase
498// FitQuality* fqual, // FitQuality(float chi2, short ndf)
499// Segment::Author author = Segment::AuthorUnknown
500// );
502
503void InDet::TRT_TrackSegmentsMaker_BarrelCosmics::linearRegressionParabolaFit(double *mean, double *a) { // input: elsewhere calculated mean, output: result a
504
505 double A[6], discriminant;
506
507 A[0] = mean[8]*mean[8]-mean[4]*mean[9]; // parabola fit
508 A[1] = mean[1]*mean[9]-mean[4]*mean[8];
509 A[2] = mean[4]*mean[4]-mean[1]*mean[8];
510 A[3] = mean[4]*mean[4]-mean[9];
511 A[4] = mean[8]-mean[1]*mean[4];
512 A[5] = mean[1]*mean[1]-mean[4];
513 discriminant = std::abs( std::pow(mean[4], 3) + mean[8]*mean[8] + mean[1]*mean[1]*mean[9] - mean[4]*(2.*mean[1]*mean[8]+mean[9]) );
514
515 a[0] = A[0] * mean[0] + A[1] * mean[2] + A[2] * mean[6];
516 a[1] = A[1] * mean[0] + A[3] * mean[2] + A[4] * mean[6];
517 a[2] = A[2] * mean[0] + A[4] * mean[2] + A[5] * mean[6];
518 discriminant *= -1.;
519 for (int i=0; i<3; i++) { a[i] /= discriminant; }
520
521 double inverseSin3phi = 1. + a[1]*a[1]; inverseSin3phi *= std::sqrt(inverseSin3phi); // 1/sin^3phi
522 a[2] *= 2. / inverseSin3phi; }
523
525
526 double y1 = (dc1->detectorElement())->surface(dc1->identify()).center().y();
527 double y2 = (dc2->detectorElement())->surface(dc2->identify()).center().y();
528 return ( y1 > y2 );
529}
530
531void InDet::TRT_TrackSegmentsMaker_BarrelCosmics::convert(std::vector<const InDet::TRT_DriftCircle *> &hits,
532 double *trackpar,
534//Track Segment production - based on TRT_TrackSegmentsMaker_ECcosmics Tool
535
536 ATH_MSG_DEBUG("InDet::TRT_TrackSegmentsMaker_BarrelCosmics::convert() segment " << event_data.m_segments.size()
537 << ", N hits = " << hits.size());
538
539 if (hits.size()<5) { msg(MSG::ERROR) << "convert(): empty list of hits! size: " << hits.size() << endmsg; return; }
540
541 // sort the vector of hits
543
544 // remove hits from the opposite side
545 int countOppositeSide(0);
546 for (std::vector<const InDet::TRT_DriftCircle *>::iterator it = hits.begin(); it != hits.end(); ) {
547 int side = m_trtid->barrel_ec( (*it)->identify() );
548 int previous = ( it == hits.begin() ) ? m_trtid->barrel_ec( (*(it+2))->identify() ) : m_trtid->barrel_ec( (*(it-1))->identify() );
549 int next = ( it == hits.end()-1 ) ? m_trtid->barrel_ec( (*(it-2))->identify() ) : m_trtid->barrel_ec( (*(it+1))->identify() );
550 if ( previous==next && side*next==-1 ) { it = hits.erase( it ); countOppositeSide++; }
551 else { ++it; }
552 }
553 if (countOppositeSide>5) {
554 ATH_MSG_DEBUG( "convert(): removed " << countOppositeSide << " hits from the other side, N remaining hits: " << hits.size() );
555 }
556 if (hits.size()<5) {
557 ATH_MSG_WARNING( "convert(): not enough hits after opposite side removal: " << hits.size() << ", removed: "
558 << countOppositeSide );
559 return;
560 }
561
562 // make the linear regression fit
563 double mean[10]; for (double & i : mean) i = 0.;
564 for (auto & hit : hits) {
565
566 mean[0] += (hit->detectorElement())->surface(hit->identify()).center().x();
567 mean[1] += (hit->detectorElement())->surface(hit->identify()).center().y();
568 }
569 for (int i=0; i<2; i++) mean[i] /= (double) hits.size();
570 int iPivot(-1); double yPivot(10000000.); // choose pivot closest to the mean
571 for (unsigned int i=0; i<hits.size(); i++) {
572 double x = (hits[i]->detectorElement())->surface(hits[i]->identify()).center().x() - mean[0];
573 double y = (hits[i]->detectorElement())->surface(hits[i]->identify()).center().y() - mean[1];
574 double x2 = x * x;
575 double y2 = y * y;
576 mean[2] += x * y;
577 mean[3] += x2;
578 mean[4] += y2;
579 mean[5] += x2 * y;
580 mean[6] += y2 * x;
581 mean[7] += x2 * x;
582 mean[8] += y2 * y;
583 mean[9] += y2 * y2;
584 if (x2+y2<yPivot) { yPivot = x2+y2; iPivot = i; }
585 }
586
587 if (iPivot==-1) {
588 msg(MSG::ERROR) << "SF pivot index not found!!! " << yPivot << " " << iPivot << endmsg;
589 iPivot = 0;
590 }
591
592 double cotanPhi = mean[2] / mean[4];
593 double phi = std::atan(1./cotanPhi); if (phi<0.) phi += M_PI;
594 ATH_MSG_DEBUG("compare parameters X : " << trackpar[0] << " vs. " << mean[0]-mean[1]*cotanPhi );
595 ATH_MSG_DEBUG("compare parameters phi: " << trackpar[1] << " vs. " << phi );
596
597 double qOverp = 0.; // units q / MeV, set only if there is magnetic field
598
599 if (m_magneticField) { // fit parabola instead of line, add circle if neccessary
600
601 double A[6], discriminant, a[3];
602
603 A[0] = mean[8]*mean[8]-mean[4]*mean[9]; // parabola fit
604 A[1] = mean[1]*mean[9]-mean[4]*mean[8];
605 A[2] = mean[4]*mean[4]-mean[1]*mean[8];
606 A[3] = mean[4]*mean[4]-mean[9];
607 A[4] = mean[8]-mean[1]*mean[4];
608 A[5] = mean[1]*mean[1]-mean[4];
609 discriminant = std::abs( std::pow(mean[4], 3) + mean[8]*mean[8] + mean[1]*mean[1]*mean[9] - mean[4]*(2.*mean[1]*mean[8]+mean[9]) );
610
611 a[0] = A[0] * mean[0] + A[1] * mean[2] + A[2] * mean[6];
612 a[1] = A[1] * mean[0] + A[3] * mean[2] + A[4] * mean[6];
613 a[2] = A[2] * mean[0] + A[4] * mean[2] + A[5] * mean[6];
614 discriminant *= -1.;
615 for (double & i : a) { i /= discriminant; }
616
617 double Xparabola = a[0];
618 double cotanParabola = a[1];
619 double inverseSin3phi = 1. + a[1]*a[1]; inverseSin3phi *= std::sqrt(inverseSin3phi); // 1/sin^3phi
620 double inverseR = 2. * a[2] / inverseSin3phi;
621
622 if (msgLvl(MSG::DEBUG)) {
623 msg(MSG::DEBUG) << "TRT_TrackSegmentsMaker_BarrelCosmics: parabola fit, X: " << Xparabola << endmsg;
624 msg(MSG::DEBUG) << " parabola fit, cotan: " << cotanParabola << " compare to " << cotanPhi << endmsg;
625 msg(MSG::DEBUG) << " parabola fit, 1/R: " << inverseR << "/ mm " << endmsg;
626 }
627
628 qOverp = inverseR / 0.6; // [1/MeV]; 1 / (eBRC) = 1 / (e 2T 3 100 M m/s mm) / R[mm] = (1/R[mm]) / 0.6
629 }
630
631// errors from the fit:
632// - for d: sigma_r / sqrt(N hits) -> sigma_d^2 =
633// - for cotanPhi: sigma_r / (sqrt(N hits * <y^2>) * sinPhi)
634// - for phi: 1 / (1 + cotanPhi^2) * DcotanPhi // note that this mthod gives a biased estimate for phi!
635
636 Amg::MatrixX cov(5,5);
637
638 double f33 = 4. / ( 3. * (1.+cotanPhi*cotanPhi) * mean[4] ); // see above, mean[4] = N hits * <y^2>
639 cov<<
640 4. , 0., 0., 0., 0.,
641 0. ,45633., 0., 0., 0.,
642 0. , 0., f33, 0., 0.,
643 0. , 0., 0., 0.2, 0.,
644 0. , 0., 0., 0., 1.;
645
646 Amg::Vector3D hepVec( mean[0], mean[1],
647 (hits[iPivot]->detectorElement())->surface(hits[iPivot]->identify()).center().z() );
648
649
650 Amg::Transform3D htrans = Amg::Transform3D(Amg::Translation3D(hepVec)*Amg::RotationMatrix3D::Identity());
651
652
653
655
656 if (phi>0.) phi -= M_PI;
657 if (phi<-M_PI || phi>0.) msg(MSG::ERROR) << "phi value problem: " << phi << endmsg;
658
659 Trk::LocalParameters par(0., 0., phi, M_PI_2, qOverp);
660
661 ATH_MSG_DEBUG("pivot: " << sur->center() << ", cross-check: " << mean[0] << " " << mean[1] );
662
663 // calculate TrackParameters so that one can calculate input for TRT_DriftCircleOnTrack
664 // Global direction of the track parameters
665 Amg::Vector3D dir(std::cos(phi), std::sin(phi), 0.);
666
668
669 for (const auto *DC : hits) { // rewrite: InDet::TRT_DriftCircle -> InDet::TRT_DriftCircleOnTrack
670 // based on http://atlas-sw.cern.ch/cgi-bin/viewcvs-atlas.cgi/offline/InnerDetector/InDetRecTools/TRT_DriftCircleOnTrackTool/src/TRT_DriftCircleOnTrackNoDriftTimeTool.cxx?revision=1.11&view=markup
671
672 // Straw identification
673 const InDetDD::TRT_BaseElement* pE = DC->detectorElement();
674 if (!pE) {
675 ATH_MSG_ERROR("convert(): no detectorElement info!");
676 continue;
677 }
678 Identifier iD = DC->identify();
679 IdentifierHash iH = m_trtid->straw_layer_hash(m_trtid->layer_id(iD));
680
681 // TRT_DriftCircleOnTrack production
683 Trk::LocalParameters lp(radius);
684 Amg::MatrixX cov(1,1); cov<<1.33333;
685 const InDet::TRT_DriftCircleOnTrack *element = new InDet::TRT_DriftCircleOnTrack(DC, std::move(lp), std::move(cov),
686 iH, 1., dir, Trk::NODRIFTTIME);
687 rio.push_back( dynamic_cast<const Trk::MeasurementBase *>(element) );
688 } // end fill rio
689
690 // add pseudo measurement - follow up https://savannah.cern.ch/bugs/?41079
691if (1) { // limit the scope of all these variables
692
693 Trk::DefinedParameter pseudoz( 0., Trk::locZ);
694 Trk::DefinedParameter pseudotheta( M_PI/2., Trk::theta);
695 Trk::LocalParameters pseudoPar( pseudoz, pseudotheta );
696
697 Amg::MatrixX pseudocov(2,2);
698 pseudocov <<
699 1. , 0.,
700 0. , .001;
701
702 const Trk::Surface &surf = (hits[hits.size()-1]->detectorElement())->surface(hits[hits.size()-1]->identify());
703
704 Amg::RotationMatrix3D linerot=surf.transform().rotation();
705 Amg::Vector3D firstcenter=surf.center();
706 Amg::Vector3D faketransl(firstcenter.x(),firstcenter.y()-5.,firstcenter.z());
707 Amg::Transform3D faketransf = Amg::Transform3D(linerot* Amg::Translation3D(faketransl));
708 Trk::StraightLineSurface *pseudoSurface = new Trk::StraightLineSurface( faketransf, 99999., 99999. );
709
711 std::move(pseudoPar),
712 std::move(pseudocov),
713 *pseudoSurface );
714 rio.push_back( pseudo );
715 delete pseudoSurface;
716}
717
718 double chi2 = mean[3] - 2.*cotanPhi*mean[2] + mean[4]*cotanPhi*cotanPhi;
719 chi2 /= ( 1. + cotanPhi*cotanPhi );
720 int ndf = (int) hits.size() - 2;
721 ATH_MSG_DEBUG( "Chi2 = " << chi2 << ", ndf = " << ndf << ", chi2/ndf = " << chi2/ndf );
722 Trk::FitQuality * fqu = new Trk::FitQuality(chi2, ndf);
723
725 std::move(par), std::move(cov), sur, std::move(rio), fqu, Trk::Segment::TRT_SegmentMaker);
726
727 //add segment to list of segments
728 ATH_MSG_DEBUG( "Add " << event_data.m_segments.size() << "th segment to list" );
729 event_data.m_segments.push_back(segment);
730}
731
733//
734// DEBUG FUNCTIONS
735//
737
738void InDet::TRT_TrackSegmentsMaker_BarrelCosmics::segFit(double *measx, double *measy, int nhits, double *residuals, double *result ) {
739
740 if (nhits<10) return;
741
742 double shift[] = {0., 0.}; // shift the points so that mean x = 0 and mean y = 0
743 for (int i=0; i<nhits; i++) { shift[0] += measx[i]; shift[1] += measy[i]; }
744 for (double & i : shift) i /= (double) nhits;
745 for (int i=0; i<nhits; i++) { measx[i] -= shift[0]; measy[i] -= shift[1]; }
746
747 double mean[10]; for (double & i : mean) i = 0.; // calculate momenta needed for the line, circle and parabola fit
748 for (int i=0; i<nhits; i++) {
749 double x2 = measx[i]*measx[i];
750 double y2 = measy[i]*measy[i];
751 mean[0] += measx[i];
752 mean[1] += measy[i];
753 mean[2] += measx[i]*measy[i];
754 mean[3] += x2;
755 mean[4] += y2;
756 mean[5] += x2 * measy[i];
757 mean[6] += y2 * measx[i];
758 mean[7] += x2 * measx[i];
759 mean[8] += y2 * measy[i];
760 mean[9] += y2 * y2;
761 }
762 for (double & i : mean) i /= (double) nhits;
763
764 double A[6], discriminant, a[3];
765
766 A[0] = mean[8]*mean[8]-mean[4]*mean[9]; // parabola fit, no need for mean[3], [5], [7], se mi zdi
767 A[1] = mean[1]*mean[9]-mean[4]*mean[8];
768 A[2] = mean[4]*mean[4]-mean[1]*mean[8];
769 A[3] = mean[4]*mean[4]-mean[9];
770 A[4] = mean[8]-mean[1]*mean[4];
771 A[5] = mean[1]*mean[1]-mean[4];
772 discriminant = std::abs( std::pow(mean[4], 3) + mean[8]*mean[8] + mean[1]*mean[1]*mean[9] - mean[4]*(2.*mean[1]*mean[8]+mean[9]) );
773
774 a[0] = A[0] * mean[0] + A[1] * mean[2] + A[2] * mean[6];
775 a[1] = A[1] * mean[0] + A[3] * mean[2] + A[4] * mean[6];
776 a[2] = A[2] * mean[0] + A[4] * mean[2] + A[5] * mean[6];
777 discriminant *= -1.;
778 for (double & i : a) i /= discriminant;
779 double Xparabola = a[0];
780 double cotanParabola = a[1];
781 double inverseSin3phi = 1. + a[1]*a[1]; inverseSin3phi *= std::sqrt(inverseSin3phi); // 1/sin^3phi
782 double inverseR = 2. * a[2] / inverseSin3phi;
783 double sinphi = std::sqrt(1./(1+a[1]*a[1]));
784 double cosphi = std::sqrt(1.-sinphi*sinphi); if (cotanParabola<0.) cosphi *= -1.;
785
786 if (result) {
787 for (int i=0; i<2; i++) result[i] = shift[i]; // pivot
788 result[2] = Xparabola;
789 result[3] = cotanParabola;
790 result[4] = inverseR;
791 }
792
793 if (residuals) {
794 for (int i=0; i<nhits; i++) {
795 double tmp = (Xparabola-measx[i])*sinphi+measy[i]*cosphi;
796 double res = (tmp+0.5*inverseR*(std::pow(measx[i]-Xparabola, 2.) + std::pow(measy[i], 2.) - std::pow(tmp, 2)));
797 residuals[i] = res;
798 }
799 }
800
801 for (int i=0; i<nhits; i++) { measx[i] += shift[0]; measy[i] += shift[1]; }
802
803 }
804
805
807//
808// DEBUG cpu consumption: old find method
809//
811
812
814
815 ATH_MSG_DEBUG( "InDet::TRT_TrackSegmentsMaker_BarrelCosmics::find()");
816
817 if ((int)event_data.m_listHitCenter.size()<m_minHitsAboveTOT) return;
818
819 if (!event_data.m_segmentDriftCircles.empty()) { // debug only
820 ATH_MSG_WARNING("find probably called twice per event? or newEvent / newRegion have not been called. check your program" );
821 event_data.clear(); return;
822 }
823
824 std::vector<double> x0, phi;
825 std::vector<int> nHitsPosY, nHitsNegY;
826
827 // first find seeds
828 int nIterations(15), countCalls(0);
829 while (countCalls<150) {
830
831 double xRange(1000.), phiRange(M_PI_4);
832 double par[] = {0., M_PI_2};
833 int foundSeed(0);
834
835 for (int i=0; i<nIterations; i++) {
836
837 countCalls++;
838 foundSeed = findSeed( par[0]-xRange, par[0]+xRange, par[1]-phiRange, par[1]+phiRange, par, event_data);
839 if ( !foundSeed ) break;
840
841 xRange /= 3.;
842 phiRange /= 2.;
843 if ( xRange < 0.01 || phiRange < 0.00001) break;
844
845 }
846 if (!foundSeed) break;
847
848// remove the hits associated with this region, save this region, search again
849 int countAssociatedHits[] = {0, 0};
850 double cosphi = std::cos(par[1]);
851 double sinphi = std::sqrt(1.-cosphi*cosphi);
852
853 for (std::vector< Amg::Vector3D >::iterator it = event_data.m_listHitCenter.begin(); it!=event_data.m_listHitCenter.end(); ) {
854
855 if ( std::abs( (par[0]-it->x())*sinphi + it->y()*cosphi ) > 2. ) { ++it; continue; }
856 countAssociatedHits[(it->y()>0)?0:1]++;
857 it = event_data.m_listHitCenter.erase( it );
858 }
859
860 if ( countAssociatedHits[0]+countAssociatedHits[1] < m_minHitsAboveTOT ) continue;
861
862 x0.push_back( par[0] ); phi.push_back( par[1] ); nHitsPosY.push_back( countAssociatedHits[0] ); nHitsNegY.push_back( countAssociatedHits[1] );
863
864 } // end for (int i=0; i<nIterations; i++) loop
865// end finding seeds
866
867// save all TRT hits for found segments
868 int nFoundSegments = x0.size();
869 if (nFoundSegments>10) nFoundSegments = 10; // number of found segments
870 for (int i=0; i<nFoundSegments; i++) {
871 event_data.m_segmentDriftCircles.emplace_back();
872 }
873
874 for (unsigned int i=0; i<event_data.m_listHits.size(); i++) {
875 const Amg::Vector3D &
876 center = event_data.m_listHits[i]->detectorElement()->surface(event_data.m_listHits[i]->identify()).center();
877
878 for (int j=0; j<nFoundSegments; j++) {
879 if (nHitsPosY[j]<5 && center.y()>0.) continue;
880 if (nHitsNegY[j]<5 && center.y()<0.) continue;
881 double cosphi = std::cos(phi[j]);
882 double sinphi = std::sqrt(1.-cosphi*cosphi);
883 if ( std::abs((x0[j]-center.x())*sinphi+center.y()*cosphi) < 2.) {
884 event_data.m_segmentDriftCircles[j].push_back(event_data.m_listHits[i]);
885 break;
886 }
887 }
888 }
889 // end saving TRT hits
890
891 if (msgLvl(MSG::DEBUG)) {
892 // debug: check how many hits per found segments
893 msg(MSG::DEBUG) << "find() debug (" << nFoundSegments << ")" ;
894 for (unsigned int i=0; i<event_data.m_segmentDriftCircles.size(); i++) {
895 msg(MSG::DEBUG) << " " << i << " " \
896 << event_data.m_segmentDriftCircles[i].size() ;
897 }
898 msg(MSG::DEBUG) << endmsg;
899 }
900
901 for (unsigned int i=0; i<event_data.m_segmentDriftCircles.size(); i++) { // convert to segments
902 if ((int)event_data.m_segmentDriftCircles[i].size()<m_minHitsForSegment) continue;
903 double trackpar[] = {x0[i], phi[i]};
904 convert(event_data.m_segmentDriftCircles[i], trackpar, event_data);
905 }
906 ATH_MSG_DEBUG( "find(), number of converted segments: " << event_data.m_segments.size() );
907}
908
909
910
911
912
913
#define M_PI
Scalar phi() const
phi method
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
std::pair< std::vector< unsigned int >, bool > res
static Double_t a
static Double_t sc
Handle class for reading from StoreGate.
This is an Identifier helper class for the TRT subdetector.
#define y
#define x
#define z
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
const ServiceHandle< StoreGateSvc > & detStore() const
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
Derived DataVector<T>.
Definition DataVector.h:795
This is a "hash" representation of an Identifier.
Virtual base class of TRT readout elements.
Represents 'corrected' measurements from the TRT (for example, corrected for wire sag).
virtual const InDetDD::TRT_BaseElement * detectorElement() const override final
return the detector element corresponding to this PRD
std::vector< std::vector< const InDet::TRT_DriftCircle * > > m_segmentDriftCircles
std::vector< Trk::TrackSegment * > m_segments
List of found segments.
void findSeedInverseR(double *par, TRT_TrackSegmentsMaker_BarrelCosmics::EventData &event_data) const
TRT_TrackSegmentsMaker_BarrelCosmics(const std::string &, const std::string &, const IInterface *)
virtual Trk::TrackSegment * next(InDet::ITRT_TrackSegmentsMaker::IEventData &event_data) const override
static void segFit(double *measx, double *measy, int nhits, double *residuals=0, double *result=0)
virtual std::unique_ptr< InDet::ITRT_TrackSegmentsMaker::IEventData > newRegion(const EventContext &ctx, const std::vector< IdentifierHash > &) const override
SG::ReadHandleKey< InDet::TRT_DriftCircleContainer > m_driftCirclesName
Container with TRT clusters.
int findSeed(double xmin, double xmax, double phimin, double phimax, double *bestParameters, TRT_TrackSegmentsMaker_BarrelCosmics::EventData &event_data) const
void findOld(TRT_TrackSegmentsMaker_BarrelCosmics::EventData &event_data) const
virtual void find(const EventContext &ctx, InDet::ITRT_TrackSegmentsMaker::IEventData &event_data, InDet::TRT_DetElementLink_xk::TRT_DetElemUsedMap &used) const override
void convert(std::vector< const InDet::TRT_DriftCircle * > &hits, double *trackpar, TRT_TrackSegmentsMaker_BarrelCosmics::EventData &event_data) const
void endEvent(InDet::ITRT_TrackSegmentsMaker::IEventData &event_data) const override
virtual std::unique_ptr< InDet::ITRT_TrackSegmentsMaker::IEventData > newEvent(const EventContext &ctx) const override
static bool sortHits(const InDet::TRT_DriftCircle *, const InDet::TRT_DriftCircle *)
sort hits on segment such that they are ordered from larger y to smaller (from top down)
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
static EventData & getPrivateEventData(InDet::ITRT_TrackSegmentsMaker::IEventData &virt_event_data)
Class to represent and store fit qualities from track reconstruction in terms of and number of degre...
Definition FitQuality.h:97
This class is the pure abstract base class for all fittable tracking measurements.
Identifier identify() const
return the identifier
Class to handle pseudo-measurements in fitters and on track objects.
Class for a StraightLineSurface in the ATLAS detector to describe dirft tube and straw like detectors...
Abstract Base Class for tracking surfaces.
const Amg::Transform3D & transform() const
Returns HepGeom::Transform3D by reference.
const Amg::Vector3D & center() const
Returns the center position of the Surface.
Class for a generic track segment that holdes polymorphic Trk::MeasurementBase objects,...
double chi2(TH1 *h0, TH1 *h1)
void mean(std::vector< double > &bins, std::vector< double > &values, const std::vector< std::string > &files, const std::string &histname, const std::string &tplotname, const std::string &label="")
double xmax
Definition listroot.cxx:61
double xmin
Definition listroot.cxx:60
Eigen::Matrix< double, 3, 3 > RotationMatrix3D
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Affine3d Transform3D
Eigen::Matrix< double, 3, 1 > Vector3D
Eigen::Translation< double, 3 > Translation3D
@ NODRIFTTIME
drift time was not used - drift radius is 0.
@ locX
Definition ParamDefs.h:37
@ theta
Definition ParamDefs.h:66
@ locZ
local cylindrical
Definition ParamDefs.h:42
std::pair< double, ParamDefs > DefinedParameter
Typedef to of a std::pair<double, ParamDefs> to identify a passed-through double as a specific type o...
Definition index.py:1
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
hold the test vectors and ease the comparison