Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/eventgenerator/HepMC/MCTruth/src/MCTruthManager.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 /// \file eventgenerator/HepMC/MCTruth/src/MCTruthManager.cc
 27 /// \brief Implementation of the MCTruthManager class
 28 //
 29 //
 30 //
 31 //
 32 // --------------------------------------------------------------
 33 //      GEANT 4 - MCTruthManager class
 34 // --------------------------------------------------------------
 35 //
 36 // Author: Witold POKORSKI (Witold.Pokorski@cern.ch)
 37 //
 38 // --------------------------------------------------------------
 39 //
 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 41 
 42 #include "MCTruthManager.hh"
 43 
 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 45 
 46 static MCTruthManager* instance = 0;
 47 
 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 49 
 50 MCTruthManager::MCTruthManager() : fEvent(0), fConfig(0) {}
 51 
 52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 53 
 54 MCTruthManager::~MCTruthManager() {}
 55 
 56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 57 
 58 MCTruthManager* MCTruthManager::GetInstance()
 59 {
 60   if (instance == 0) {
 61     instance = new MCTruthManager();
 62   }
 63   return instance;
 64 }
 65 
 66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 67 
 68 void MCTruthManager::NewEvent()
 69 {
 70   // first delete the old event
 71   delete fEvent;
 72   // and now instaciate a new one
 73   fEvent = new HepMC::GenEvent();
 74 }
 75 
 76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 77 
 78 void MCTruthManager::AddParticle(G4LorentzVector& momentum, G4LorentzVector& prodpos,
 79                                  G4LorentzVector& endpos, G4int pdg_id, G4int partID,
 80                                  G4int motherID, G4bool directParent)
 81 {
 82   // we create a new particle with barcode = partID
 83   HepMC::GenParticle* particle = new HepMC::GenParticle(momentum, pdg_id);
 84   particle->suggest_barcode(partID);
 85   // we initialize the 'segmentations' map
 86   // for the moment particle is not 'segmented'
 87   fSegmentations[partID] = 1;
 88 
 89   // we create the GenVertex corresponding to the end point of the track
 90   HepMC::GenVertex* endvertex = new HepMC::GenVertex(endpos);
 91 
 92   // barcode of the endvertex = - barcode of the track
 93   endvertex->suggest_barcode(-partID);
 94   endvertex->add_particle_in(particle);
 95   fEvent->add_vertex(endvertex);
 96 
 97   if (motherID)  // not a primary
 98   {
 99     // here we could try to improve speed by searching only through particles which
100     // belong to the given primary tree
101     HepMC::GenParticle* mother = fEvent->barcode_to_particle(motherID);
102     //
103     if (mother) {
104       // we first check whether the mother's end vertex corresponds to the particle's
105       // production vertex
106       HepMC::GenVertex* motherendvtx = mother->end_vertex();
107       HepMC::FourVector mp0 = motherendvtx->position();
108       G4LorentzVector motherendpos(mp0.x(), mp0.y(), mp0.z(), mp0.t());
109 
110       if (motherendpos.x() == prodpos.x() && motherendpos.y() == prodpos.y()
111           && motherendpos.z() == prodpos.z())  // if yes, we attach the particle
112       {
113         motherendvtx->add_particle_out(particle);
114       }
115       else  // if not, we check whether the mother is biological or adopted
116       {
117         if (!directParent)  // adopted
118         {
119           G4bool found = false;
120 
121           // first check if any of the dummy particles
122           // has the end vertex at the right place
123           //
124           for (HepMC::GenVertex::particles_out_const_iterator it =
125                  motherendvtx->particles_out_const_begin();
126                it != motherendvtx->particles_out_const_end(); it++)
127           {
128             if ((*it)->pdg_id() == -999999) {
129               HepMC::FourVector dp0 = (*it)->end_vertex()->position();
130               G4LorentzVector dummypos(dp0.x(), dp0.y(), dp0.z(), dp0.t());
131               ;
132 
133               if (dummypos.x() == prodpos.x() && dummypos.y() == prodpos.y()
134                   && dummypos.z() == prodpos.z())
135               {
136                 (*it)->end_vertex()->add_particle_out(particle);
137                 found = true;
138                 break;
139               }
140             }
141           }
142 
143           // and if not, create a dummy particle connecting
144           // to the end vertex of the mother
145           //
146           if (!found) {
147             HepMC::GenVertex* childvtx = new HepMC::GenVertex(prodpos);
148             childvtx->add_particle_out(particle);
149 
150             // the dummy vertex gets the barcode -500000
151             // minus the daughter particle barcode
152             //
153             childvtx->suggest_barcode(-500000 - partID);
154             fEvent->add_vertex(childvtx);
155 
156             HepMC::GenParticle* dummypart = new HepMC::GenParticle(G4LorentzVector(), -999999);
157 
158             // the dummy particle gets the barcode 500000
159             // plus the daughter particle barcode
160             //
161             dummypart->suggest_barcode(500000 + partID);
162             childvtx->add_particle_in(dummypart);
163             motherendvtx->add_particle_out(dummypart);
164           }
165         }
166         else  // biological
167         {
168           // in case mother was already 'split' we need to look for
169           // the right 'segment' to add the new daugther.
170           // We use Time coordinate to locate the place for the new vertex
171 
172           G4int number_of_segments = fSegmentations[motherID];
173           G4int segment = 0;
174 
175           // we loop through the segments
176           //
177           while (!((mother->end_vertex()->position().t() > prodpos.t())
178                    && (mother->production_vertex()->position().t() < prodpos.t())))
179           {
180             segment++;
181             if (segment == number_of_segments)
182               G4cerr << "Problem!!!! Time coordinates incompatible!" << G4endl;
183 
184             mother = fEvent->barcode_to_particle(segment * 10000000 + motherID);
185           }
186 
187           // now, we 'split' the appropriate 'segment' of the mother particle
188           // into two particles and create a new vertex
189           //
190           HepMC::GenVertex* childvtx = new HepMC::GenVertex(prodpos);
191           childvtx->add_particle_out(particle);
192           fEvent->add_vertex(childvtx);
193 
194           // we first detach the mother from its original vertex
195           //
196           HepMC::GenVertex* orig_mother_end_vtx = mother->end_vertex();
197           orig_mother_end_vtx->remove_particle(mother);
198 
199           // and attach it to the new vertex
200           //
201           childvtx->add_particle_in(mother);
202 
203           // now we create a new particle representing the mother after
204           // interaction the barcode of the new particle is 10000000 + the
205           // original barcode
206           //
207           HepMC::GenParticle* mothertwo = new HepMC::GenParticle(*mother);
208           mothertwo->suggest_barcode(fSegmentations[motherID] * 10000000 + mother->barcode());
209 
210           // we also reset the barcodes of the vertices
211           //
212           orig_mother_end_vtx->suggest_barcode(-fSegmentations[motherID] * 10000000
213                                                - mother->barcode());
214           childvtx->suggest_barcode(-mother->barcode());
215 
216           // we attach it to the new vertex where interaction took place
217           //
218           childvtx->add_particle_out(mothertwo);
219 
220           // and we attach it to the original endvertex
221           //
222           orig_mother_end_vtx->add_particle_in(mothertwo);
223 
224           // and finally ... the increase the 'segmentation counter'
225           //
226           fSegmentations[motherID] = fSegmentations[motherID] + 1;
227         }
228       }
229     }
230     else
231     // mother GenParticle is not there for some reason...
232     // if this happens, we need to revise the philosophy...
233     // a solution would be to create HepMC particles
234     // at the begining of each track
235     {
236       G4cerr << "barcode " << motherID << " mother not there! " << G4endl;
237     }
238   }
239   else  // primary
240   {
241     HepMC::GenVertex* primaryvtx = new HepMC::GenVertex(prodpos);
242     primaryvtx->add_particle_out(particle);
243     fEvent->add_vertex(primaryvtx);
244 
245     // add id to the list of primaries
246     //
247     fPrimarybarcodes.push_back(partID);
248   }
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252 
253 void MCTruthManager::PrintEvent()
254 {
255   fEvent->print();
256 
257   // looping over primaries and print the decay tree for each of them
258   //
259   for (std::vector<int>::const_iterator primarybar = fPrimarybarcodes.begin();
260        primarybar != fPrimarybarcodes.end(); primarybar++)
261   {
262     PrintTree(fEvent->barcode_to_particle(*primarybar), " | ");
263   }
264 }
265 
266 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
267 
268 void MCTruthManager::PrintTree(HepMC::GenParticle* particle, G4String offset)
269 {
270   G4cout << offset << "---  barcode: " << particle->barcode() << " pdg: " << particle->pdg_id()
271          << " energy: " << particle->momentum().e()
272          << " production vertex: " << particle->production_vertex()->position().x() << ", "
273          << particle->production_vertex()->position().y() << ", "
274          << particle->production_vertex()->position().z() << ", "
275          << particle->production_vertex()->position().t() << G4endl;
276 
277   for (HepMC::GenVertex::particles_out_const_iterator it =
278          particle->end_vertex()->particles_out_const_begin();
279        it != particle->end_vertex()->particles_out_const_end(); it++)
280   {
281     G4String deltaoffset = "";
282 
283     G4int curr = std::fmod(double((*it)->barcode()), 10000000.);
284     G4int part = std::fmod(double(particle->barcode()), 10000000.);
285     if (curr != part) {
286       deltaoffset = " | ";
287     }
288 
289     PrintTree((*it), offset + deltaoffset);
290   }
291 }
292 
293 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
294