CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/SimMuon/Neutron/src/SubsystemNeutronWriter.cc

Go to the documentation of this file.
00001 #include "SimMuon/Neutron/interface/SubsystemNeutronWriter.h"
00002 #include "FWCore/Framework/interface/Event.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "FWCore/Utilities/interface/EDMException.h"
00006 #include "SimMuon/Neutron/src/AsciiNeutronWriter.h"
00007 #include "SimMuon/Neutron/src/RootNeutronWriter.h"
00008 #include "SimMuon/Neutron/src/NeutronWriter.h"
00009 #include "SimMuon/Neutron/src/EDMNeutronWriter.h"
00010 #include "CLHEP/Random/RandomEngine.h"
00011 #include "FWCore/ServiceRegistry/interface/Service.h"
00012 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00013 
00014 #include <algorithm>
00015 
00016 using namespace std;
00017 
00018 class SortByTime {
00019 public:
00020   bool operator()(const PSimHit & h1, const PSimHit & h2) {
00021    return (h1.tof() < h2.tof());
00022   }
00023 };
00024 
00025 
00026 SubsystemNeutronWriter::SubsystemNeutronWriter(edm::ParameterSet const& pset) 
00027 : theHitWriter(0),
00028   theRandFlat(0),
00029   theInputTag(pset.getParameter<edm::InputTag>("input")),
00030   theNeutronTimeCut(pset.getParameter<double>("neutronTimeCut")),
00031   theTimeWindow(pset.getParameter<double>("timeWindow")),
00032   theT0(pset.getParameter<double>("t0")),
00033   theNEvents(0),
00034   initialized(false),
00035   useLocalDetId_(true)
00036 {
00037   string writer = pset.getParameter<string>("writer");
00038   string output = pset.getParameter<string>("output");
00039   if(writer == "ASCII")
00040   {
00041     theHitWriter = new AsciiNeutronWriter(output);
00042   }
00043   else if (writer == "ROOT")
00044   {
00045     theHitWriter = new RootNeutronWriter(output);
00046   }
00047   else if (writer == "EDM")
00048   {
00049     produces<edm::PSimHitContainer>();
00050     theHitWriter = new EDMNeutronWriter();
00051     // write out the real DetId, not the local one
00052     useLocalDetId_ = false;
00053     // smear the times
00054     edm::Service<edm::RandomNumberGenerator> rng;
00055     if ( ! rng.isAvailable()) {
00056      throw cms::Exception("Configuration")
00057        << "SubsystemNeutronWriter requires the RandomNumberGeneratorService\n"
00058           "which is not present in the configuration file.  You must add the service\n"
00059           "in the configuration file or remove the modules that require it.";
00060     }
00061     CLHEP::HepRandomEngine& engine = rng->getEngine();
00062     theRandFlat = new CLHEP::RandFlat(engine);
00063   }
00064   else 
00065   {
00066     throw cms::Exception("NeutronWriter") << "Bad writer: "
00067       << writer;
00068   } 
00069 }
00070 
00071 
00072 SubsystemNeutronWriter::~SubsystemNeutronWriter()
00073 {
00074   printStats();
00075   delete theHitWriter;
00076   delete theRandFlat;
00077 }
00078 
00079 
00080 
00081 void SubsystemNeutronWriter::printStats() {
00082   edm::LogInfo("SubsystemNeutronWriter") << "SubsystemNeutronWriter Statistics:\n";
00083   for(map<int,int>::iterator mapItr = theCountPerChamberType.begin();
00084       mapItr != theCountPerChamberType.end();  ++mapItr) {
00085      edm::LogInfo("SubsystemNeutronWriter") << "   theEventOccupancy[" << mapItr->first << "] = "
00086          << mapItr->second << " / " << theNEvents << " / NCT \n";
00087   }
00088 }
00089 
00090 
00091 void SubsystemNeutronWriter::produce(edm::Event & e, edm::EventSetup const& c)
00092 {
00093   theHitWriter->beginEvent(e,c);
00094   ++theNEvents;
00095   edm::Handle<edm::PSimHitContainer> hits;
00096   e.getByLabel(theInputTag, hits);
00097 
00098   // sort hits by chamber
00099   map<int, edm::PSimHitContainer> hitsByChamber;
00100   for(edm::PSimHitContainer::const_iterator hitItr = hits->begin();
00101       hitItr != hits->end(); ++hitItr)
00102   {
00103     int chamberIndex = chamberId(hitItr->detUnitId());
00104     hitsByChamber[chamberIndex].push_back(*hitItr);
00105   }
00106 
00107   // now write out each chamber's contents
00108   for(map<int, edm::PSimHitContainer>::iterator hitsByChamberItr = hitsByChamber.begin();
00109       hitsByChamberItr != hitsByChamber.end(); ++hitsByChamberItr)
00110   {
00111     int chambertype = chamberType(hitsByChamberItr->first);
00112     writeHits(chambertype, hitsByChamberItr->second);
00113   }
00114   theHitWriter->endEvent();
00115 }
00116 
00117 
00118 void SubsystemNeutronWriter::initialize(int chamberType)
00119 {
00120   // should instantiate one of every chamber type, just so
00121   // ROOT knows what file to associate them with
00122   theHitWriter->initialize(chamberType);
00123 }
00124 
00125 
00126 void SubsystemNeutronWriter::writeHits(int chamberType, edm::PSimHitContainer & chamberHits)
00127 {
00128 
00129   sort(chamberHits.begin(), chamberHits.end(), SortByTime());
00130   edm::PSimHitContainer cluster;
00131   float startTime = -1000.;
00132   float smearing = 0.;
00133   for(size_t i = 0; i < chamberHits.size(); ++i) {
00134     PSimHit hit = chamberHits[i];
00135     float tof = hit.tof();
00136     LogDebug("SubsystemNeutronWriter") << "found hit from part type " << hit.particleType()
00137                    << " at tof " << tof << " p " << hit.pabs() 
00138                    << " on det " << hit.detUnitId() 
00139                    << " chamber type " << chamberType;
00140     if(tof > theNeutronTimeCut) {
00141       if(tof > (startTime + theTimeWindow) ) { // 1st in cluster
00142         startTime = tof;
00143         // set the time to be [t0, t0+25] at start of event
00144         smearing = theT0;
00145         if(theRandFlat) {
00146           smearing += theRandFlat->fire(25.);
00147         }
00148         if(!cluster.empty()) {
00149           LogDebug("SubsystemNeutronWriter") << "filling old cluster";
00150           writeCluster(chamberType, cluster);
00151           cluster.clear();
00152         }
00153         LogDebug("SubsystemNeutronWriter") << "starting neutron cluster at time " << startTime 
00154           << " on detType " << chamberType;
00155       }
00156       adjust(hit, -1.*startTime, smearing);
00157       cluster.push_back( hit );
00158     }
00159   }
00160   // get any leftover cluster
00161   if(!cluster.empty()) {
00162     writeCluster(chamberType, cluster);
00163   }
00164 }
00165 
00166 
00167 void SubsystemNeutronWriter::writeCluster(int chamberType, const edm::PSimHitContainer & cluster)
00168 {
00169   if(accept(cluster))
00170   {
00171     theHitWriter->writeCluster(chamberType, cluster);
00172     updateCount(chamberType);
00173   }
00174 }
00175 
00176 
00177 void SubsystemNeutronWriter::adjust(PSimHit & h, float timeOffset, float smearing) {
00178   unsigned int detId = useLocalDetId_ ? localDetId(h.detUnitId()) : h.detUnitId();
00179   float htime = h.timeOfFlight() + timeOffset + smearing;
00180   // prevent float precision loss
00181   if (h.timeOfFlight() > 1.E+6) {
00182     htime = smearing;
00183   }
00184   h = PSimHit( h.entryPoint(), h.exitPoint(), h.pabs(), htime,
00185                h.energyLoss(), h.particleType(), 
00186                detId, h.trackId(),
00187                h.momentumAtEntry().theta(),
00188                h.momentumAtEntry().phi() );
00189 }
00190 
00191 
00192 void SubsystemNeutronWriter::updateCount(int chamberType) {
00193   map<int,int>::iterator entry = theCountPerChamberType.find(chamberType);
00194   if(entry == theCountPerChamberType.end()) {
00195     theCountPerChamberType.insert( pair<int,int>(chamberType, 1) );
00196   } else {
00197     ++(entry->second);
00198   }
00199 }
00200 
00201