CMS 3D CMS Logo

VisDTDigi.cc

Go to the documentation of this file.
00001 #include "VisReco/Analyzer/interface/VisDTDigi.h"
00002 #include "VisReco/Analyzer/interface/IguanaService.h"
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "FWCore/Framework/interface/EventSetup.h"
00006 #include "FWCore/Framework/interface/MakerMacros.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 #include "FWCore/ServiceRegistry/interface/Service.h"
00009 #include "FWCore/Utilities/interface/Exception.h"
00010 #include "Iguana/Framework/interface/IgCollection.h"
00011 
00012 #include "DataFormats/DTDigi/interface/DTDigi.h"
00013 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
00014 #include "Geometry/DTGeometry/interface/DTLayer.h"
00015 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00016 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00017 #include "VisReco/VisRecoGeometry/interface/VisTrackingGeometry.h"
00018 
00019 using namespace edm::service;
00020 
00021 VisDTDigi::VisDTDigi(const edm::ParameterSet& iConfig)
00022     : inputTag_(iConfig.getParameter<edm::InputTag>("visDTDigiTag"))
00023 {}
00024 
00025 void VisDTDigi::analyze (const edm::Event& event, const edm::EventSetup& eventSetup)
00026 {
00027     edm::Service<IguanaService> config;
00028     if (! config.isAvailable ()) 
00029     {
00030         throw cms::Exception ("Configuration")
00031             << "VisDTDigi requires the IguanaService\n"
00032             "which is not present in the configuration file.\n"
00033             "You must add the service in the configuration file\n"
00034             "or remove the module that requires it";
00035     }
00036 
00037     edm::ESHandle<DTGeometry> geom;
00038     eventSetup.get<MuonGeometryRecord>().get(geom);
00039 
00040     edm::Handle<DTDigiCollection> collection;
00041     event.getByLabel (inputTag_, collection);
00042 
00043     if ( collection.isValid () && geom.isValid ())
00044     {
00045         IgDataStorage* storage = config->storage();
00046         IgCollection& digis = storage->getCollection("DTDigis_V1");
00047 
00048         IgProperty WIREN = digis.addProperty("wireNumber", int(0));
00049         IgProperty LAYER_ID = digis.addProperty("layerId", int(0));
00050         IgProperty SUPERLAYER_ID = digis.addProperty("superLayerId", int(0));
00051         IgProperty SECTOR_ID = digis.addProperty("sectorId", int(0));
00052         IgProperty STATION_ID = digis.addProperty("stationId", int(0));
00053         IgProperty WHEEL_ID = digis.addProperty("wheelId", int(0));
00054 
00055         IgProperty POS = digis.addProperty("pos", IgV3d());
00056         //IgProperty AXIS = digis.addProperty("axis", IgV3d());
00057         //IgProperty ANGLE = digis.addProperty("angle", 0.0);
00058         IgProperty COUNT = digis.addProperty("countsTDC", int(0));
00059         IgProperty NUMBER = digis.addProperty("number", int(0));
00060 
00061         IgProperty CELL_WIDTH = digis.addProperty("cellWidth", 0.0);
00062         IgProperty CELL_LENGTH = digis.addProperty("cellLength", 0.0);
00063         IgProperty CELL_HEIGHT = digis.addProperty("cellHeight", 0.0);
00064         
00065         for ( DTDigiCollection::DigiRangeIterator dri = collection->begin();
00066               dri != collection->end(); ++dri )
00067         {
00068             const DTLayerId& dtlayerId = (*dri).first;
00069             const DTDigiCollection::Range& range = (*dri).second;
00070 
00071             const DTLayer* layer = geom->layer(dtlayerId);
00072 
00073             const DTChamber* chamber = layer->chamber();
00074             const DTChamberId chamberId = chamber->id();
00075 
00076             //const GeomDetUnit* det = geom->idToDetUnit((*dri).first);
00077 
00078             int layerId = dtlayerId.layer();
00079             int superLayerId = dtlayerId.superlayerId().superLayer();
00080             int sectorId = dtlayerId.superlayerId().chamberId().sector();                 
00081             int stationId = dtlayerId.superlayerId().chamberId().station();
00082             int wheelId = dtlayerId.superlayerId().chamberId().wheel();
00083 
00084             for ( DTDigiCollection::const_iterator dit = range.first;
00085                   dit != range.second; ++dit )
00086             {
00087                 IgCollectionItem digi = digis.create();
00088 
00089                 digi[LAYER_ID] = static_cast<int>(layerId);
00090                 digi[SUPERLAYER_ID] = static_cast<int>(superLayerId);
00091                 digi[SECTOR_ID] = static_cast<int>(sectorId);
00092                 digi[STATION_ID] = static_cast<int>(stationId);
00093                 digi[WHEEL_ID] = static_cast<int>(wheelId);
00094            
00095                 int wireNumber = (*dit).wire();
00096                 digi[WIREN] = static_cast<int>(wireNumber);
00097 
00098                 double localXPos = layer->specificTopology().wirePosition(wireNumber);
00099                 LocalPoint localPos(localXPos, 0.0, 0.0);
00100                 GlobalPoint pos = layer->surface().toGlobal(localPos);
00101             
00102                 digi[POS] = IgV3d(static_cast<double>(pos.x()/100.0),
00103                                   static_cast<double>(pos.y()/100.0),
00104                                   static_cast<double>(pos.z()/100.0));
00105 
00106                 //NOTE: May need to hand-code transformation by hand.
00107                 // Leave this for next version.
00108                 /*
00109                 float angle;
00110                 SbVec3f axis;
00111                 VisTrackingGeometry::createRotation(det, axis, angle);
00112                 */
00113 
00114                 int countsTDC = (*dit).countsTDC();
00115                 digi[COUNT] = static_cast<int>(countsTDC);
00116                 int number = (*dit).number();
00117                 digi[NUMBER] = static_cast<int>(number);
00118                 
00119                 digi[CELL_WIDTH] =      
00120                     static_cast<double>(layer->specificTopology().cellWidth()/100.0);
00121                 digi[CELL_LENGTH] = 
00122                     static_cast<double>(layer->specificTopology().cellLenght()/100.0);  
00123                 digi[CELL_HEIGHT] = 
00124                     static_cast<double>(layer->specificTopology().cellHeight()/100.0);
00125             }
00126         }
00127         
00128     }
00129     
00130     else 
00131     {
00132         // friendlyName:moduleLabel:instanceName:processName
00133         std::string error = "### Error: DTDigis "
00134                             + edm::TypeID (typeid (DTDigiCollection)).friendlyClassName () + ":" 
00135                             + inputTag_.label() + ":"
00136                             + inputTag_.instance() + ":" 
00137                             + inputTag_.process() + " are not found.";
00138 
00139         IgDataStorage *storage = config->storage ();
00140         IgCollection &collection = storage->getCollection ("Errors_V1");
00141         IgProperty ERROR_MSG = collection.addProperty ("Error", std::string ());
00142         IgCollectionItem item = collection.create ();
00143         item [ERROR_MSG] = error;
00144     }
00145 }
00146 
00147 DEFINE_FWK_MODULE(VisDTDigi);

Generated on Tue Jun 9 17:50:07 2009 for CMSSW by  doxygen 1.5.4