1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // Particles stack class //
21 // Implements the TMCVirtualStack of the Virtual Monte Carlo //
22 // Holds the particles transported during simulation //
23 // Is used to compare results of reconstruction with simulation //
26 ///////////////////////////////////////////////////////////////////////////////
29 #include <TClonesArray.h>
30 #include <TObjArray.h>
32 #include <TParticle.h>
33 #include <TParticlePDG.h>
35 #include <TDirectory.h>
42 //_______________________________________________________________________
44 fParticles("TParticle", 1000),
59 // Default constructor
63 //_______________________________________________________________________
64 AliStack::AliStack(Int_t size, const char* /*evfoldname*/):
65 fParticles("TParticle",1000),
84 //_______________________________________________________________________
85 AliStack::AliStack(const AliStack& st):
87 fParticles("TParticle",1000),
88 fParticleMap(*(st.Particles())),
89 fParticleFileMap(st.fParticleFileMap),
92 fTreeK((TTree*)(st.fTreeK->Clone())),
93 fNtrack(st.GetNtrack()),
94 fNprimary(st.GetNprimary()),
105 //_______________________________________________________________________
106 void AliStack::Copy(TObject&) const
108 AliFatal("Not implemented!");
111 //_______________________________________________________________________
112 AliStack::~AliStack()
125 //_____________________________________________________________________________
126 void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg, Float_t *pmom,
127 Float_t *vpos, Float_t *polar, Float_t tof,
128 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is)
131 // Load a track on the stack
133 // done 1 if the track has to be transported
135 // parent identifier of the parent track. -1 for a primary
137 // pmom momentum GeV/c
139 // polar polarisation
140 // tof time of flight in seconds
141 // mecha production mechanism
142 // ntr on output the number of the track stored
145 // const Float_t tlife=0;
148 // Here we get the static mass
149 // For MC is ok, but a more sophisticated method could be necessary
150 // if the calculated mass is required
151 // also, this method is potentially dangerous if the mass
152 // used in the MC is not the same of the PDG database
154 TParticlePDG* pmc = TDatabasePDG::Instance()->GetParticle(pdg);
156 Float_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
157 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
158 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
160 // printf("Loading mass %f ene %f No %d ip %d parent %d done %d pos %f %f %f mom %f %f %f kS %d m \n",
161 // mass,e,fNtrack,pdg,parent,done,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],kS);
164 PushTrack(done, parent, pdg, pmom[0], pmom[1], pmom[2], e,
165 vpos[0], vpos[1], vpos[2], tof, polar[0], polar[1], polar[2],
166 mech, ntr, weight, is);
168 AliWarning(Form("Particle type %d not defined in PDG Database !", pdg));
169 AliWarning("Particle skipped !");
173 //_____________________________________________________________________________
174 void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg,
175 Double_t px, Double_t py, Double_t pz, Double_t e,
176 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
177 Double_t polx, Double_t poly, Double_t polz,
178 TMCProcess mech, Int_t &ntr, Double_t weight, Int_t is)
181 // Load a track on the stack
183 // done 1 if the track has to be transported
185 // parent identifier of the parent track. -1 for a primary
187 // kS generation status code
188 // px, py, pz momentum GeV/c
189 // vx, vy, vz position
190 // polar polarisation
191 // tof time of flight in seconds
192 // mech production mechanism
193 // ntr on output the number of the track stored
195 // New method interface:
196 // arguments were changed to be in correspondence with TParticle
198 // Note: the energy is not calculated from the static mass but
199 // it is passed by argument e.
201 const Int_t kFirstDaughter=-1;
202 const Int_t kLastDaughter=-1;
206 = new(fParticles[fLoadPoint++])
207 TParticle(pdg, is, parent, -1, kFirstDaughter, kLastDaughter,
208 px, py, pz, e, vx, vy, vz, tof);
210 particle->SetPolarisation(polx, poly, polz);
211 particle->SetWeight(weight);
212 particle->SetUniqueID(mech);
217 particle->SetBit(kDoneBit);
219 particle->SetBit(kTransportBit);
224 // Declare that the daughter information is valid
225 particle->SetBit(kDaughtersBit);
226 // Add the particle to the stack
229 fParticleMap.AddAtAndExpand(particle, fNtrack);//CHECK!!
232 particle = dynamic_cast<TParticle*>(fParticleMap.At(parent));
234 particle->SetLastDaughter(fNtrack);
235 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
238 AliError(Form("Parent %d does not exist",parent));
242 // This is a primary track. Set high water mark for this event
245 // Set also number if primary tracks
246 fNprimary = fHgwmk+1;
252 //_____________________________________________________________________________
253 TParticle* AliStack::PopNextTrack(Int_t& itrack)
256 // Returns next track from stack of particles
260 TParticle* track = GetNextParticle();
264 track->SetBit(kDoneBit);
269 fCurrentTrack = track;
273 //_____________________________________________________________________________
274 TParticle* AliStack::PopPrimaryForTracking(Int_t i)
277 // Returns i-th primary particle if it is flagged to be tracked,
281 TParticle* particle = Particle(i);
283 if (!particle->TestBit(kDoneBit))
289 //_____________________________________________________________________________
290 Bool_t AliStack::PurifyKine()
293 // Compress kinematic tree keeping only flagged particles
294 // and renaming the particle id's in all the hits
297 int nkeep = fHgwmk + 1, parent, i;
298 TParticle *part, *father;
299 fTrackLabelMap.Set(fParticleMap.GetLast()+1);
301 // Save in Header total number of tracks before compression
302 // If no tracks generated return now
303 if(fHgwmk+1 == fNtrack) return kFALSE;
305 // First pass, invalid Daughter information
306 for(i=0; i<fNtrack; i++) {
307 // Preset map, to be removed later
308 if(i<=fHgwmk) fTrackLabelMap[i]=i ;
310 fTrackLabelMap[i] = -99;
311 if((part=dynamic_cast<TParticle*>(fParticleMap.At(i)))) {
313 // Check of this track should be kept for physics reasons
314 if (KeepPhysics(part)) KeepTrack(i);
316 part->ResetBit(kDaughtersBit);
317 part->SetFirstDaughter(-1);
318 part->SetLastDaughter(-1);
322 // Invalid daughter information for the parent of the first particle
323 // generated. This may or may not be the current primary according to
324 // whether decays have been recorded among the primaries
325 part = dynamic_cast<TParticle*>(fParticleMap.At(fHgwmk+1));
326 fParticleMap.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
327 // Second pass, build map between old and new numbering
328 for(i=fHgwmk+1; i<fNtrack; i++) {
329 if(fParticleMap.At(i)->TestBit(kKeepBit)) {
330 // This particle has to be kept
331 fTrackLabelMap[i]=nkeep;
332 // If old and new are different, have to move the pointer
333 if(i!=nkeep) fParticleMap[nkeep]=fParticleMap.At(i);
334 part = dynamic_cast<TParticle*>(fParticleMap.At(nkeep));
335 // as the parent is always *before*, it must be already
336 // in place. This is what we are checking anyway!
337 if((parent=part->GetFirstMother())>fHgwmk)
338 if(fTrackLabelMap[parent]==-99) Fatal("PurifyKine","fTrackLabelMap[%d] = -99!\n",parent);
339 else part->SetFirstMother(fTrackLabelMap[parent]);
344 // Fix daughters information
345 for (i=fHgwmk+1; i<nkeep; i++) {
346 part = dynamic_cast<TParticle*>(fParticleMap.At(i));
347 parent = part->GetFirstMother();
349 father = dynamic_cast<TParticle*>(fParticleMap.At(parent));
350 if(father->TestBit(kDaughtersBit)) {
352 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
353 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
355 // Initialise daughters info for first pass
356 father->SetFirstDaughter(i);
357 father->SetLastDaughter(i);
358 father->SetBit(kDaughtersBit);
363 // Now the output bit, from fHgwmk to nkeep we write everything and we erase
364 if(nkeep > fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
365 for (i=fHgwmk+1; i<nkeep; ++i) {
366 fParticleBuffer = dynamic_cast<TParticle*>(fParticleMap.At(i));
367 fParticleFileMap[i]=static_cast<Int_t>(TreeK()->GetEntries());
369 fParticleMap[i]=fParticleBuffer=0;
372 for (i = nkeep; i < fNtrack; ++i) fParticleMap[i]=0;
374 Int_t toshrink = fNtrack-fHgwmk-1;
375 fLoadPoint-=toshrink;
377 for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles.RemoveAt(i);
384 Bool_t AliStack::ReorderKine()
387 // In some transport code children might not come in a continuous sequence.
388 // In this case the stack has to be reordered in order to establish the
389 // mother daughter relation using index ranges.
391 if(fHgwmk+1 == fNtrack) return kFALSE;
394 // Howmany secondaries have been produced ?
395 Int_t nNew = fNtrack - fHgwmk - 1;
401 // Copy pointers to temporary array
402 TParticle** tmp = new TParticle*[nNew];
404 for (i = 0; i < nNew; i++) {
405 if (fParticleMap.At(fHgwmk + 1 + i)) {
406 tmp[i] = (TParticle*) (fParticleMap.At(fHgwmk + 1 + i));
417 Int_t loadPoint = fHgwmk + 1;
419 // Re-Push particles into stack
420 // The outer loop is over parents, the inner over children.
421 // -1 refers to the primary particle
423 for (i = -1; i < nNew-1; i++) {
427 ipa = tmp[0]->GetFirstMother();
428 parP =dynamic_cast<TParticle*>(fParticleMap.At(ipa));
430 ipa = (fHgwmk + 1 + i);
431 // Skip deleted particles
432 if (!tmp[i]) continue;
433 // Skip particles without children
434 if (tmp[i]->GetFirstDaughter() == -1) continue;
437 // Reset daughter information
439 Int_t idaumin = parP->GetFirstDaughter() - fHgwmk - 1;
440 Int_t idaumax = parP->GetLastDaughter() - fHgwmk - 1;
441 parP->SetFirstDaughter(-1);
442 parP->SetLastDaughter(-1);
443 for (j = idaumin; j <= idaumax; j++) {
444 // Skip deleted particles
445 if (!tmp[j]) continue;
446 // Skip particles already handled
447 if (map1[j] != -99) continue;
448 Int_t jpa = tmp[j]->GetFirstMother();
449 // Check if daughter of current parent
451 fParticleMap[loadPoint] = tmp[j];
452 // Re-establish daughter information
453 parP->SetLastDaughter(loadPoint);
454 if (parP->GetFirstDaughter() == -1) parP->SetFirstDaughter(loadPoint);
455 // Set Mother information
457 tmp[j]->SetFirstMother(map1[i]);
461 // Increase load point
470 // Build map for remapping of hits
472 fTrackLabelMap.Set(fNtrack);
473 for (i = 0; i < fNtrack; i ++) {
475 fTrackLabelMap[i] = i;
477 fTrackLabelMap[i] = map1[i - fHgwmk -1];
480 } // new particles poduced
485 Bool_t AliStack::KeepPhysics(TParticle* part)
488 // Some particles have to kept on the stack for reasons motivated
489 // by physics analysis. Decision is put here.
491 Bool_t keep = kFALSE;
493 // Keep first-generation daughter from primaries with heavy flavor
495 Int_t parent = part->GetFirstMother();
496 if (parent >= 0 && parent <= fHgwmk) {
497 TParticle* father = dynamic_cast<TParticle*>(Particles()->At(parent));
498 Int_t kf = father->GetPdgCode();
502 if (kfl > 10) kfl/=100;
504 if (kfl > 10) kfl/=10;
505 if (kfl > 10) kfl/=10;
513 //_____________________________________________________________________________
514 void AliStack::FinishEvent()
517 // Write out the kinematics that was not yet filled
520 // Update event header
523 // Fatal("FinishEvent", "No kinematics tree is defined.");
524 // Don't panic this is a probably a lego run
529 if(TreeK()->GetEntries() ==0) {
530 // set the fParticleFileMap size for the first time
531 fParticleFileMap.Set(fHgwmk+1);
534 Bool_t allFilled = kFALSE;
536 for(Int_t i=0; i<fHgwmk+1; ++i)
537 if((part=fParticleMap.At(i))) {
538 fParticleBuffer = dynamic_cast<TParticle*>(part);
539 fParticleFileMap[i]= static_cast<Int_t>(TreeK()->GetEntries());
542 fParticleMap.AddAt(0,i);
544 // When all primaries were filled no particle!=0
545 // should be left => to be removed later.
546 if (allFilled) AliWarning(Form("Why != 0 part # %d?\n",i));
550 // // printf("Why = 0 part # %d?\n",i); => We know.
552 // we don't break now in order to be sure there is no
553 // particle !=0 left.
554 // To be removed later and replaced with break.
555 if(!allFilled) allFilled = kTRUE;
558 //_____________________________________________________________________________
560 void AliStack::FlagTrack(Int_t track)
563 // Flags a track and all its family tree to be kept
570 particle=dynamic_cast<TParticle*>(fParticleMap.At(curr));
572 // If the particle is flagged the three from here upward is saved already
573 if(particle->TestBit(kKeepBit)) return;
575 // Save this particle
576 particle->SetBit(kKeepBit);
578 // Move to father if any
579 if((curr=particle->GetFirstMother())==-1) return;
583 //_____________________________________________________________________________
584 void AliStack::KeepTrack(Int_t track)
587 // Flags a track to be kept
590 fParticleMap.At(track)->SetBit(kKeepBit);
593 //_____________________________________________________________________________
594 void AliStack::Clean(Int_t size)
597 // Reset stack data except for fTreeK
608 //_____________________________________________________________________________
609 void AliStack::Reset(Int_t size)
612 // Reset stack data including fTreeK
616 delete fParticleBuffer; fParticleBuffer = 0;
620 //_____________________________________________________________________________
621 void AliStack::ResetArrays(Int_t size)
624 // Resets stack arrays
627 fParticleMap.Clear();
628 if (size>0) fParticleMap.Expand(size);
631 //_____________________________________________________________________________
632 void AliStack::SetHighWaterMark(Int_t)
635 // Set high water mark for last track in event
639 fCurrentPrimary=fHgwmk;
640 // Set also number of primary tracks
641 fNprimary = fHgwmk+1;
644 //_____________________________________________________________________________
645 TParticle* AliStack::Particle(Int_t i)
648 // Return particle with specified ID
650 if(!fParticleMap.At(i)) {
651 Int_t nentries = fParticles.GetEntriesFast();
652 // algorithmic way of getting entry index
653 // (primary particles are filled after secondaries)
654 Int_t entry = TreeKEntry(i);
655 // check whether algorithmic way and
656 // and the fParticleFileMap[i] give the same;
657 // give the fatal error if not
658 if (entry != fParticleFileMap[i]) {
660 "!! The algorithmic way and map are different: !!\n entry: %d map: %d",
661 entry, fParticleFileMap[i]));
663 // Load particle at entry into fParticleBuffer
664 TreeK()->GetEntry(entry);
665 // Add to the TClonesarray
666 new (fParticles[nentries]) TParticle(*fParticleBuffer);
667 // Store a pointer in the TObjArray
668 fParticleMap.AddAt(fParticles[nentries],i);
670 return dynamic_cast<TParticle*>(fParticleMap.At(i));
673 //_____________________________________________________________________________
674 TParticle* AliStack::ParticleFromTreeK(Int_t id) const
677 // return pointer to TParticle with label id
680 if ((entry = TreeKEntry(id)) < 0) return 0;
681 if (fTreeK->GetEntry(entry)<=0) return 0;
682 return fParticleBuffer;
685 //_____________________________________________________________________________
686 Int_t AliStack::TreeKEntry(Int_t id) const
689 // Return entry number in the TreeK for particle with label id
690 // Return negative number if label>fNtrack
692 // The order of particles in TreeK reflects the order of the transport of primaries and production of secondaries:
694 // Before transport there are fNprimary particles on the stack.
695 // They are transported one by one and secondaries (fNtrack - fNprimary) are produced.
696 // After the transport of each particles secondaries are written to the TreeK
697 // They occupy the entries 0 ... fNtrack - fNprimary - 1
698 // The primaries are written after they have been transported and occupy
699 // fNtrack - fNprimary .. fNtrack - 1
703 entry = id+fNtrack-fNprimary;
705 entry = id-fNprimary;
709 //_____________________________________________________________________________
710 Int_t AliStack::GetCurrentParentTrackNumber() const
713 // Return number of the parent of the current track
716 TParticle* current = (TParticle*)fParticleMap.At(fCurrent);
719 return current->GetFirstMother();
721 AliWarning("Current track not found in the stack");
726 //_____________________________________________________________________________
727 Int_t AliStack::GetPrimary(Int_t id)
730 // Return number of primary that has generated track
738 parent=Particle(current)->GetFirstMother();
739 if(parent<0) return current;
743 //_____________________________________________________________________________
744 void AliStack::DumpPart (Int_t i) const
747 // Dumps particle i in the stack
749 dynamic_cast<TParticle*>(fParticleMap.At(i))->Print();
752 //_____________________________________________________________________________
753 void AliStack::DumpPStack ()
756 // Dumps the particle stack
761 printf("\n\n=======================================================================\n");
762 for (i=0;i<fNtrack;i++)
764 TParticle* particle = Particle(i);
766 printf("-> %d ",i); particle->Print();
767 printf("--------------------------------------------------------------\n");
770 Warning("DumpPStack", "No particle with id %d.", i);
773 printf("\n=======================================================================\n\n");
775 // print particle file map
776 printf("\nParticle file map: \n");
777 for (i=0; i<fNtrack; i++)
778 printf(" %d th entry: %d \n",i,fParticleFileMap[i]);
782 //_____________________________________________________________________________
783 void AliStack::DumpLoadedStack() const
786 // Dumps the particle in the stack
787 // that are loaded in memory.
791 "\n\n=======================================================================\n");
792 for (Int_t i=0;i<fNtrack;i++)
794 TParticle* particle = dynamic_cast<TParticle*>(fParticleMap[i]);
796 printf("-> %d ",i); particle->Print();
797 printf("--------------------------------------------------------------\n");
800 printf("-> %d Particle not loaded.\n",i);
801 printf("--------------------------------------------------------------\n");
805 "\n=======================================================================\n\n");
812 //_____________________________________________________________________________
813 void AliStack::CleanParents()
816 // Clean particles stack
817 // Set parent/daughter relations
822 for(i=0; i<fHgwmk+1; i++) {
823 part = dynamic_cast<TParticle*>(fParticleMap.At(i));
824 if(part) if(!part->TestBit(kDaughtersBit)) {
825 part->SetFirstDaughter(-1);
826 part->SetLastDaughter(-1);
831 //_____________________________________________________________________________
832 TParticle* AliStack::GetNextParticle()
835 // Return next particle from stack of particles
838 TParticle* particle = 0;
840 // search secondaries
841 //for(Int_t i=fNtrack-1; i>=0; i--) {
842 for(Int_t i=fNtrack-1; i>fHgwmk; i--) {
843 particle = dynamic_cast<TParticle*>(fParticleMap.At(i));
844 if ((particle) && (!particle->TestBit(kDoneBit))) {
850 // take next primary if all secondaries were done
851 while (fCurrentPrimary>=0) {
852 fCurrent = fCurrentPrimary;
853 particle = dynamic_cast<TParticle*>(fParticleMap.At(fCurrentPrimary--));
854 if ((particle) && (!particle->TestBit(kDoneBit))) {
859 // nothing to be tracked
865 //__________________________________________________________________________________________
867 TTree* AliStack::TreeK()
872 //__________________________________________________________________________________________
874 void AliStack::ConnectTree(TTree* tree)
877 // Creates branch for writing particles
882 AliDebug(1, "Connecting TreeK");
887 AliFatal("Parameter is NULL");//we don't like such a jokes
890 return;//in this case TreeK() calls back this method (ConnectTree)
891 //tree after setting fTreeK, the rest was already executed
892 //it is safe to return now
895 // Create a branch for particles
897 AliDebug(2, Form("Tree name is %s",fTreeK->GetName()));
899 if (fTreeK->GetDirectory())
901 AliDebug(2, Form("and dir is %s",fTreeK->GetDirectory()->GetName()));
904 AliWarning("DIR IS NOT SET !!!");
906 TBranch *branch=fTreeK->GetBranch("Particles");
909 branch = fTreeK->Branch("Particles", "TParticle", &fParticleBuffer, 4000);
910 AliDebug(2, "Creating Branch in Tree");
914 AliDebug(2, "Branch Found in Tree");
915 branch->SetAddress(&fParticleBuffer);
917 if (branch->GetDirectory())
919 AliDebug(1, Form("Branch Dir Name is %s",branch->GetDirectory()->GetName()));
922 AliWarning("Branch Dir is NOT SET");
925 //_____________________________________________________________________________
927 Bool_t AliStack::GetEvent()
930 // Get new event from TreeK
932 // Reset/Create the particle stack
933 Int_t size = (Int_t)TreeK()->GetEntries();
937 //_____________________________________________________________________________
939 Bool_t AliStack::IsStable(Int_t pdg) const
942 // Decide whether particle (pdg) is stable
946 // All ions/nucleons are considered as stable
947 // Nuclear code is 10LZZZAAAI
948 if(pdg>1000000000)return kTRUE;
950 const Int_t kNstable = 15;
953 Int_t pdgStable[kNstable] = {
955 kElectron, // Electron
963 kLambda0, // Lambda_0
964 kSigmaMinus, // Sigma Minus
965 kSigmaPlus, // Sigma Plus
971 Bool_t isStable = kFALSE;
972 for (i = 0; i < kNstable; i++) {
973 if (pdg == TMath::Abs(pdgStable[i])) {
982 Bool_t AliStack::IsPhysicalPrimary(Int_t index)
985 // Test if a particle is a physical primary according to the following definition:
986 // Particles produced in the collision including products of strong and
987 // electromagnetic decay and excluding feed-down from weak decays of strange
990 TParticle* p = Particle(index);
991 Int_t ist = p->GetStatusCode();
994 // Initial state particle
995 if (ist > 1) return kFALSE;
997 Int_t pdg = TMath::Abs(p->GetPdgCode());
999 if (!IsStable(pdg)) return kFALSE;
1001 if (index < GetNprimary()) {
1003 // Particle produced by generator
1007 // Particle produced during transport
1010 Int_t imo = p->GetFirstMother();
1011 TParticle* pm = Particle(imo);
1012 Int_t mpdg = TMath::Abs(pm->GetPdgCode());
1013 // Check if it comes from a pi0 decay
1015 // What about the pi0 Dalitz ??
1016 // if ((mpdg == kPi0) && (imo < GetNprimary())) return kTRUE;
1018 // Check if this is a heavy flavor decay product
1019 Int_t mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1022 if (mfl < 4) return kFALSE;
1025 // Heavy flavor hadron produced by generator
1026 if (imo < GetNprimary()) {
1030 // To be sure that heavy flavor has not been produced in a secondary interaction
1031 // Loop back to the generated mother
1032 while (imo >= GetNprimary()) {
1033 imo = pm->GetFirstMother();
1036 mpdg = TMath::Abs(pm->GetPdgCode());
1037 mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1044 } // produced by generator ?