Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-22 04:02:34

0001 #include "FWCore/Framework/interface/EventSetup.h"
0002 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0003 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0004 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeomBuilderFromGeometricDet.h"
0005 #include "CondFormats/GeometryObjects/interface/PTrackerParameters.h"
0006 #include "Geometry/Records/interface/PTrackerParametersRcd.h"
0007 #include "Geometry/Records/interface/PTrackerAdditionalParametersPerDetRcd.h"
0008 
0009 #include "Alignment/CommonAlignment/interface/SurveyDet.h"
0010 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
0011 
0012 #include "CondFormats/Alignment/interface/Alignments.h"
0013 #include "CondFormats/AlignmentRecord/interface/TrackerAlignmentRcd.h"
0014 
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 
0018 #include "Alignment/SurveyAnalysis/plugins/SurveyMisalignmentInput.h"
0019 
0020 SurveyMisalignmentInput::SurveyMisalignmentInput(const edm::ParameterSet& cfg)
0021     : tTopoToken_(esConsumes()),
0022       geomDetToken_(esConsumes()),
0023       ptpToken_(esConsumes()),
0024       ptitpToken_(esConsumes()),
0025       aliToken_(esConsumes()),
0026       textFileName(cfg.getParameter<std::string>("textFileName")) {}
0027 
0028 void SurveyMisalignmentInput::analyze(const edm::Event&, const edm::EventSetup& setup) {
0029   if (theFirstEvent) {
0030     //Retrieve tracker topology from geometry
0031     const TrackerTopology* const tTopo = &setup.getData(tTopoToken_);
0032     const GeometricDet* geom = &setup.getData(geomDetToken_);
0033     const PTrackerParameters& ptp = setup.getData(ptpToken_);
0034     const PTrackerAdditionalParametersPerDet* ptitp = &setup.getData(ptitpToken_);
0035     TrackerGeometry* tracker = TrackerGeomBuilderFromGeometricDet().build(geom, ptitp, ptp, tTopo);
0036 
0037     addComponent(new AlignableTracker(tracker, tTopo));
0038 
0039     edm::LogInfo("SurveyMisalignmentInput") << "Starting!";
0040     // Retrieve alignment[Error]s from DBase
0041     alignments = setup.getHandle(aliToken_);
0042 
0043     //Get map from textreader
0044     SurveyInputTextReader dataReader;
0045     dataReader.readFile(textFileName);
0046     uIdMap = dataReader.UniqueIdMap();
0047 
0048     addSurveyInfo(detector());
0049 
0050     theFirstEvent = false;
0051   }
0052 }
0053 
0054 void SurveyMisalignmentInput::addSurveyInfo(Alignable* ali) {
0055   const align::Alignables& comp = ali->components();
0056   unsigned int nComp = comp.size();
0057   for (unsigned int i = 0; i < nComp; ++i)
0058     addSurveyInfo(comp[i]);
0059 
0060   SurveyInputTextReader::MapType::const_iterator it = uIdMap.find(std::make_pair(ali->id(), ali->alignableObjectId()));
0061 
0062   align::ErrorMatrix error;
0063 
0064   if (it != uIdMap.end()) {
0065     //survey error values
0066     const align::Scalars& parameters = (it)->second;
0067     //sets the errors for the hierarchy level
0068     double* errorData = error.Array();
0069     for (unsigned int i = 0; i < 21; ++i) {
0070       errorData[i] = parameters[i + 6];
0071     }
0072 
0073     //because record only needs global value of modules
0074     if (ali->alignableObjectId() == align::AlignableDetUnit) {
0075       // fill survey values
0076       ali->setSurvey(new SurveyDet(getAlignableSurface(ali->id()), error));
0077     } else {
0078       ali->setSurvey(new SurveyDet(ali->surface(), error));
0079     }
0080   } else {
0081     //fill
0082     error = ROOT::Math::SMatrixIdentity();
0083     ali->setSurvey(new SurveyDet(ali->surface(), error * (1e-6)));
0084   }
0085   //std::cout << "UniqueId: " << id.first << ", " << id.second << std::endl;
0086   //std::cout << error << std::endl;
0087 }
0088 
0089 AlignableSurface SurveyMisalignmentInput::getAlignableSurface(align::ID id) {
0090   std::vector<AlignTransform>::const_iterator it;
0091 
0092   for (it = alignments->m_align.begin(); it != alignments->m_align.end(); ++it) {
0093     if (id == (*it).rawId()) {
0094       align::PositionType position((*it).translation().x(), (*it).translation().y(), (*it).translation().z());
0095       CLHEP::HepRotation rot((*it).rotation());
0096       align::RotationType rotation(
0097           rot.xx(), rot.xy(), rot.xz(), rot.yx(), rot.yy(), rot.yz(), rot.zx(), rot.zy(), rot.zz());
0098       return AlignableSurface(position, rotation);
0099     }
0100   }
0101 
0102   return AlignableSurface();
0103 }
0104 
0105 // Plug in to framework
0106 
0107 #include "FWCore/Framework/interface/MakerMacros.h"
0108 
0109 DEFINE_FWK_MODULE(SurveyMisalignmentInput);