ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeFiller/src/FillerMCParticles.cc
Revision: 1.1
Committed: Fri Jul 25 11:33:58 2008 UTC (16 years, 9 months ago) by bendavid
Content type: text/plain
Branch: MAIN
Log Message:
Merged gen and sim particles into new MCParticle class

File Contents

# User Rev Content
1 bendavid 1.1 // $Id: FillerMCParticles.cc,v 1.13 2008/07/14 21:01:00 loizides Exp $
2    
3     #include "MitProd/TreeFiller/interface/FillerMCParticles.h"
4     #include "FWCore/MessageLogger/interface/MessageLogger.h"
5     #include "DataFormats/Common/interface/Handle.h"
6     #include "SimDataFormats/HepMCProduct/interface/HepMCProduct.h"
7     #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
8     #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
9     #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertex.h"
10     #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertexContainer.h"
11     #include "MitAna/DataTree/interface/Names.h"
12    
13     using namespace std;
14     using namespace edm;
15     using namespace mithep;
16    
17     //--------------------------------------------------------------------------------------------------
18     FillerMCParticles::FillerMCParticles(const ParameterSet &cfg, const char *name, bool active) :
19     BaseFiller(cfg, name, active),
20     genActive_(Conf().getUntrackedParameter<bool>("genActive",true)),
21     simActive_(Conf().getUntrackedParameter<bool>("simActive",true)),
22     genEdmName_(Conf().getUntrackedParameter<string>("genEdmName","source")),
23     simEdmName_(Conf().getUntrackedParameter<string>("simEdmName","mergedtruth:MergedTrackTruth")),
24     mitName_(Conf().getUntrackedParameter<string>("mitName",Names::gkMCPartBrn)),
25     mcParticles_(new mithep::MCParticleArr(250)),
26     genMap_(new mithep::GenParticleMap),
27     simMap_(new mithep::SimParticleMap)
28     {
29     // Constructor.
30     }
31    
32     //--------------------------------------------------------------------------------------------------
33     FillerMCParticles::~FillerMCParticles()
34     {
35     // Destructor.
36    
37     delete mcParticles_;
38     delete genMap_;
39     delete simMap_;
40     }
41    
42     //--------------------------------------------------------------------------------------------------
43     void FillerMCParticles::BookDataBlock(TreeWriter &tws)
44     {
45     // Add branch to tree.
46    
47     tws.AddBranch(mitName_.c_str(),&mcParticles_);
48    
49     }
50    
51     //--------------------------------------------------------------------------------------------------
52     void FillerMCParticles::FillDataBlock(const edm::Event &event,
53     const edm::EventSetup &setup)
54     {
55     // Loop over HepMC particle and fill their information.
56    
57     genMap_->Reset();
58     simMap_->Reset();
59     mcParticles_->Reset();
60    
61     if (genActive_) {
62    
63     Handle<edm::HepMCProduct> hHepMCProduct;
64     GetProduct(genEdmName_, hHepMCProduct, event);
65    
66     const HepMC::GenEvent &GenEvent = hHepMCProduct->getHepMCData();
67    
68     //loop over all hepmc particles and copy their information
69     for (HepMC::GenEvent::particle_const_iterator pgen = GenEvent.particles_begin();
70     pgen != GenEvent.particles_end(); ++pgen) {
71    
72     HepMC::GenParticle *mcPart = (*pgen);
73     if(!mcPart) continue;
74    
75     mithep::MCParticle *genParticle = mcParticles_->Allocate();
76     new (genParticle) mithep::MCParticle(mcPart->momentum().x(),mcPart->momentum().y(),
77     mcPart->momentum().z(),mcPart->momentum().e(),
78     mcPart->pdg_id(),
79     mcPart->status());
80    
81     genParticle->SetIsGenerated();
82    
83     genMap_->Add(mcPart->barcode(), genParticle);
84     }
85    
86     }
87    
88     if (simActive_) {
89    
90     Handle<TrackingParticleCollection> hTrackingParticleProduct;
91     GetProduct(simEdmName_, hTrackingParticleProduct, event);
92    
93     const TrackingParticleCollection trackingParticles = *(hTrackingParticleProduct.product());
94    
95     // loop through all simParticles
96     for (TrackingParticleCollection::const_iterator iM = trackingParticles.begin();
97     iM != trackingParticles.end(); ++iM) {
98    
99     TrackingParticleRef theRef(hTrackingParticleProduct, iM-trackingParticles.begin());
100     mithep::MCParticle *outSimParticle;
101    
102     if (genActive_ && iM->genParticle().size()) {
103     const HepMC::GenParticle *mcPart = iM->genParticle_begin()->get();
104     outSimParticle=genMap_->GetMit(mcPart->barcode());
105     outSimParticle->SetStatus(iM->status());
106     }
107     else {
108     outSimParticle = mcParticles_->Allocate();
109     new (outSimParticle) mithep::MCParticle(iM->px(),iM->py(),iM->pz(),iM->energy(),iM->pdgId(), iM->status());
110     }
111    
112     outSimParticle->SetIsSimulated();
113     simMap_->Add(theRef, outSimParticle);
114    
115     }
116    
117     }
118    
119    
120     mcParticles_->Trim();
121     }
122    
123     //--------------------------------------------------------------------------------------------------
124     void FillerMCParticles::ResolveLinks(const edm::Event &event,
125     const edm::EventSetup &setup)
126     {
127     // Loop over HepMC particle and resolve their links.
128    
129     if (genActive_) {
130    
131     Handle<edm::HepMCProduct> hHepMCProduct;
132     GetProduct(genEdmName_, hHepMCProduct, event);
133    
134     const HepMC::GenEvent &GenEvent = hHepMCProduct->getHepMCData();
135    
136     for (HepMC::GenEvent::particle_const_iterator pgen = GenEvent.particles_begin();
137     pgen != GenEvent.particles_end(); ++pgen) {
138    
139     HepMC::GenParticle *mcPart = (*pgen);
140     if(!mcPart) continue;
141    
142     //check if genpart has a decay vertex
143     HepMC::GenVertex *dVertex = mcPart->end_vertex();
144     if (!dVertex) continue;
145    
146     //find corresponding mithep genparticle parent in association table
147     mithep::MCParticle *genParent = genMap_->GetMit(mcPart->barcode());
148    
149     if (genParent->IsSimulated()) continue;
150    
151     //set decay vertex
152     //division by 10.0 is needed due to HepMC use of mm instead of cm for distance units
153     genParent->SetVertex(dVertex->point3d().x()/10.0,
154     dVertex->point3d().y()/10.0,
155     dVertex->point3d().z()/10.0);
156    
157     //loop through daugthers
158     for (HepMC::GenVertex::particles_out_const_iterator pgenD = dVertex->particles_out_const_begin();
159     pgenD != dVertex->particles_out_const_end(); ++pgenD) {
160     HepMC::GenParticle *mcDaughter = (*pgenD);
161     mithep::MCParticle *genDaughter = genMap_->GetMit(mcDaughter->barcode());
162     genParent->AddDaughter(genDaughter);
163     if (!(genDaughter->HasMother()))
164     genDaughter->SetMother(genParent);
165     }
166     }
167    
168     }
169    
170     if (simActive_) {
171    
172     Handle<TrackingParticleCollection> hTrackingParticleProduct;
173     GetProduct(simEdmName_, hTrackingParticleProduct, event);
174    
175     const TrackingParticleCollection trackingParticles = *(hTrackingParticleProduct.product());
176    
177     // loop through all simParticles
178     for (TrackingParticleCollection::const_iterator iM = trackingParticles.begin();
179     iM != trackingParticles.end(); ++iM) {
180    
181     if (iM->decayVertices().size() <= 0)
182     continue;
183    
184     TrackingParticleRef theRef(hTrackingParticleProduct, iM-trackingParticles.begin());
185     mithep::MCParticle *simParent = simMap_->GetMit(theRef);
186     for (TrackingVertexRefVector::iterator v= iM->decayVertices().begin();
187     v != iM->decayVertices().end(); ++v) {
188     for (TrackingParticleRefVector::iterator pd = v->get()->daughterTracks().begin();
189     pd != v->get()->daughterTracks().end(); ++pd) {
190     mithep::MCParticle *simDaughter = simMap_->GetMit(*pd);
191     simParent->AddDaughter(simDaughter);
192     simDaughter->SetMother(simParent);
193     }
194     if (v == iM->decayVertices().end()-1) {
195     simParent->SetVertex(v->get()->position().x(),
196     v->get()->position().y(),
197     v->get()->position().z());
198     }
199     }
200     }
201    
202     }
203    
204     }