CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/Alignment/SurveyAnalysis/plugins/SurveyMisalignmentInput.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/EventSetup.h"
00002 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00003 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeomBuilderFromGeometricDet.h"
00004 
00005 #include "Alignment/CommonAlignment/interface/SurveyDet.h"
00006 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
00007 
00008 #include "CondFormats/Alignment/interface/Alignments.h"
00009 #include "CondFormats/AlignmentRecord/interface/TrackerAlignmentRcd.h"
00010 
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 
00014 #include "Alignment/SurveyAnalysis/plugins/SurveyMisalignmentInput.h"
00015 
00016 SurveyMisalignmentInput::SurveyMisalignmentInput(const edm::ParameterSet& cfg):
00017   textFileName( cfg.getParameter<std::string>("textFileName") ),
00018   theParameterSet( cfg )
00019 {}
00020 
00021 void SurveyMisalignmentInput::analyze(const edm::Event&, const edm::EventSetup& setup)
00022 {
00023   if (theFirstEvent) {
00024     edm::ESHandle<GeometricDet> geom;
00025     setup.get<IdealGeometryRecord>().get(geom);  
00026     TrackerGeometry* tracker = TrackerGeomBuilderFromGeometricDet().build(&*geom, theParameterSet);
00027     
00028     addComponent(new AlignableTracker( tracker ) );
00029 
00030     edm::LogInfo("SurveyMisalignmentInput") << "Starting!";
00031     // Retrieve alignment[Error]s from DBase
00032     setup.get<TrackerAlignmentRcd>().get( alignments );
00033     
00034     //Get map from textreader
00035     SurveyInputTextReader dataReader;
00036     dataReader.readFile( textFileName );
00037     uIdMap = dataReader.UniqueIdMap();
00038     
00039     addSurveyInfo( detector() );
00040 
00041     theFirstEvent = false;
00042   }
00043 }
00044 
00045 
00046 void SurveyMisalignmentInput::addSurveyInfo(Alignable* ali)
00047 {
00048 
00049   const align::Alignables& comp = ali->components();
00050   unsigned int nComp = comp.size();
00051   for (unsigned int i = 0; i < nComp; ++i) addSurveyInfo(comp[i]);
00052         
00053   SurveyInputTextReader::MapType::const_iterator it
00054     = uIdMap.find(std::make_pair(ali->id(), ali->alignableObjectId()));
00055 
00056   align::ErrorMatrix error;
00057 
00058   if (it != uIdMap.end()){
00059     //survey error values
00060     const align::Scalars& parameters = (it)->second;
00061     //sets the errors for the hierarchy level
00062     double* errorData = error.Array();
00063     for (unsigned int i = 0; i < 21; ++i){errorData[i] = parameters[i+6];}
00064                 
00065     //because record only needs global value of modules
00066     if (ali->alignableObjectId() == align::AlignableDetUnit){
00067       // fill survey values
00068       ali->setSurvey( new SurveyDet(getAlignableSurface(ali->id()), error) );
00069     }
00070     else{
00071       ali->setSurvey( new SurveyDet(ali->surface(), error) );
00072     }
00073   }
00074   else{
00075     //fill
00076     error = ROOT::Math::SMatrixIdentity();
00077     ali->setSurvey( new SurveyDet(ali->surface(), error*(1e-6)) );
00078   }
00079   //std::cout << "UniqueId: " << id.first << ", " << id.second << std::endl;
00080   //std::cout << error << std::endl;
00081         
00082 }
00083 
00084 AlignableSurface SurveyMisalignmentInput::getAlignableSurface(align::ID id)
00085 {
00086   std::vector<AlignTransform>::const_iterator it;
00087 
00088   for (it = alignments->m_align.begin(); it != alignments->m_align.end(); ++it)
00089   {
00090     if (id == (*it).rawId())
00091     {
00092       align::PositionType position( (*it).translation().x(), (*it).translation().y(), (*it).translation().z() );
00093       CLHEP::HepRotation rot( (*it).rotation() );
00094       align::RotationType rotation( rot.xx(), rot.xy(), rot.xz(),
00095                                     rot.yx(), rot.yy(), rot.yz(),
00096                                     rot.zx(), rot.zy(), rot.zz() );
00097       return AlignableSurface(position,rotation);
00098     }
00099   }
00100         
00101   return AlignableSurface();
00102 }
00103 
00104 // Plug in to framework
00105 
00106 #include "FWCore/Framework/interface/MakerMacros.h"
00107 
00108 DEFINE_FWK_MODULE(SurveyMisalignmentInput);