CMS 3D CMS Logo

LaserAlignment.cc
Go to the documentation of this file.
1 
15 
20  theEvents(0),
21  theDoPedestalSubtraction( theConf.getUntrackedParameter<bool>( "SubtractPedestals", true ) ),
22  theUseMinuitAlgorithm( theConf.getUntrackedParameter<bool>( "RunMinuitAlignmentTubeAlgorithm", false ) ),
23  theApplyBeamKinkCorrections( theConf.getUntrackedParameter<bool>( "ApplyBeamKinkCorrections", true ) ),
24  peakFinderThreshold( theConf.getUntrackedParameter<double>( "PeakFinderThreshold", 10. ) ),
25  enableJudgeZeroFilter( theConf.getUntrackedParameter<bool>( "EnableJudgeZeroFilter", true ) ),
26  judgeOverdriveThreshold( theConf.getUntrackedParameter<unsigned int>( "JudgeOverdriveThreshold", 220 ) ),
27  updateFromInputGeometry( theConf.getUntrackedParameter<bool>( "UpdateFromInputGeometry", false ) ),
28  misalignedByRefGeometry( theConf.getUntrackedParameter<bool>( "MisalignedByRefGeometry", false ) ),
29  theStoreToDB ( theConf.getUntrackedParameter<bool>( "SaveToDbase", false ) ),
30  theDigiProducersList( theConf.getParameter<std::vector<edm::ParameterSet> >( "DigiProducersList" ) ),
31  theSaveHistograms( theConf.getUntrackedParameter<bool>( "SaveHistograms", false ) ),
32  theCompression( theConf.getUntrackedParameter<int>( "ROOTFileCompression", 1 ) ),
33  theFileName( theConf.getUntrackedParameter<std::string>( "ROOTFileName", "test.root" ) ),
34  theMaskTecModules( theConf.getUntrackedParameter<std::vector<unsigned int> >( "MaskTECModules" ) ),
35  theMaskAtModules( theConf.getUntrackedParameter<std::vector<unsigned int> >( "MaskATModules" ) ),
36  theSetNominalStrips( theConf.getUntrackedParameter<bool>( "ForceFitterToNominalStrips", false ) ),
37  theLasConstants( theConf.getUntrackedParameter<std::vector<edm::ParameterSet> >( "LaserAlignmentConstants" ) ),
38  theFile(),
39  theAlignableTracker(),
40  theAlignRecordName( "TrackerAlignmentRcd" ),
41  theErrorRecordName( "TrackerAlignmentErrorExtendedRcd" ),
42  firstEvent_(true)
43 {
44 
45 
46  std::cout << std::endl;
47  std::cout << "=============================================================="
48  << "\n=== LaserAlignment module configuration ==="
49  << "\n"
50  << "\n Write histograms to file = " << (theSaveHistograms?"true":"false")
51  << "\n Histogram file name = " << theFileName
52  << "\n Histogram file compression = " << theCompression
53  << "\n Subtract pedestals = " << (theDoPedestalSubtraction?"true":"false")
54  << "\n Run Minuit AT algorithm = " << (theUseMinuitAlgorithm?"true":"false")
55  << "\n Apply beam kink corrections = " << (theApplyBeamKinkCorrections?"true":"false")
56  << "\n Peak Finder Threshold = " << peakFinderThreshold
57  << "\n EnableJudgeZeroFilter = " << (enableJudgeZeroFilter?"true":"false")
58  << "\n JudgeOverdriveThreshold = " << judgeOverdriveThreshold
59  << "\n Update from input geometry = " << (updateFromInputGeometry?"true":"false")
60  << "\n Misalignment from ref geometry = " << (misalignedByRefGeometry?"true":"false")
61  << "\n Number of TEC modules masked = " << theMaskTecModules.size() << " (s. below list if > 0)"
62  << "\n Number of AT modules masked = " << theMaskAtModules.size() << " (s. below list if > 0)"
63  << "\n Store to database = " << (theStoreToDB?"true":"false")
64  << "\n ----------------------------------------------- ----------"
65  << (theSetNominalStrips?"\n Set strips to nominal = true":"\n")
66  << "\n=============================================================" << std::endl;
67 
68  // tell about masked modules
69  if( !theMaskTecModules.empty() ) {
70  std::cout << " ===============================================================================================\n" << std::flush;
71  std::cout << " The following " << theMaskTecModules.size() << " TEC modules have been masked out and will not be considered by the TEC algorithm:\n " << std::flush;
72  for( std::vector<unsigned int>::iterator moduleIt = theMaskTecModules.begin(); moduleIt != theMaskTecModules.end(); ++moduleIt ) {
73  std::cout << *moduleIt << (moduleIt!=--theMaskTecModules.end()?", ":"") << std::flush;
74  }
75  std::cout << std::endl << std::flush;
76  std::cout << " ===============================================================================================\n\n" << std::flush;
77  }
78  if( !theMaskAtModules.empty() ) {
79  std::cout << " ===============================================================================================\n" << std::flush;
80  std::cout << " The following " << theMaskAtModules.size() << " AT modules have been masked out and will not be considered by the AT algorithm:\n " << std::flush;
81  for( std::vector<unsigned int>::iterator moduleIt = theMaskAtModules.begin(); moduleIt != theMaskAtModules.end(); ++moduleIt ) {
82  std::cout << *moduleIt << (moduleIt!=--theMaskAtModules.end()?", ":"") << std::flush;
83  }
84  std::cout << std::endl << std::flush;
85  std::cout << " ===============================================================================================\n\n" << std::flush;
86  }
87 
88 
89 
90  // alias for the Branches in the root files
91  std::string alias ( theConf.getParameter<std::string>("@module_label") );
92 
93  // declare the product to produce
94  produces<TkLasBeamCollection, edm::Transition::EndRun>( "tkLaserBeams" ).setBranchAlias( alias + "TkLasBeamCollection" );
95 
96  // switch judge's zero filter depending on cfg
98 
99  // set the upper threshold for zero suppressed data
101 
102 }
103 
104 
105 
106 
107 
112 
113  if ( theSaveHistograms ) theFile->Write();
114  if ( theFile ) { delete theFile; }
115  if ( theAlignableTracker ) { delete theAlignableTracker; }
116 
117 }
118 
119 
120 
121 
122 
127 
128 
129  // write sumed histograms to file (if selected in cfg)
130  if( theSaveHistograms ) {
131 
132  // creating a new file
133  theFile = new TFile( theFileName.c_str(), "RECREATE", "CMS ROOT file" );
134 
135  // initialize the histograms
136  if ( theFile ) {
137  theFile->SetCompressionLevel(theCompression);
138  singleModulesDir = theFile->mkdir( "single modules" );
139  } else
140  throw cms::Exception( " [LaserAlignment::beginJob]") << " ** ERROR: could not open file:"
141  << theFileName.c_str() << " for writing." << std::endl;
142 
143  }
144 
145  // detector id maps (hard coded)
146  fillDetectorId();
147 
148  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
149  // PROFILE, HISTOGRAM & FITFUNCTION INITIALIZATION
150  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
151 
152  // object used to build various strings for names and labels
153  std::stringstream nameBuilder;
154 
155  // loop variables for use with LASGlobalLoop object
156  int det, ring, beam, disk, pos;
157 
158  // loop TEC modules
159  det = 0; ring = 0; beam = 0; disk = 0;
160  do { // loop using LASGlobalLoop functionality
161  // init the profiles
162  pedestalProfiles.GetTECEntry( det, ring, beam, disk ).SetAllValuesTo( 0. );
163  currentDataProfiles.GetTECEntry( det, ring, beam, disk ).SetAllValuesTo( 0. );
164  collectedDataProfiles.GetTECEntry( det, ring, beam, disk ).SetAllValuesTo( 0. );
165 
166  // init the hit maps
167  isAcceptedProfile.SetTECEntry( det, ring, beam, disk, 0 );
168  numberOfAcceptedProfiles.SetTECEntry( det, ring, beam, disk, 0 );
169 
170  // create strings for histo names
171  nameBuilder.clear();
172  nameBuilder.str( "" );
173  nameBuilder << "TEC";
174  if( det == 0 ) nameBuilder << "+"; else nameBuilder << "-";
175  nameBuilder << "_Ring";
176  if( ring == 0 ) nameBuilder << "4"; else nameBuilder << "6";
177  nameBuilder << "_Beam" << beam;
178  nameBuilder << "_Disk" << disk;
179  theProfileNames.SetTECEntry( det, ring, beam, disk, nameBuilder.str() );
180 
181  // init the histograms
182  if( theSaveHistograms ) {
183  nameBuilder << "_Histo";
184  summedHistograms.SetTECEntry( det, ring, beam, disk, new TH1D( nameBuilder.str().c_str(), nameBuilder.str().c_str(), 512, 0, 512 ) );
185  summedHistograms.GetTECEntry( det, ring, beam, disk )->SetDirectory( singleModulesDir );
186  }
187 
188  } while ( moduleLoop.TECLoop( det, ring, beam, disk ) );
189 
190 
191  // TIB & TOB section
192  det = 2; beam = 0; pos = 0;
193  do { // loop using LASGlobalLoop functionality
194  // init the profiles
195  pedestalProfiles.GetTIBTOBEntry( det, beam, pos ).SetAllValuesTo( 0. );
196  currentDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetAllValuesTo( 0. );
197  collectedDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetAllValuesTo( 0. );
198 
199  // init the hit maps
200  isAcceptedProfile.SetTIBTOBEntry( det, beam, pos, 0 );
201  numberOfAcceptedProfiles.SetTIBTOBEntry( det, beam, pos, 0 );
202 
203  // create strings for histo names
204  nameBuilder.clear();
205  nameBuilder.str( "" );
206  if( det == 2 ) nameBuilder << "TIB"; else nameBuilder << "TOB";
207  nameBuilder << "_Beam" << beam;
208  nameBuilder << "_Zpos" << pos;
209 
210  theProfileNames.SetTIBTOBEntry( det, beam, pos, nameBuilder.str() );
211 
212  // init the histograms
213  if( theSaveHistograms ) {
214  nameBuilder << "_Histo";
215  summedHistograms.SetTIBTOBEntry( det, beam, pos, new TH1D( nameBuilder.str().c_str(), nameBuilder.str().c_str(), 512, 0, 512 ) );
216  summedHistograms.GetTIBTOBEntry( det, beam, pos )->SetDirectory( singleModulesDir );
217  }
218 
219  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
220 
221 
222  // TEC2TEC AT section
223  det = 0; beam = 0; disk = 0;
224  do { // loop using LASGlobalLoop functionality
225  // init the profiles
226  pedestalProfiles.GetTEC2TECEntry( det, beam, disk ).SetAllValuesTo( 0. );
227  currentDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetAllValuesTo( 0. );
228  collectedDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetAllValuesTo( 0. );
229 
230  // init the hit maps
231  isAcceptedProfile.SetTEC2TECEntry( det, beam, disk, 0 );
232  numberOfAcceptedProfiles.SetTEC2TECEntry( det, beam, disk, 0 );
233 
234  // create strings for histo names
235  nameBuilder.clear();
236  nameBuilder.str( "" );
237  nameBuilder << "TEC(AT)";
238  if( det == 0 ) nameBuilder << "+"; else nameBuilder << "-";
239  nameBuilder << "_Beam" << beam;
240  nameBuilder << "_Disk" << disk;
241  theProfileNames.SetTEC2TECEntry( det, beam, disk, nameBuilder.str() );
242 
243  // init the histograms
244  if( theSaveHistograms ) {
245  nameBuilder << "_Histo";
246  summedHistograms.SetTEC2TECEntry( det, beam, disk, new TH1D( nameBuilder.str().c_str(), nameBuilder.str().c_str(), 512, 0, 512 ) );
247  summedHistograms.GetTEC2TECEntry( det, beam, disk )->SetDirectory( singleModulesDir );
248  }
249 
250  } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
251 
252  firstEvent_ = true;
253 }
254 
255 
256 
257 
258 
259 
263 void LaserAlignment::produce(edm::Event& theEvent, edm::EventSetup const& theSetup) {
264 
265  if (firstEvent_) {
266 
267  //Retrieve tracker topology from geometry
268  edm::ESHandle<TrackerTopology> tTopoHandle;
269  theSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
270  const TrackerTopology* const tTopo = tTopoHandle.product();
271 
272  // access the tracker
274  theSetup.get<IdealGeometryRecord>().get( gD );
275 
276  // access pedestals (from db..) if desired
277  edm::ESHandle<SiStripPedestals> pedestalsHandle;
279  theSetup.get<SiStripPedestalsRcd>().get( pedestalsHandle );
280  fillPedestalProfiles( pedestalsHandle );
281  }
282 
283  // global positions
284  // edm::ESHandle<Alignments> theGlobalPositionRcd;
285  theSetup.get<TrackerDigiGeometryRecord>().getRecord<GlobalPositionRcd>().get( theGlobalPositionRcd );
286 
287  // select the reference geometry
288  if( !updateFromInputGeometry ) {
289  // the AlignableTracker object is initialized with the ideal geometry
290  edm::ESHandle<GeometricDet> theGeometricDet;
291  theSetup.get<IdealGeometryRecord>().get(theGeometricDet);
293  theSetup.get<PTrackerParametersRcd>().get( ptp );
294  TrackerGeomBuilderFromGeometricDet trackerBuilder;
295  TrackerGeometry* theRefTracker = trackerBuilder.build(&*theGeometricDet, *ptp, tTopo );
296 
297  theAlignableTracker = new AlignableTracker(&(*theRefTracker), tTopo);
298  }
299  else {
300  // the AlignableTracker object is initialized with the input geometry from DB
302  }
303 
304  firstEvent_ = false;
305  }
306 
307  LogDebug("LaserAlignment") << "==========================================================="
308  << "\n Private analysis of event #"<< theEvent.id().event()
309  << " in run #" << theEvent.id().run();
310 
311 
312  // do the Tracker Statistics to retrieve the current profiles
313  fillDataProfiles( theEvent, theSetup );
314 
315  // index variables for the LASGlobalLoop object
316  int det, ring, beam, disk, pos;
317 
318  //
319  // first pre-loop on selected entries to find out
320  // whether the TEC or the AT beams have fired
321  // (pedestal profiles are left empty if false in cfg)
322  //
323 
324 
325  // TEC+- (only ring 6)
326  ring = 1;
327  for( det = 0; det < 2; ++det ) {
328  for( beam = 0; beam < 8; ++ beam ) {
329  for( disk = 0; disk < 9; ++disk ) {
330  if( judge.IsSignalIn( currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk ), 0 ) ) {
331  isAcceptedProfile.SetTECEntry( det, ring, beam, disk, 1 );
332  }
333  else { // assume no initialization
334  isAcceptedProfile.SetTECEntry( det, ring, beam, disk, 0 );
335  }
336  }
337  }
338  }
339 
340  // TIBTOB
341  det = 2; beam = 0; pos = 0;
342  do {
343  // add current event's data and subtract pedestals
344  if( judge.IsSignalIn( currentDataProfiles.GetTIBTOBEntry( det, beam, pos ) - pedestalProfiles.GetTIBTOBEntry( det, beam, pos ), getTIBTOBNominalBeamOffset( det, beam, pos ) ) ) {
345  isAcceptedProfile.SetTIBTOBEntry( det, beam, pos, 1 );
346  }
347  else { // dto.
348  isAcceptedProfile.SetTIBTOBEntry( det, beam, pos, 0 );
349  }
350 
351  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
352 
353 
354 
355 
356  // now come the beam finders
357  bool isTECMode = isTECBeam();
358  // LogDebug( " [LaserAlignment::produce]" ) << "LaserAlignment::isTECBeam declares this event " << ( isTECMode ? "" : "NOT " ) << "a TEC event." << std::endl;
359  std::cout << " [LaserAlignment::produce] -- LaserAlignment::isTECBeam declares this event " << ( isTECMode ? "" : "NOT " ) << "a TEC event." << std::endl;
360 
361  bool isATMode = isATBeam();
362  // LogDebug( " [LaserAlignment::produce]" ) << "LaserAlignment::isATBeam declares this event " << ( isATMode ? "" : "NOT " ) << "an AT event." << std::endl;
363  std::cout << " [LaserAlignment::produce] -- LaserAlignment::isATBeam declares this event " << ( isATMode ? "" : "NOT " ) << "an AT event." << std::endl;
364 
365 
366 
367 
368  //
369  // now pass the pedestal subtracted profiles to the judge
370  // if they're accepted, add them on the collectedDataProfiles
371  // (pedestal profiles are left empty if false in cfg)
372  //
373 
374 
375  // loop TEC+- modules
376  det = 0; ring = 0; beam = 0; disk = 0;
377  do {
378 
379  LogDebug( "[LaserAlignment::produce]" ) << "Profile is: " << theProfileNames.GetTECEntry( det, ring, beam, disk ) << "." << std::endl;
380 
381  // this now depends on the AT/TEC mode, is this a doubly hit module? -> look for it in vector<int> tecDoubleHitDetId
382  // (ring == 0 not necessary but makes it a little faster)
383  if( ring == 0 && find( tecDoubleHitDetId.begin(), tecDoubleHitDetId.end(), detectorId.GetTECEntry( det, ring, beam, disk ) ) != tecDoubleHitDetId.end() ) {
384 
385  if( isTECMode ) { // add profile to TEC collection
386  // add current event's data and subtract pedestals
387  if( judge.JudgeProfile( currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk ), 0 ) ) {
388  collectedDataProfiles.GetTECEntry( det, ring, beam, disk ) += currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk );
389  numberOfAcceptedProfiles.GetTECEntry( det, ring, beam, disk )++;
390  }
391  }
392  }
393 
394  else { // not a doubly hit module, don't care about the mode
395  // add current event's data and subtract pedestals
396  if( judge.JudgeProfile( currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk ), 0 ) ) {
397  collectedDataProfiles.GetTECEntry( det, ring, beam, disk ) += currentDataProfiles.GetTECEntry( det, ring, beam, disk ) - pedestalProfiles.GetTECEntry( det, ring, beam, disk );
398  numberOfAcceptedProfiles.GetTECEntry( det, ring, beam, disk )++;
399  }
400  }
401 
402  } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
403 
404 
405 
406 
407 
408  // loop TIB/TOB modules
409  det = 2; beam = 0; pos = 0;
410  do {
411 
412  LogDebug( "[LaserAlignment::produce]" ) << "Profile is: " << theProfileNames.GetTIBTOBEntry( det, beam, pos ) << "." << std::endl;
413 
414  // add current event's data and subtract pedestals
415  if( judge.JudgeProfile( currentDataProfiles.GetTIBTOBEntry( det, beam, pos ) - pedestalProfiles.GetTIBTOBEntry( det, beam, pos ), getTIBTOBNominalBeamOffset( det, beam, pos ) ) ) {
416  collectedDataProfiles.GetTIBTOBEntry( det, beam, pos ) += currentDataProfiles.GetTIBTOBEntry( det, beam, pos ) - pedestalProfiles.GetTIBTOBEntry( det, beam, pos );
417  numberOfAcceptedProfiles.GetTIBTOBEntry( det, beam, pos )++;
418  }
419 
420  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
421 
422 
423 
424  // loop TEC2TEC modules
425  det = 0; beam = 0; disk = 0;
426  do {
427 
428  LogDebug( "[LaserAlignment::produce]" ) << "Profile is: " << theProfileNames.GetTEC2TECEntry( det, beam, disk ) << "." << std::endl;
429 
430  // this again depends on the AT/TEC mode, is this a doubly hit module?
431  // (ring == 0 not necessary but makes it a little faster)
432  if( ring == 0 && find( tecDoubleHitDetId.begin(), tecDoubleHitDetId.end(), detectorId.GetTECEntry( det, ring, beam, disk ) ) != tecDoubleHitDetId.end() ) {
433 
434  if( isATMode ) { // add profile to TEC2TEC collection
435  // add current event's data and subtract pedestals
436  if( judge.JudgeProfile( currentDataProfiles.GetTEC2TECEntry( det, beam, disk ) - pedestalProfiles.GetTEC2TECEntry( det, beam, disk ), 0 ) ) {
437  collectedDataProfiles.GetTEC2TECEntry( det, beam, disk ) += currentDataProfiles.GetTEC2TECEntry( det, beam, disk ) - pedestalProfiles.GetTEC2TECEntry( det, beam, disk );
438  numberOfAcceptedProfiles.GetTEC2TECEntry( det, beam, disk )++;
439  }
440  }
441 
442  }
443 
444  else { // not a doubly hit module, don't care about the mode
445  // add current event's data and subtract pedestals
446  if( judge.JudgeProfile( currentDataProfiles.GetTEC2TECEntry( det, beam, disk ) - pedestalProfiles.GetTEC2TECEntry( det, beam, disk ), 0 ) ) {
447  collectedDataProfiles.GetTEC2TECEntry( det, beam, disk ) += currentDataProfiles.GetTEC2TECEntry( det, beam, disk ) - pedestalProfiles.GetTEC2TECEntry( det, beam, disk );
448  numberOfAcceptedProfiles.GetTEC2TECEntry( det, beam, disk )++;
449  }
450  }
451 
452 
453  } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
454 
455 
456 
457  // total event number counter
458  theEvents++;
459 
460 }
461 
462 
463 
464 
465 
469 void LaserAlignment::endRunProduce( edm::Run& theRun, const edm::EventSetup& theSetup ) {
470 
471 
472  std::cout << " [LaserAlignment::endRun] -- Total number of events processed: " << theEvents << std::endl;
473 
474  // for debugging only..
476 
477  // index variables for the LASGlobalLoop objects
478  int det, ring, beam, disk, pos;
479 
480  // measured positions container for the algorithms
481  LASGlobalData<LASCoordinateSet> measuredCoordinates;
482 
483  // fitted peak positions in units of strips (pair for value,error)
484  LASGlobalData<std::pair<float,float> > measuredStripPositions;
485 
486  // the peak finder, a pair (pos/posErr in units of strips) for its results, and the success confirmation
487  LASPeakFinder peakFinder;
489  std::pair<double,double> peakFinderResults;
490  bool isGoodFit;
491 
492  // tracker geom. object for calculating the global beam positions
493  const TrackerGeometry& theTracker( *theTrackerGeometry );
494 
495  // fill LASGlobalData<LASCoordinateSet> nominalCoordinates
497 
498  // for determining the phi errors
499  // ErrorFrameTransformer errorTransformer; // later...
500 
501 
502 
503 
504 
505  // do the fits for TEC+- internal
506  det = 0; ring = 0; beam = 0; disk = 0;
507  do {
508 
509  // do the fit
510  isGoodFit = peakFinder.FindPeakIn( collectedDataProfiles.GetTECEntry( det, ring, beam, disk ), peakFinderResults,
511  summedHistograms.GetTECEntry( det, ring, beam, disk ), 0 ); // offset is 0 for TEC
512 
513  // now we have the measured positions in units of strips.
514  if( !isGoodFit ) std::cout << " [LaserAlignment::endRun] ** WARNING: Fit failed for TEC det: "
515  << det << ", ring: " << ring << ", beam: " << beam << ", disk: " << disk
516  << " (id: " << detectorId.GetTECEntry( det, ring, beam, disk ) << ")." << std::endl;
517 
518 
519 
520  // <- here we will later implement the kink corrections
521 
522  // access the tracker geometry for this module
523  const DetId theDetId( detectorId.GetTECEntry( det, ring, beam, disk ) );
524  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
525 
526  if (theStripDet) {
527  // first, set the measured coordinates to their nominal values
528  measuredCoordinates.SetTECEntry( det, ring, beam, disk, nominalCoordinates.GetTECEntry( det, ring, beam, disk ) );
529 
530  if( isGoodFit ) { // convert strip position to global phi and replace the nominal phi value/error
531 
532  measuredStripPositions.GetTECEntry( det, ring, beam, disk ) = peakFinderResults;
533  const float positionInStrips = theSetNominalStrips ? 256. : peakFinderResults.first; // implementation of "ForceFitterToNominalStrips" config parameter
534  const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( positionInStrips ) );
535  measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
536 
537  // const GlobalError& globalError = errorTransformer.transform( theStripDet->specificTopology().localError( peakFinderResults.first, pow( peakFinderResults.second, 2 ) ), theStripDet->surface() );
538  // measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhiError( globalError.phierr( globalPoint ) );
539  measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhiError( 0.00046 ); // PRELIMINARY ESTIMATE
540 
541  }
542  else { // keep nominal position (middle-of-module) but set a giant phi error so that the module can be ignored by the alignment algorithm
543  measuredStripPositions.GetTECEntry( det, ring, beam, disk ) = std::pair<float,float>( 256., 1000. );
544  const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( 256. ) );
545  measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
546  measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhiError( 1000. );
547  }
548  }
549 
550  } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
551 
552 
553 
554 
555  // do the fits for TIB/TOB
556  det = 2; beam = 0; pos = 0;
557  do {
558 
559  // do the fit
560  isGoodFit = peakFinder.FindPeakIn( collectedDataProfiles.GetTIBTOBEntry( det, beam, pos ), peakFinderResults,
561  summedHistograms.GetTIBTOBEntry( det, beam, pos ), getTIBTOBNominalBeamOffset( det, beam, pos ) );
562 
563  // now we have the measured positions in units of strips.
564  if( !isGoodFit ) std::cout << " [LaserAlignment::endJob] ** WARNING: Fit failed for TIB/TOB det: "
565  << det << ", beam: " << beam << ", pos: " << pos
566  << " (id: " << detectorId.GetTIBTOBEntry( det, beam, pos ) << ")." << std::endl;
567 
568 
569  // <- here we will later implement the kink corrections
570 
571  // access the tracker geometry for this module
572  const DetId theDetId( detectorId.GetTIBTOBEntry( det, beam, pos ) );
573  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
574 
575  if (theStripDet) {
576  // first, set the measured coordinates to their nominal values
577  measuredCoordinates.SetTIBTOBEntry( det, beam, pos, nominalCoordinates.GetTIBTOBEntry( det, beam, pos ) );
578 
579  if( isGoodFit ) { // convert strip position to global phi and replace the nominal phi value/error
580  measuredStripPositions.GetTIBTOBEntry( det, beam, pos ) = peakFinderResults;
581  const float positionInStrips = theSetNominalStrips ? 256. + getTIBTOBNominalBeamOffset( det, beam, pos ) : peakFinderResults.first; // implementation of "ForceFitterToNominalStrips" config parameter
582  const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( positionInStrips ) );
583  measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
584  measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhiError( 0.00028 ); // PRELIMINARY ESTIMATE
585  }
586  else { // keep nominal position but set a giant phi error so that the module can be ignored by the alignment algorithm
587  measuredStripPositions.GetTIBTOBEntry( det, beam, pos ) = std::pair<float,float>( 256. + getTIBTOBNominalBeamOffset( det, beam, pos ), 1000. );
588  const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( 256. + getTIBTOBNominalBeamOffset( det, beam, pos ) ) );
589  measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
590  measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhiError( 1000. );
591  }
592  }
593 
594  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
595 
596 
597 
598 
599  // do the fits for TEC AT
600  det = 0; beam = 0; disk = 0;
601  do {
602 
603  // do the fit
604  isGoodFit = peakFinder.FindPeakIn( collectedDataProfiles.GetTEC2TECEntry( det, beam, disk ), peakFinderResults,
605  summedHistograms.GetTEC2TECEntry( det, beam, disk ), getTEC2TECNominalBeamOffset( det, beam, disk ) );
606  // now we have the positions in units of strips.
607  if( !isGoodFit ) std::cout << " [LaserAlignment::endRun] ** WARNING: Fit failed for TEC2TEC det: "
608  << det << ", beam: " << beam << ", disk: " << disk
609  << " (id: " << detectorId.GetTEC2TECEntry( det, beam, disk ) << ")." << std::endl;
610 
611 
612  // <- here we will later implement the kink corrections
613 
614  // access the tracker geometry for this module
615  const DetId theDetId( detectorId.GetTEC2TECEntry( det, beam, disk ) );
616  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
617 
618  if (theStripDet) {
619  // first, set the measured coordinates to their nominal values
620  measuredCoordinates.SetTEC2TECEntry( det, beam, disk, nominalCoordinates.GetTEC2TECEntry( det, beam, disk ) );
621 
622  if( isGoodFit ) { // convert strip position to global phi and replace the nominal phi value/error
623  measuredStripPositions.GetTEC2TECEntry( det, beam, disk ) = peakFinderResults;
624  const float positionInStrips = theSetNominalStrips ? 256. + getTEC2TECNominalBeamOffset( det, beam, disk ) : peakFinderResults.first; // implementation of "ForceFitterToNominalStrips" config parameter
625  const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( positionInStrips ) );
626  measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
627  measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhiError( 0.00047 ); // PRELIMINARY ESTIMATE
628  }
629  else { // keep nominal position but set a giant phi error so that the module can be ignored by the alignment algorithm
630  measuredStripPositions.GetTEC2TECEntry( det, beam, disk ) = std::pair<float,float>( 256. + getTEC2TECNominalBeamOffset( det, beam, disk ), 1000. );
631  const GlobalPoint& globalPoint = theStripDet->surface().toGlobal( theStripDet->specificTopology().localPosition( 256. + getTEC2TECNominalBeamOffset( det, beam, disk ) ) );
632  measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhi( ConvertAngle( globalPoint.barePhi() ) );
633  measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhiError( 1000. );
634  }
635  }
636 
637  } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
638 
639 
640 
641 
642 
643 
644 
645  // see what we got (for debugging)
646  // DumpStripFileSet( measuredStripPositions );
647  // DumpPosFileSet( measuredCoordinates );
648 
649 
650 
651 
652 
653 
654  // CALCULATE PARAMETERS AND UPDATE DB OBJECT
655  // for beam kink corrections, reconstructing the geometry and updating the db object
657 
658  // apply all beam corrections
659  if( theApplyBeamKinkCorrections ) geometryUpdater.ApplyBeamKinkCorrections( measuredCoordinates );
660 
661  // if we start with input geometry instead of IDEAL,
662  // reverse the adjustments in the AlignableTracker object
663  if( updateFromInputGeometry ) geometryUpdater.SetReverseDirection( true );
664 
665  // if we have "virtual" misalignment which is introduced via the reference geometry,
666  // tell the LASGeometryUpdater to reverse x & y adjustments
667  if( misalignedByRefGeometry ) geometryUpdater.SetMisalignmentFromRefGeometry( true );
668 
669  // run the endcap algorithm
670  LASEndcapAlgorithm endcapAlgorithm;
671  LASEndcapAlignmentParameterSet endcapParameters;
672 
673 
674  // this basically sets all the endcap modules to be masked
675  // to their nominal positions (since endcapParameters is overall zero)
676  if( !theMaskTecModules.empty() ) {
677  ApplyEndcapMaskingCorrections( measuredCoordinates, nominalCoordinates, endcapParameters );
678  }
679 
680  // run the algorithm
681  endcapParameters = endcapAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
682 
683  //
684  // loop to mask out events
685  // DESCRIPTION:
686  //
687 
688  // do this only if there are modules to be masked..
689  if( !theMaskTecModules.empty() ) {
690 
691  const unsigned int nIterations = 30;
692  for( unsigned int iteration = 0; iteration < nIterations; ++iteration ) {
693 
694  // set the endcap modules to be masked to their positions
695  // according to the reconstructed parameters
696  ApplyEndcapMaskingCorrections( measuredCoordinates, nominalCoordinates, endcapParameters );
697 
698  // modifications applied, so re-run the algorithm
699  endcapParameters = endcapAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
700 
701  }
702 
703  }
704 
705  // these are now final, so:
706  endcapParameters.Print();
707 
708 
709 
710 
711  // do a pre-alignment of the endcaps (TEC2TEC only)
712  // so that the alignment tube algorithms finds orderly disks
713  geometryUpdater.EndcapUpdate( endcapParameters, measuredCoordinates );
714 
715 
716  // the alignment tube algorithms, choose from config
717  LASBarrelAlignmentParameterSet alignmentTubeParameters;
718  // the MINUIT-BASED alignment tube algorithm
719  LASBarrelAlgorithm barrelAlgorithm;
720  // the ANALYTICAL alignment tube algorithm
721  LASAlignmentTubeAlgorithm alignmentTubeAlgorithm;
722 
723 
724  // this basically sets all the modules to be masked
725  // to their nominal positions (since alignmentTubeParameters is overall zero)
726  if( !theMaskAtModules.empty() ) {
727  ApplyATMaskingCorrections( measuredCoordinates, nominalCoordinates, alignmentTubeParameters );
728  }
729 
730  if( theUseMinuitAlgorithm ) {
731  // run the MINUIT-BASED alignment tube algorithm
732  alignmentTubeParameters = barrelAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
733  }
734  else {
735  // the ANALYTICAL alignment tube algorithm
736  alignmentTubeParameters = alignmentTubeAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
737  }
738 
739 
740 
741  //
742  // loop to mask out events
743  // DESCRIPTION:
744  //
745 
746  // do this only if there are modules to be masked..
747  if( !theMaskAtModules.empty() ) {
748 
749  const unsigned int nIterations = 30;
750  for( unsigned int iteration = 0; iteration < nIterations; ++iteration ) {
751 
752  // set the AT modules to be masked to their positions
753  // according to the reconstructed parameters
754  ApplyATMaskingCorrections( measuredCoordinates, nominalCoordinates, alignmentTubeParameters );
755 
756  // modifications applied, so re-run the algorithm
757  if( theUseMinuitAlgorithm ) {
758  alignmentTubeParameters = barrelAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
759  }
760  else {
761  alignmentTubeParameters = alignmentTubeAlgorithm.CalculateParameters( measuredCoordinates, nominalCoordinates );
762  }
763 
764  }
765 
766  }
767 
768 
769  // these are now final, so:
770  alignmentTubeParameters.Print();
771 
772 
773  // combine the results and update the db object
774  geometryUpdater.TrackerUpdate( endcapParameters, alignmentTubeParameters, *theAlignableTracker );
775 
776 
777 
778 
779 
780 
781 
786 
787 
788  // the collection container
789  auto laserBeams = std::make_unique<TkLasBeamCollection>();
790 
791 
792  // first for the endcap internal beams
793  for( det = 0; det < 2; ++det ) {
794  for( ring = 0; ring < 2; ++ring ) {
795  for( beam = 0; beam < 8; ++beam ) {
796 
797  // the beam and its identifier (see TkLasTrackBasedInterface TWiki)
798  TkLasBeam currentBeam( 100 * det + 10 * beam + ring );
799 
800  // order the hits in the beam by increasing z
801  const int firstDisk = det==0 ? 0 : 8;
802  const int lastDisk = det==0 ? 8 : 0;
803 
804  // count upwards or downwards
805  for( disk = firstDisk; det==0 ? disk <= lastDisk : disk >= lastDisk; det==0 ? ++disk : --disk ) {
806 
807  // detId for the SiStripLaserRecHit2D
808  const SiStripDetId theDetId( detectorId.GetTECEntry( det, ring, beam, disk ) );
809 
810  // need this to calculate the localPosition and its error
811  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
812 
813  // the hit container
814  const SiStripLaserRecHit2D currentHit(
815  theStripDet->specificTopology().localPosition( measuredStripPositions.GetTECEntry( det, ring, beam, disk ).first ),
816  theStripDet->specificTopology().localError( measuredStripPositions.GetTECEntry( det, ring, beam, disk ).first, measuredStripPositions.GetTECEntry( det, ring, beam, disk ).second ),
817  theDetId
818  );
819 
820  currentBeam.push_back( currentHit );
821 
822  }
823 
824  laserBeams->push_back( currentBeam );
825 
826  }
827  }
828  }
829 
830 
831 
832  // then, following the convention in TkLasTrackBasedInterface TWiki, the alignment tube beams;
833  // they comprise hits in TIBTOB & TEC2TEC
834 
835  for( beam = 0; beam < 8; ++beam ) {
836 
837  // the beam and its identifier (see TkLasTrackBasedInterface TWiki)
838  TkLasBeam currentBeam( 100 * 2 /*beamGroup=AT=2*/ + 10 * beam + 0 /*ring=0*/);
839 
840 
841  // first: tec-
842  det = 1;
843  for( disk = 4; disk >= 0; --disk ) {
844 
845  // detId for the SiStripLaserRecHit2D
846  const SiStripDetId theDetId( detectorId.GetTEC2TECEntry( det, beam, disk ) );
847 
848  // need this to calculate the localPosition and its error
849  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
850 
851  // the hit container
852  const SiStripLaserRecHit2D currentHit(
853  theStripDet->specificTopology().localPosition( measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first ),
854  theStripDet->specificTopology().localError( measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first, measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).second ),
855  theDetId
856  );
857 
858  currentBeam.push_back( currentHit );
859 
860  }
861 
862 
863  // now TIB and TOB in one go
864  for( det = 2; det < 4; ++det ) {
865  for( pos = 5; pos >= 0; --pos ) { // stupidly, pos is defined from +z to -z in LASGlobalLoop
866 
867  // detId for the SiStripLaserRecHit2D
868  const SiStripDetId theDetId( detectorId.GetTIBTOBEntry( det, beam, pos ) );
869 
870  // need this to calculate the localPosition and its error
871  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
872 
873  // the hit container
874  const SiStripLaserRecHit2D currentHit(
875  theStripDet->specificTopology().localPosition( measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).first ),
876  theStripDet->specificTopology().localError( measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).first, measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).second ),
877  theDetId
878  );
879 
880  currentBeam.push_back( currentHit );
881 
882  }
883  }
884 
885 
886  // then: tec+
887  det = 0;
888  for( disk = 0; disk < 5; ++disk ) {
889 
890  // detId for the SiStripLaserRecHit2D
891  const SiStripDetId theDetId( detectorId.GetTEC2TECEntry( det, beam, disk ) );
892 
893  // need this to calculate the localPosition and its error
894  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
895 
896  // the hit container
897  const SiStripLaserRecHit2D currentHit(
898  theStripDet->specificTopology().localPosition( measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first ),
899  theStripDet->specificTopology().localError( measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first, measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).second ),
900  theDetId
901  );
902 
903  currentBeam.push_back( currentHit );
904 
905  }
906 
907 
908 
909  // save this beam to the beamCollection
910  laserBeams->push_back( currentBeam );
911 
912  } // (close beam loop)
913 
914 
915  // now attach the collection to the run
916  theRun.put(std::move(laserBeams), "tkLaserBeams" );
917 
918 
919 
920 
921 
922  // store the estimated alignment parameters into the DB
923  // first get them
924  Alignments* alignments = theAlignableTracker->alignments();
926 
927  if ( theStoreToDB ) {
928 
929  std::cout << " [LaserAlignment::endRun] -- Storing the calculated alignment parameters to the DataBase:" << std::endl;
930 
931  // Call service
933  if( !poolDbService.isAvailable() ) // Die if not available
934  throw cms::Exception( "NotAvailable" ) << "PoolDBOutputService not available";
935 
936  // Store
937 
938  // if ( poolDbService->isNewTagRequest(theAlignRecordName) ) {
939  // poolDbService->createNewIOV<Alignments>( alignments, poolDbService->currentTime(), poolDbService->endOfTime(), theAlignRecordName );
940  // }
941  // else {
942  // poolDbService->appendSinceTime<Alignments>( alignments, poolDbService->currentTime(), theAlignRecordName );
943  // }
944  poolDbService->writeOne<Alignments>( alignments, poolDbService->beginOfTime(), theAlignRecordName );
945 
946  // if ( poolDbService->isNewTagRequest(theErrorRecordName) ) {
947  // poolDbService->createNewIOV<AlignmentErrorsExtended>( alignmentErrors, poolDbService->currentTime(), poolDbService->endOfTime(), theErrorRecordName );
948  // }
949  // else {
950  // poolDbService->appendSinceTime<AlignmentErrorsExtended>( alignmentErrors, poolDbService->currentTime(), theErrorRecordName );
951  // }
952  poolDbService->writeOne<AlignmentErrorsExtended>( alignmentErrors, poolDbService->beginOfTime(), theErrorRecordName );
953 
954  std::cout << " [LaserAlignment::endRun] -- Storing done." << std::endl;
955 
956  }
957 
958 }
959 
960 
961 
962 
963 
968 }
969 
970 
971 
972 
973 
978 void LaserAlignment::fillDataProfiles( edm::Event const& theEvent, edm::EventSetup const& theSetup ) {
979 
980  // two handles for the two different kinds of digis
983 
984  bool isRawDigi = false;
985 
986  // indices for the LASGlobalLoop object
987  int det = 0, ring = 0, beam = 0, disk = 0, pos = 0;
988 
989  // query config set and loop over all PSets in the VPSet
990  for ( std::vector<edm::ParameterSet>::iterator itDigiProducersList = theDigiProducersList.begin(); itDigiProducersList != theDigiProducersList.end(); ++itDigiProducersList ) {
991 
992  std::string digiProducer = itDigiProducersList->getParameter<std::string>( "DigiProducer" );
993  std::string digiLabel = itDigiProducersList->getParameter<std::string>( "DigiLabel" );
994  std::string digiType = itDigiProducersList->getParameter<std::string>( "DigiType" );
995 
996  // now branch according to digi type (raw or processed);
997  // first we go for raw digis => SiStripRawDigi
998  if( digiType == "Raw" ) {
999  theEvent.getByLabel( digiProducer, digiLabel, theStripRawDigis );
1000  isRawDigi = true;
1001  }
1002  else if( digiType == "Processed" ) {
1003  theEvent.getByLabel( digiProducer, digiLabel, theStripDigis );
1004  isRawDigi = false;
1005  }
1006  else {
1007  throw cms::Exception( " [LaserAlignment::fillDataProfiles]") << " ** ERROR: Invalid digi type: \"" << digiType << "\" specified in configuration." << std::endl;
1008  }
1009 
1010 
1011 
1012  // loop TEC internal modules
1013  det = 0; ring = 0; beam = 0; disk = 0;
1014  do {
1015 
1016  // first clear the profile
1017  currentDataProfiles.GetTECEntry( det, ring, beam, disk ).SetAllValuesTo( 0. );
1018 
1019  // retrieve the raw id of that module
1020  const int detRawId = detectorId.GetTECEntry( det, ring, beam, disk );
1021 
1022  if( isRawDigi ) { // we have raw SiStripRawDigis
1023 
1024  // search the digis for the raw id
1025  edm::DetSetVector<SiStripRawDigi>::const_iterator detSetIter = theStripRawDigis->find( detRawId );
1026  if( detSetIter == theStripRawDigis->end() ) {
1027  throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: No raw DetSet found for det: " << detRawId << "." << std::endl;
1028  }
1029 
1030  // fill the digis to the profiles
1031  edm::DetSet<SiStripRawDigi>::const_iterator digiRangeIterator = detSetIter->data.begin(); // for the loop
1032  edm::DetSet<SiStripRawDigi>::const_iterator digiRangeStart = digiRangeIterator; // save starting positions
1033 
1034  // loop all digis
1035  for (; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
1036  const SiStripRawDigi& digi = *digiRangeIterator;
1037  const int channel = distance( digiRangeStart, digiRangeIterator );
1038  if ( channel >= 0 && channel < 512 ) currentDataProfiles.GetTECEntry( det, ring, beam, disk ).SetValue( channel, digi.adc() );
1039  else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: raw digi channel: " << channel << " out of range for det: " << detRawId << "." << std::endl;
1040  }
1041 
1042  }
1043 
1044  else { // we have zero suppressed SiStripDigis
1045 
1046  // search the digis for the raw id
1047  edm::DetSetVector<SiStripDigi>::const_iterator detSetIter = theStripDigis->find( detRawId );
1048 
1049  // processed DetSets may be missing, just skip
1050  if( detSetIter == theStripDigis->end() ) continue;
1051 
1052  // fill the digis to the profiles
1053  edm::DetSet<SiStripDigi>::const_iterator digiRangeIterator = detSetIter->data.begin(); // for the loop
1054 
1055  for(; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
1056  const SiStripDigi& digi = *digiRangeIterator;
1057  if ( digi.strip() < 512 ) currentDataProfiles.GetTECEntry( det, ring, beam, disk ).SetValue( digi.strip(), digi.adc() );
1058  else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: digi strip: " << digi.strip() << " out of range for det: " << detRawId << "." << std::endl;
1059  }
1060 
1061  }
1062 
1063 
1064  } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
1065 
1066 
1067 
1068 
1069 
1070  // loop TIBTOB modules
1071  det = 2; beam = 0; pos = 0;
1072  do {
1073 
1074  // first clear the profile
1076 
1077  // retrieve the raw id of that module
1078  const int detRawId = detectorId.GetTIBTOBEntry( det, beam, pos );
1079 
1080  if( isRawDigi ) { // we have raw SiStripRawDigis
1081 
1082  // search the digis for the raw id
1083  edm::DetSetVector<SiStripRawDigi>::const_iterator detSetIter = theStripRawDigis->find( detRawId );
1084  if( detSetIter == theStripRawDigis->end() ) {
1085  throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: No raw DetSet found for det: " << detRawId << "." << std::endl;
1086  }
1087 
1088  // fill the digis to the profiles
1089  edm::DetSet<SiStripRawDigi>::const_iterator digiRangeIterator = detSetIter->data.begin(); // for the loop
1090  edm::DetSet<SiStripRawDigi>::const_iterator digiRangeStart = digiRangeIterator; // save starting positions
1091 
1092  // loop all digis
1093  for (; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
1094  const SiStripRawDigi& digi = *digiRangeIterator;
1095  const int channel = distance( digiRangeStart, digiRangeIterator );
1096  if ( channel >= 0 && channel < 512 ) currentDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetValue( channel, digi.adc() );
1097  else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: raw digi channel: " << channel << " out of range for det: " << detRawId << "." << std::endl;
1098  }
1099 
1100  }
1101 
1102  else { // we have zero suppressed SiStripDigis
1103 
1104  // search the digis for the raw id
1105  edm::DetSetVector<SiStripDigi>::const_iterator detSetIter = theStripDigis->find( detRawId );
1106 
1107  // processed DetSets may be missing, just skip
1108  if( detSetIter == theStripDigis->end() ) continue;
1109 
1110  // fill the digis to the profiles
1111  edm::DetSet<SiStripDigi>::const_iterator digiRangeIterator = detSetIter->data.begin(); // for the loop
1112 
1113  for(; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
1114  const SiStripDigi& digi = *digiRangeIterator;
1115  if ( digi.strip() < 512 ) currentDataProfiles.GetTIBTOBEntry( det, beam, pos ).SetValue( digi.strip(), digi.adc() );
1116  else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: digi strip: " << digi.strip() << " out of range for det: " << detRawId << "." << std::endl;
1117  }
1118 
1119  }
1120 
1121  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
1122 
1123 
1124 
1125  // loop TEC AT modules
1126  det = 0; beam = 0; disk = 0;
1127  do {
1128 
1129  // first clear the profile
1130  currentDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetAllValuesTo( 0. );
1131 
1132  // retrieve the raw id of that module
1133  const int detRawId = detectorId.GetTEC2TECEntry( det, beam, disk );
1134 
1135  if( isRawDigi ) { // we have raw SiStripRawDigis
1136 
1137  // search the digis for the raw id
1138  edm::DetSetVector<SiStripRawDigi>::const_iterator detSetIter = theStripRawDigis->find( detRawId );
1139  if( detSetIter == theStripRawDigis->end() ) {
1140  throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: No raw DetSet found for det: " << detRawId << "." << std::endl;
1141  }
1142 
1143  // fill the digis to the profiles
1144  edm::DetSet<SiStripRawDigi>::const_iterator digiRangeIterator = detSetIter->data.begin(); // for the loop
1145  edm::DetSet<SiStripRawDigi>::const_iterator digiRangeStart = digiRangeIterator; // save starting positions
1146 
1147  // loop all digis
1148  for (; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
1149  const SiStripRawDigi& digi = *digiRangeIterator;
1150  const int channel = distance( digiRangeStart, digiRangeIterator );
1151  if ( channel >= 0 && channel < 512 ) currentDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetValue( channel, digi.adc() );
1152  else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: raw digi channel: " << channel << " out of range for det: " << detRawId << "." << std::endl;
1153  }
1154 
1155  }
1156 
1157  else { // we have zero suppressed SiStripDigis
1158 
1159  // search the digis for the raw id
1160  edm::DetSetVector<SiStripDigi>::const_iterator detSetIter = theStripDigis->find( detRawId );
1161 
1162  // processed DetSets may be missing, just skip
1163  if( detSetIter == theStripDigis->end() ) continue;
1164 
1165  // fill the digis to the profiles
1166  edm::DetSet<SiStripDigi>::const_iterator digiRangeIterator = detSetIter->data.begin(); // for the loop
1167 
1168  for(; digiRangeIterator != detSetIter->data.end(); ++digiRangeIterator ) {
1169  const SiStripDigi& digi = *digiRangeIterator;
1170  if ( digi.strip() < 512 ) currentDataProfiles.GetTEC2TECEntry( det, beam, disk ).SetValue( digi.strip(), digi.adc() );
1171  else throw cms::Exception( "[Laser Alignment::fillDataProfiles]" ) << " ** ERROR: digi strip: " << digi.strip() << " out of range for det: " << detRawId << "." << std::endl;
1172  }
1173 
1174  }
1175 
1176  } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
1177 
1178  } // theDigiProducersList loop
1179 
1180 }
1181 
1182 
1183 
1184 
1185 
1195 
1196  int det, ring, beam, disk, pos;
1197 
1198  // loop TEC modules (yet without AT)
1199  det = 0; ring = 0; beam = 0; disk = 0;
1200  do { // loop using LASGlobalLoop functionality
1201  SiStripPedestals::Range pedRange = pedestalsHandle->getRange( detectorId.GetTECEntry( det, ring, beam, disk ) );
1202  for( int strip = 0; strip < 512; ++strip ) {
1203  int thePedestal = int( pedestalsHandle->getPed( strip, pedRange ) );
1204  if( thePedestal > 895 ) thePedestal -= 1024;
1205  pedestalProfiles.GetTECEntry( det, ring, beam, disk ).SetValue( strip, thePedestal );
1206  }
1207  } while ( moduleLoop.TECLoop( det, ring, beam, disk ) );
1208 
1209 
1210  // TIB & TOB section
1211  det = 2; beam = 0; pos = 0;
1212  do { // loop using LASGlobalLoop functionality
1213  SiStripPedestals::Range pedRange = pedestalsHandle->getRange( detectorId.GetTIBTOBEntry( det, beam, pos ) );
1214  for( int strip = 0; strip < 512; ++strip ) {
1215  int thePedestal = int( pedestalsHandle->getPed( strip, pedRange ) );
1216  if( thePedestal > 895 ) thePedestal -= 1024;
1217  pedestalProfiles.GetTIBTOBEntry( det, beam, pos ).SetValue( strip, thePedestal );
1218  }
1219  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
1220 
1221 
1222  // TEC2TEC AT section
1223  det = 0; beam = 0; disk = 0;
1224  do { // loop using LASGlobalLoop functionality
1225  SiStripPedestals::Range pedRange = pedestalsHandle->getRange( detectorId.GetTEC2TECEntry( det, beam, disk ) );
1226  for( int strip = 0; strip < 512; ++strip ) {
1227  int thePedestal = int( pedestalsHandle->getPed( strip, pedRange ) );
1228  if( thePedestal > 895 ) thePedestal -= 1024;
1229  pedestalProfiles.GetTEC2TECEntry( det, beam, disk ).SetValue( strip, thePedestal );
1230  }
1231  } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
1232 
1233 }
1234 
1235 
1236 
1237 
1238 
1245 
1246  int numberOfProfiles = 0;
1247 
1248  int ring = 1; // search all ring6 modules for signals
1249  for( int det = 0; det < 2; ++det ) {
1250  for( int beam = 0; beam < 8; ++ beam ) {
1251  for( int disk = 0; disk < 9; ++disk ) {
1252  if( isAcceptedProfile.GetTECEntry( det, ring, beam, disk ) == 1 ) numberOfProfiles++;
1253  }
1254  }
1255  }
1256 
1257  LogDebug( "[LaserAlignment::isTECBeam]" ) << " Found: " << numberOfProfiles << "hits." << std::endl;
1258  std::cout << " [LaserAlignment::isTECBeam] -- Found: " << numberOfProfiles << " hits." << std::endl;
1259 
1260  if( numberOfProfiles > 10 ) return( true );
1261  return( false );
1262 
1263 }
1264 
1265 
1266 
1267 
1268 
1274 
1276 
1277  int numberOfProfiles = 0;
1278 
1279  int det = 2; int beam = 0; int pos = 0; // search all TIB/TOB for signals
1280  do {
1281  if( isAcceptedProfile.GetTIBTOBEntry( det, beam, pos ) == 1 ) numberOfProfiles++;
1282  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
1283 
1284  LogDebug( "[LaserAlignment::isATBeam]" ) << " Found: " << numberOfProfiles << "hits." << std::endl;
1285  std::cout << " [LaserAlignment::isATBeam] -- Found: " << numberOfProfiles << " hits." << std::endl;
1286 
1287  if( numberOfProfiles > 10 ) return( true );
1288  return( false );
1289 
1290 }
1291 
1292 
1293 
1294 
1295 
1304 double LaserAlignment::getTIBTOBNominalBeamOffset( unsigned int det, unsigned int beam, unsigned int pos ) {
1305 
1306  if( det < 2 || det > 3 || beam > 7 || pos > 5 ) {
1307  throw cms::Exception( "[LaserAlignment::getTIBTOBNominalBeamOffset]" ) << " ERROR ** Called with nonexisting parameter set: det " << det << " beam " << beam << " pos " << pos << "." << std::endl;
1308  }
1309 
1310  const double nominalOffsetsTIB[8] = { 0.00035, 2.10687, -2.10827, -0.00173446, 2.10072, -0.00135114, 2.10105, -2.10401 };
1311 
1312  // in tob, modules have alternating orientations along the rods.
1313  // this is described by the following pattern.
1314  // (even more confusing, this pattern is inversed for beams 0, 5, 6, 7)
1315  const int orientationPattern[6] = { -1, 1, 1, -1, -1, 1 };
1316  const double nominalOffsetsTOB[8] = { 0.00217408, 1.58678, 117.733, 119.321, 120.906, 119.328, 117.743, 1.58947 };
1317 
1318 
1319  if( det == 2 ) return( -1. * nominalOffsetsTIB[beam] );
1320 
1321  else {
1322  if( beam == 0 or beam > 4 ) return( nominalOffsetsTOB[beam] * orientationPattern[pos] );
1323  else return( -1. * nominalOffsetsTOB[beam] * orientationPattern[pos] );
1324  }
1325 
1326 }
1327 
1328 
1329 
1330 
1339 double LaserAlignment::getTEC2TECNominalBeamOffset( unsigned int det, unsigned int beam, unsigned int disk ) {
1340 
1341  if( det > 1 || beam > 7 || disk > 5 ) {
1342  throw cms::Exception( "[LaserAlignment::getTEC2TECNominalBeamOffset]" ) << " ERROR ** Called with nonexisting parameter set: det " << det << " beam " << beam << " disk " << disk << "." << std::endl;
1343  }
1344 
1345  const double nominalOffsets[8] = { 0., 2.220, -2.221, 0., 2.214, 0., 2.214, -2.217 };
1346 
1347  if( det == 0 ) return -1. * nominalOffsets[beam];
1348  else return nominalOffsets[beam];
1349 
1350 }
1351 
1352 
1353 
1354 
1355 
1360 
1361  //
1362  // hard coded data yet...
1363  //
1364 
1365  // nominal phi values of tec beam / alignment tube hits (parameter is beam 0-7)
1366  const double tecPhiPositions[8] = { 0.392699, 1.178097, 1.963495, 2.748894, 3.534292, 4.319690, 5.105088, 5.890486 }; // new values calculated by maple
1367  const double atPhiPositions[8] = { 0.392699, 1.289799, 1.851794, 2.748894, 3.645995, 4.319690, 5.216791, 5.778784 }; // new values calculated by maple
1368 
1369  // nominal r values (mm) of hits
1370  const double tobRPosition = 600.;
1371  const double tibRPosition = 514.;
1372  const double tecRPosition[2] = { 564., 840. }; // ring 4,6
1373 
1374  // nominal z values (mm) of hits in barrel (parameter is pos 0-6)
1375  const double tobZPosition[6] = { 1040., 580., 220., -140., -500., -860. };
1376  const double tibZPosition[6] = { 620., 380., 180., -100., -340., -540. };
1377 
1378  // nominal z values (mm) of hits in tec (parameter is disk 0-8); FOR TEC-: (* -1.)
1379  const double tecZPosition[9] = { 1322.5, 1462.5, 1602.5, 1742.5, 1882.5, 2057.5, 2247.5, 2452.5, 2667.5 };
1380 
1381 
1382  //
1383  // now we fill these into the nominalCoordinates container;
1384  // errors are zero for nominal values..
1385  //
1386 
1387  // loop object and its variables
1389  int det, ring, beam, disk, pos;
1390 
1391 
1392  // TEC+- section
1393  det = 0; ring = 0, beam = 0; disk = 0;
1394  do {
1395 
1396  if( det == 0 ) { // this is TEC+
1397  nominalCoordinates.SetTECEntry( det, ring, beam, disk, LASCoordinateSet( tecPhiPositions[beam], 0., tecRPosition[ring], 0., tecZPosition[disk], 0. ) );
1398  }
1399  else { // now TEC-
1400  nominalCoordinates.SetTECEntry( det, ring, beam, disk, LASCoordinateSet( tecPhiPositions[beam], 0., tecRPosition[ring], 0., -1. * tecZPosition[disk], 0. ) ); // just * -1.
1401  }
1402 
1403  } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
1404 
1405 
1406 
1407  // TIB & TOB section
1408  det = 2; beam = 0; pos = 0;
1409  do {
1410  if( det == 2 ) { // this is TIB
1411  nominalCoordinates.SetTIBTOBEntry( det, beam, pos, LASCoordinateSet( atPhiPositions[beam], 0., tibRPosition, 0., tibZPosition[pos], 0. ) );
1412  }
1413  else { // now TOB
1414  nominalCoordinates.SetTIBTOBEntry( det, beam, pos, LASCoordinateSet( atPhiPositions[beam], 0., tobRPosition, 0., tobZPosition[pos], 0. ) );
1415  }
1416 
1417  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
1418 
1419 
1420 
1421 
1422  // TEC2TEC AT section
1423  det = 0; beam = 0; disk = 0;
1424  do {
1425 
1426  if( det == 0 ) { // this is TEC+, ring4 only
1427  nominalCoordinates.SetTEC2TECEntry( det, beam, disk, LASCoordinateSet( atPhiPositions[beam], 0., tecRPosition[0], 0., tecZPosition[disk], 0. ) );
1428  }
1429  else { // now TEC-
1430  nominalCoordinates.SetTEC2TECEntry( det, beam, disk, LASCoordinateSet( atPhiPositions[beam], 0., tecRPosition[0], 0., -1. * tecZPosition[disk], 0. ) ); // just * -1.
1431  }
1432 
1433  } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
1434 
1435 
1436 }
1437 
1438 
1439 
1440 
1441 
1447 
1448  if( angle < -1. * M_PI || angle > M_PI ) {
1449  throw cms::Exception(" [LaserAlignment::ConvertAngle] ") << "** ERROR: Called with illegal input angle: " << angle << "." << std::endl;
1450  }
1451 
1452  if( angle >= 0. ) return angle;
1453  else return( angle + 2. * M_PI );
1454 
1455 }
1456 
1457 
1458 
1459 
1460 
1465 
1467  int det, ring, beam, disk, pos;
1468 
1469  std:: cout << std::endl << " [LaserAlignment::DumpPosFileSet] -- Dump: " << std::endl;
1470 
1471  // TEC INTERNAL
1472  det = 0; ring = 0; beam = 0; disk = 0;
1473  do {
1474  std::cout << "POS " << det << "\t" << beam << "\t" << disk << "\t" << ring << "\t" << coordinates.GetTECEntry( det, ring, beam, disk ).GetPhi() << "\t" << coordinates.GetTECEntry( det, ring, beam, disk ).GetPhiError() << std::endl;
1475  } while ( loop.TECLoop( det, ring, beam, disk ) );
1476 
1477  // TIBTOB
1478  det = 2; beam = 0; pos = 0;
1479  do {
1480  std::cout << "POS " << det << "\t" << beam << "\t" << pos << "\t" << "-1" << "\t" << coordinates.GetTIBTOBEntry( det, beam, pos ).GetPhi() << "\t" << coordinates.GetTIBTOBEntry( det, beam, pos ).GetPhiError() << std::endl;
1481  } while( loop.TIBTOBLoop( det, beam, pos ) );
1482 
1483  // TEC2TEC
1484  det = 0; beam = 0; disk = 0;
1485  do {
1486  std::cout << "POS " << det << "\t" << beam << "\t" << disk << "\t" << "-1" << "\t" << coordinates.GetTEC2TECEntry( det, beam, disk ).GetPhi() << "\t" << coordinates.GetTEC2TECEntry( det, beam, disk ).GetPhiError() << std::endl;
1487  } while( loop.TEC2TECLoop( det, beam, disk ) );
1488 
1489  std:: cout << std::endl << " [LaserAlignment::DumpPosFileSet] -- End dump: " << std::endl;
1490 
1491 }
1492 
1493 
1494 
1495 
1496 
1500 void LaserAlignment::DumpStripFileSet( LASGlobalData<std::pair<float,float> >& measuredStripPositions ) {
1501 
1503  int det, ring, beam, disk, pos;
1504 
1505  std:: cout << std::endl << " [LaserAlignment::DumpStripFileSet] -- Dump: " << std::endl;
1506 
1507  // TEC INTERNAL
1508  det = 0; ring = 0; beam = 0; disk = 0;
1509  do {
1510  std::cout << "STRIP " << det << "\t" << beam << "\t" << disk << "\t" << ring << "\t" << measuredStripPositions.GetTECEntry( det, ring, beam, disk ).first
1511  << "\t" << measuredStripPositions.GetTECEntry( det, ring, beam, disk ).second << std::endl;
1512  } while ( loop.TECLoop( det, ring, beam, disk ) );
1513 
1514  // TIBTOB
1515  det = 2; beam = 0; pos = 0;
1516  do {
1517  std::cout << "STRIP " << det << "\t" << beam << "\t" << pos << "\t" << "-1" << "\t" << measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).first
1518  << "\t" << measuredStripPositions.GetTIBTOBEntry( det, beam, pos ).second << std::endl;
1519  } while( loop.TIBTOBLoop( det, beam, pos ) );
1520 
1521  // TEC2TEC
1522  det = 0; beam = 0; disk = 0;
1523  do {
1524  std::cout << "STRIP " << det << "\t" << beam << "\t" << disk << "\t" << "-1" << "\t" << measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).first
1525  << "\t" << measuredStripPositions.GetTEC2TECEntry( det, beam, disk ).second << std::endl;
1526  } while( loop.TEC2TECLoop( det, beam, disk ) );
1527 
1528  std:: cout << std::endl << " [LaserAlignment::DumpStripFileSet] -- End dump: " << std::endl;
1529 
1530 
1531 }
1532 
1533 
1534 
1535 
1536 
1541 
1542  std::cout << " [LaserAlignment::DumpHitmaps] -- Dumping hitmap for TEC+:" << std::endl;
1543  std::cout << " [LaserAlignment::DumpHitmaps] -- Ring4:" << std::endl;
1544  std::cout << " disk0 disk1 disk2 disk3 disk4 disk5 disk6 disk7 disk8" << std::endl;
1545 
1546  for( int beam = 0; beam < 8; ++beam ) {
1547  std::cout << " beam" << beam << ":";
1548  for( int disk = 0; disk < 9; ++disk ) {
1549  std::cout << "\t" << numberOfAcceptedProfiles.GetTECEntry( 0, 0, beam, disk );
1550  }
1551  std::cout << std::endl;
1552  }
1553 
1554  std::cout << " [LaserAlignment::DumpHitmaps] -- Ring6:" << std::endl;
1555  std::cout << " disk0 disk1 disk2 disk3 disk4 disk5 disk6 disk7 disk8" << std::endl;
1556 
1557  for( int beam = 0; beam < 8; ++beam ) {
1558  std::cout << " beam" << beam << ":";
1559  for( int disk = 0; disk < 9; ++disk ) {
1560  std::cout << "\t" << numberOfAcceptedProfiles.GetTECEntry( 0, 1, beam, disk );
1561  }
1562  std::cout << std::endl;
1563  }
1564 
1565  std::cout << " [LaserAlignment::DumpHitmaps] -- Dumping hitmap for TEC-:" << std::endl;
1566  std::cout << " [LaserAlignment::DumpHitmaps] -- Ring4:" << std::endl;
1567  std::cout << " disk0 disk1 disk2 disk3 disk4 disk5 disk6 disk7 disk8" << std::endl;
1568 
1569  for( int beam = 0; beam < 8; ++beam ) {
1570  std::cout << " beam" << beam << ":";
1571  for( int disk = 0; disk < 9; ++disk ) {
1572  std::cout << "\t" << numberOfAcceptedProfiles.GetTECEntry( 1, 0, beam, disk );
1573  }
1574  std::cout << std::endl;
1575  }
1576 
1577  std::cout << " [LaserAlignment::DumpHitmaps] -- Ring6:" << std::endl;
1578  std::cout << " disk0 disk1 disk2 disk3 disk4 disk5 disk6 disk7 disk8" << std::endl;
1579 
1580  for( int beam = 0; beam < 8; ++beam ) {
1581  std::cout << " beam" << beam << ":";
1582  for( int disk = 0; disk < 9; ++disk ) {
1583  std::cout << "\t" << numberOfAcceptedProfiles.GetTECEntry( 1, 1, beam, disk );
1584  }
1585  std::cout << std::endl;
1586  }
1587 
1588  std::cout << " [LaserAlignment::DumpHitmaps] -- End of dump." << std::endl << std::endl;
1589 
1590 }
1591 
1592 
1593 
1594 
1595 
1601 
1602  // loop the list of modules to be masked
1603  for( std::vector<unsigned int>::iterator moduleIt = theMaskTecModules.begin(); moduleIt != theMaskTecModules.end(); ++moduleIt ) {
1604 
1605  // loop variables
1607  int det, ring, beam, disk;
1608 
1609  // this will calculate the corrections from the alignment parameters
1610  LASEndcapAlgorithm endcapAlgorithm;
1611 
1612  // find the location of the respective module in the container with this loop
1613  det = 0; ring = 0; beam = 0; disk = 0;
1614  do {
1615 
1616  // here we got it
1617  if( detectorId.GetTECEntry( det, ring, beam, disk ) == *moduleIt ) {
1618 
1619  // the nominal phi value for this module
1620  const double nominalPhi = nominalCoordinates.GetTECEntry( det, ring, beam, disk ).GetPhi();
1621 
1622  // the offset from the alignment parameters
1623  const double phiCorrection = endcapAlgorithm.GetAlignmentParameterCorrection( det, ring, beam, disk, nominalCoordinates, endcapParameters );
1624 
1625  // apply the corrections
1626  measuredCoordinates.GetTECEntry( det, ring, beam, disk ).SetPhi( nominalPhi - phiCorrection );
1627 
1628  }
1629 
1630  } while ( moduleLoop.TECLoop( det, ring, beam, disk ) );
1631 
1632  }
1633 
1634 }
1635 
1636 
1637 
1638 
1639 
1645 
1646  // loop the list of modules to be masked
1647  for( std::vector<unsigned int>::iterator moduleIt = theMaskAtModules.begin(); moduleIt != theMaskAtModules.end(); ++moduleIt ) {
1648 
1649  // loop variables
1651  int det, beam, disk, pos;
1652 
1653  // this will calculate the corrections from the alignment parameters
1654  LASAlignmentTubeAlgorithm atAlgorithm;
1655 
1656 
1657  // find the location of the respective module in the container with these loops:
1658 
1659  // first TIB+TOB
1660  det = 2; beam = 0; pos = 0;
1661  do {
1662 
1663  // here we got it
1664  if( detectorId.GetTIBTOBEntry( det, beam, pos ) == *moduleIt ) {
1665 
1666  // the nominal phi value for this module
1667  const double nominalPhi = nominalCoordinates.GetTIBTOBEntry( det, beam, pos ).GetPhi();
1668 
1669  // the offset from the alignment parameters
1670  const double phiCorrection = atAlgorithm.GetTIBTOBAlignmentParameterCorrection( det, beam, pos, nominalCoordinates, atParameters );
1671 
1672  // apply the corrections
1673  measuredCoordinates.GetTIBTOBEntry( det, beam, pos ).SetPhi( nominalPhi - phiCorrection );
1674 
1675  }
1676 
1677  } while ( moduleLoop.TIBTOBLoop( det, beam, pos ) );
1678 
1679 
1680 
1681  // then TEC(AT)
1682  det = 0; beam = 0; disk = 0;
1683  do {
1684 
1685  // here we got it
1686  if( detectorId.GetTEC2TECEntry( det, beam, disk ) == *moduleIt ) {
1687 
1688  // the nominal phi value for this module
1689  const double nominalPhi = nominalCoordinates.GetTEC2TECEntry( det, beam, disk ).GetPhi();
1690 
1691  // the offset from the alignment parameters
1692  const double phiCorrection = atAlgorithm.GetTEC2TECAlignmentParameterCorrection( det, beam, disk, nominalCoordinates, atParameters );
1693 
1694  // apply the corrections
1695  measuredCoordinates.GetTEC2TECEntry( det, beam, disk ).SetPhi( nominalPhi - phiCorrection );
1696 
1697  }
1698 
1699  } while ( moduleLoop.TEC2TECLoop( det, beam, disk ) );
1700 
1701  }
1702 
1703 }
1704 
1705 
1706 
1707 
1708 
1714 
1715 
1716  // tracker geom. object for calculating the global beam positions
1717  const TrackerGeometry& theTracker( *theTrackerGeometry );
1718 
1719  const double atPhiPositions[8] = { 0.392699, 1.289799, 1.851794, 2.748894, 3.645995, 4.319690, 5.216791, 5.778784 };
1720  const double tecPhiPositions[8] = { 0.392699, 1.178097, 1.963495, 2.748894, 3.534292, 4.319690, 5.105088, 5.890486 };
1721  const double zPositions[9] = { 125.0, 139.0, 153.0, 167.0, 181.0, 198.5, 217.5, 238.0, 259.5 };
1722  const double zPositionsTIB[6] = { 62.0, 38.0, 18.0, -10.0, -34.0, -54.0 };
1723  const double zPositionsTOB[6] = { 104.0, 58.0, 22.0, -14.0, -50.0, -86.0 };
1724 
1725  int det, beam, disk, pos, ring;
1726 
1727  // loop TEC+- internal
1728  det = 0; ring = 0; beam = 0; disk = 0;
1729  do {
1730 
1731  const double radius = ring?84.0:56.4;
1732 
1733  // access the tracker geometry for this module
1734  const DetId theDetId( detectorId.GetTECEntry( det, ring, beam, disk ) );
1735  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
1736 
1737  if (theStripDet) {
1738  const GlobalPoint gp( GlobalPoint::Cylindrical( radius, tecPhiPositions[beam], zPositions[disk] ) );
1739 
1740  const LocalPoint lp( theStripDet->surface().toLocal( gp ) );
1741  std::cout << "__TEC: " << 256. - theStripDet->specificTopology().strip( lp ) << std::endl;
1742  }
1743 
1744  } while( moduleLoop.TECLoop( det, ring, beam, disk ) );
1745 
1746 
1747  // loop TIBTOB
1748  det = 2; beam = 0; pos = 0;
1749  do {
1750 
1751  const double radius = (det==2?51.4:58.4);
1752  const double theZ = (det==2?zPositionsTIB[pos]:zPositionsTOB[pos]);
1753 
1754  // access the tracker geometry for this module
1755  const DetId theDetId( detectorId.GetTIBTOBEntry( det, beam, pos ) );
1756  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
1757 
1758  if (theStripDet) {
1759  const GlobalPoint gp( GlobalPoint::Cylindrical( radius, atPhiPositions[beam], theZ ) );
1760 
1761  const LocalPoint lp( theStripDet->surface().toLocal( gp ) );
1762  std::cout << "__TIBTOB det " << det << " beam " << beam << " pos " << pos << " " << 256. - theStripDet->specificTopology().strip( lp );
1763  std::cout << " " << theStripDet->position().perp()<< std::endl;
1764  }
1765 
1766  } while( moduleLoop.TIBTOBLoop( det, beam, pos ) );
1767 
1768 
1769  // loop TEC2TEC
1770  det = 0; beam = 0; disk = 0;
1771  do {
1772 
1773  const double radius = 56.4;
1774 
1775  // access the tracker geometry for this module
1776  const DetId theDetId( detectorId.GetTEC2TECEntry( det, beam, disk ) );
1777  const StripGeomDetUnit* const theStripDet = dynamic_cast<const StripGeomDetUnit*>( theTracker.idToDet( theDetId ) );
1778 
1779  if (theStripDet) {
1780  const GlobalPoint gp( GlobalPoint::Cylindrical( radius, atPhiPositions[beam], zPositions[disk] ) );
1781 
1782  const LocalPoint lp( theStripDet->surface().toLocal( gp ) );
1783  std::cout << "__TEC2TEC det " << det << " beam " << beam << " disk " << disk << " " << 256. - theStripDet->specificTopology().strip( lp ) << std::endl;
1784  }
1785 
1786  } while( moduleLoop.TEC2TECLoop( det, beam, disk ) );
1787 
1788 
1789 }
1790 
1791 
1792 
1793 
1794 
1795 
1796 
1797 
1798 
1799 // define the SEAL module
1801 
1803 
1804 
1805 
1806 
1807 // the ATTIC
1808 
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
bool IsSignalIn(const LASModuleProfile &, double)
bool isATBeam(void)
void SetTEC2TECEntry(int subdetector, int beam, int tecDisk, T)
T theZ
void SetPhi(double aPhi)
const uint16_t & adc() const
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
unsigned int judgeOverdriveThreshold
config parameters for the LASProfileJudge
void fillDataProfiles(edm::Event const &, edm::EventSetup const &)
fill profiles from SiStrip(Raw)Digi container
void beginJob() override
bool theSaveHistograms
config switch
void SetPhiError(double aPhiError)
void ApplyEndcapMaskingCorrections(LASGlobalData< LASCoordinateSet > &, LASGlobalData< LASCoordinateSet > &, LASEndcapAlignmentParameterSet &)
apply endcap correction to masked modules in TEC
void produce(edm::Event &, edm::EventSetup const &) override
void TrackerUpdate(LASEndcapAlignmentParameterSet &, LASBarrelAlignmentParameterSet &, AlignableTracker &)
void SetMisalignmentFromRefGeometry(bool)
bool misalignedByRefGeometry
config switch
std::string theAlignRecordName
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool isTECBeam(void)
decide whether TEC or AT beams have fired
LaserAlignment(edm::ParameterSet const &theConf)
double getTIBTOBNominalBeamOffset(unsigned int, unsigned int, unsigned int)
returns the nominal beam position (strips) in TOB for the profileJudge
edm::ESHandle< TrackerGeometry > theTrackerGeometry
double peakFinderThreshold
config parameter
void DumpHitmaps(LASGlobalData< int > &)
for debugging only, will disappear
void DumpPosFileSet(LASGlobalData< LASCoordinateSet > &)
for debugging only, will disappear
LASGlobalData< int > numberOfAcceptedProfiles
void SetAllValuesTo(const double &)
TrackerGeometry * build(const GeometricDet *gd, const PTrackerParameters &ptp, const TrackerTopology *tTopo)
std::pair< ContainerIterator, ContainerIterator > Range
bool theDoPedestalSubtraction
config switch
LASGlobalData< LASModuleProfile > currentDataProfiles
data profiles for the current event
std::vector< unsigned int > theMaskTecModules
config parameters
bool theUseMinuitAlgorithm
config switch
TFile * theFile
Tree stuff.
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
LASGlobalData< int > isAcceptedProfile
std::vector< edm::ParameterSet > theDigiProducersList
void DumpStripFileSet(LASGlobalData< std::pair< float, float > > &)
for debugging only, will disappear
void fillPedestalProfiles(edm::ESHandle< SiStripPedestals > &)
fill pedestals from dbase
double GetPhi(void) const
~LaserAlignment() override
float getPed(const uint16_t &strip, const Range &range) const
Handle< TkLasBeamCollection > laserBeams
T barePhi() const
Definition: PV3DBase.h:68
LASGlobalData< TH1D * > summedHistograms
bool FindPeakIn(const LASModuleProfile &, std::pair< double, double > &, TH1D *, const double)
TDirectory * singleModulesDir
void push_back(const SiStripLaserRecHit2D &aHit)
insert a hit in the data vector
Definition: TkLasBeam.h:37
double ConvertAngle(double)
convert an angle in the [-pi,pi] range to the [0,2*pi] range
LASBarrelAlignmentParameterSet CalculateParameters(LASGlobalData< LASCoordinateSet > &, LASGlobalData< LASCoordinateSet > &)
edm::ESHandle< Alignments > theGlobalPositionRcd
const uint16_t & strip() const
Definition: SiStripDigi.h:40
bool TEC2TECLoop(int &, int &, int &) const
std::string theFileName
config parameter (histograms file output name)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
void get(HolderT &iHolder) const
bool isAvailable() const
Definition: Service.h:46
void endJob() override
std::vector< unsigned int > tecDoubleHitDetId
void writeOne(T *payload, Time_t time, const std::string &recordName, bool withlogging=false)
void SetAmplitudeThreshold(double)
void endRunProduce(edm::Run &, const edm::EventSetup &) override
bool theApplyBeamKinkCorrections
config switch
T & GetTIBTOBEntry(int subdetector, int beam, int tibTobPosition)
bool JudgeProfile(const LASModuleProfile &, double)
A Digi for the silicon strip detector, containing both strip and adc information, and suitable for st...
Definition: SiStripDigi.h:12
void SetTIBTOBEntry(int subdetector, int beam, int tibTobPosition, T)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:475
void CalculateNominalCoordinates(void)
fills a LASGlobalData<LASCoordinateSet> with nominal module positions
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:361
LASEndcapAlignmentParameterSet CalculateParameters(LASGlobalData< LASCoordinateSet > &, LASGlobalData< LASCoordinateSet > &)
std::string theErrorRecordName
#define M_PI
double GetTEC2TECAlignmentParameterCorrection(int, int, int, LASGlobalData< LASCoordinateSet > &, LASBarrelAlignmentParameterSet &)
void SetValue(unsigned int theStripNumber, const double &theValue)
bool theStoreToDB
config switch
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
Definition: DetId.h:18
LASGlobalData< LASModuleProfile > pedestalProfiles
T & GetTEC2TECEntry(int subdetector, int beam, int tecDisk)
LASGlobalData< LASCoordinateSet > nominalCoordinates
AlignableTracker * theAlignableTracker
int theEvents
counter for the total number of events processed
const T & get() const
Definition: EventSetup.h:58
LASBarrelAlignmentParameterSet CalculateParameters(LASGlobalData< LASCoordinateSet > &, LASGlobalData< LASCoordinateSet > &)
LASConstants theLasConstants
double getTEC2TECNominalBeamOffset(unsigned int, unsigned int, unsigned int)
returns the nominal beam position (strips) in TEC (AT) for the profileJudge
void put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Run.h:114
void SetOverdriveThreshold(unsigned int)
T & GetTECEntry(int subdetector, int tecRing, int beam, int tecDisk)
Definition: LASGlobalData.h:91
LASGlobalData< std::string > theProfileNames
double GetAlignmentParameterCorrection(int, int, int, int, LASGlobalData< LASCoordinateSet > &, LASEndcapAlignmentParameterSet &)
std::vector< unsigned int > theMaskAtModules
void ApplyATMaskingCorrections(LASGlobalData< LASCoordinateSet > &, LASGlobalData< LASCoordinateSet > &, LASBarrelAlignmentParameterSet &)
same for alignment tube modules
void testRoutine(void)
for debugging & testing only, will disappear..
LASGlobalData< LASModuleProfile > collectedDataProfiles
edm::EventID id() const
Definition: EventBase.h:60
void fillDetectorId(void)
fill hard coded detIds
HLT enums.
Alignments * alignments() const override
Return alignments, sorted by DetId.
void ApplyBeamKinkCorrections(LASGlobalData< LASCoordinateSet > &) const
int theCompression
config parameter (histograms file compression level)
AlignmentErrorsExtended * alignmentErrors() const override
Return alignment errors, sorted by DetId.
void EnableZeroFilter(bool)
const TrackerGeomDet * idToDet(DetId) const override
LASGlobalData< unsigned int > detectorId
bool TECLoop(int &, int &, int &, int &) const
bool TIBTOBLoop(int &, int &, int &) const
double GetPhiError(void) const
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:346
const Range getRange(const uint32_t &detID) const
bool enableJudgeZeroFilter
config switch
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:104
const uint16_t & adc() const
Definition: SiStripDigi.h:41
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
T const * product() const
Definition: ESHandle.h:86
double GetTIBTOBAlignmentParameterCorrection(int, int, int, LASGlobalData< LASCoordinateSet > &, LASBarrelAlignmentParameterSet &)
void SetTECEntry(int subdetector, int tecRing, int beam, int tecDisk, T)
edm::ESHandle< GeometricDet > gD
tracker geometry;
def move(src, dest)
Definition: eostools.py:510
void EndcapUpdate(LASEndcapAlignmentParameterSet &, LASGlobalData< LASCoordinateSet > &)
Definition: Run.h:43
LASGlobalLoop moduleLoop
LASProfileJudge judge
bool updateFromInputGeometry
config switch
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
bool theSetNominalStrips
config switch