ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitProd/TreeFiller/src/FillerMCParticles.cc
Revision: 1.3
Committed: Thu Jul 31 12:34:04 2008 UTC (16 years, 9 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: MITHEP_2_0_x
Changes since 1.2: +22 -13 lines
Log Message:
Consistently introduced ObjectService. Updated comments. Updated .cfg files. Switched off default handling of Stable and DecayParts (for now). ZmmFullReco.cfg shows how to use standard filler with extensions.

File Contents

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