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 <TMCProcess.h>
33 #include <TParticle.h>
34 #include <TParticlePDG.h>
35 #include <TDatabasePDG.h>
37 #include <TDirectory.h>
44 //_______________________________________________________________________
46 fParticles("TParticle", 1000),
61 // Default constructor
65 //_______________________________________________________________________
66 AliStack::AliStack(Int_t size, const char* /*evfoldname*/):
67 fParticles("TParticle",1000),
86 //_______________________________________________________________________
87 AliStack::AliStack(const AliStack& st):
89 fParticles("TParticle",1000),
90 fParticleMap(*(st.Particles())),
91 fParticleFileMap(st.fParticleFileMap),
94 fTreeK((TTree*)(st.fTreeK->Clone())),
95 fNtrack(st.GetNtrack()),
96 fNprimary(st.GetNprimary()),
107 //_______________________________________________________________________
108 void AliStack::Copy(TObject&) const
110 AliFatal("Not implemented!");
113 //_______________________________________________________________________
114 AliStack::~AliStack()
127 //_____________________________________________________________________________
128 void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg, const Float_t *pmom,
129 const Float_t *vpos, const Float_t *polar, Float_t tof,
130 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is)
133 // Load a track on the stack
135 // done 1 if the track has to be transported
137 // parent identifier of the parent track. -1 for a primary
139 // pmom momentum GeV/c
141 // polar polarisation
142 // tof time of flight in seconds
143 // mecha production mechanism
144 // ntr on output the number of the track stored
147 // const Float_t tlife=0;
150 // Here we get the static mass
151 // For MC is ok, but a more sophisticated method could be necessary
152 // if the calculated mass is required
153 // also, this method is potentially dangerous if the mass
154 // used in the MC is not the same of the PDG database
156 TParticlePDG* pmc = TDatabasePDG::Instance()->GetParticle(pdg);
158 Float_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
159 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
160 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
162 // 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",
163 // mass,e,fNtrack,pdg,parent,done,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],kS);
166 PushTrack(done, parent, pdg, pmom[0], pmom[1], pmom[2], e,
167 vpos[0], vpos[1], vpos[2], tof, polar[0], polar[1], polar[2],
168 mech, ntr, weight, is);
170 AliWarning(Form("Particle type %d not defined in PDG Database !", pdg));
171 AliWarning("Particle skipped !");
175 //_____________________________________________________________________________
176 void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg,
177 Double_t px, Double_t py, Double_t pz, Double_t e,
178 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
179 Double_t polx, Double_t poly, Double_t polz,
180 TMCProcess mech, Int_t &ntr, Double_t weight, Int_t is)
183 // Load a track on the stack
185 // done 1 if the track has to be transported
187 // parent identifier of the parent track. -1 for a primary
189 // kS generation status code
190 // px, py, pz momentum GeV/c
191 // vx, vy, vz position
192 // polar polarisation
193 // tof time of flight in seconds
194 // mech production mechanism
195 // ntr on output the number of the track stored
197 // New method interface:
198 // arguments were changed to be in correspondence with TParticle
200 // Note: the energy is not calculated from the static mass but
201 // it is passed by argument e.
203 const Int_t kFirstDaughter=-1;
204 const Int_t kLastDaughter=-1;
208 = new(fParticles[fLoadPoint++])
209 TParticle(pdg, is, parent, -1, kFirstDaughter, kLastDaughter,
210 px, py, pz, e, vx, vy, vz, tof);
212 particle->SetPolarisation(polx, poly, polz);
213 particle->SetWeight(weight);
214 particle->SetUniqueID(mech);
219 particle->SetBit(kDoneBit);
221 particle->SetBit(kTransportBit);
226 // Declare that the daughter information is valid
227 particle->SetBit(kDaughtersBit);
228 // Add the particle to the stack
230 fParticleMap.AddAtAndExpand(particle, fNtrack);//CHECK!!
233 particle = GetParticleMapEntry(parent);
235 particle->SetLastDaughter(fNtrack);
236 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
239 AliError(Form("Parent %d does not exist",parent));
243 // This is a primary track. Set high water mark for this event
246 // Set also number if primary tracks
247 fNprimary = fHgwmk+1;
253 //_____________________________________________________________________________
254 TParticle* AliStack::PopNextTrack(Int_t& itrack)
257 // Returns next track from stack of particles
261 TParticle* track = GetNextParticle();
265 track->SetBit(kDoneBit);
270 fCurrentTrack = track;
274 //_____________________________________________________________________________
275 TParticle* AliStack::PopPrimaryForTracking(Int_t i)
278 // Returns i-th primary particle if it is flagged to be tracked,
282 TParticle* particle = Particle(i);
284 if (!particle->TestBit(kDoneBit)) {
285 fCurrentTrack = particle;
292 //_____________________________________________________________________________
293 Bool_t AliStack::PurifyKine()
296 // Compress kinematic tree keeping only flagged particles
297 // and renaming the particle id's in all the hits
300 int nkeep = fHgwmk + 1, parent, i;
301 TParticle *part, *father;
302 fTrackLabelMap.Set(fParticleMap.GetLast()+1);
304 // Save in Header total number of tracks before compression
305 // If no tracks generated return now
306 if(fHgwmk+1 == fNtrack) return kFALSE;
308 // First pass, invalid Daughter information
309 for(i=0; i<fNtrack; i++) {
310 // Preset map, to be removed later
311 if(i<=fHgwmk) fTrackLabelMap[i]=i ;
313 fTrackLabelMap[i] = -99;
314 if((part=GetParticleMapEntry(i))) {
316 // Check of this track should be kept for physics reasons
317 if (KeepPhysics(part)) KeepTrack(i);
319 part->ResetBit(kDaughtersBit);
320 part->SetFirstDaughter(-1);
321 part->SetLastDaughter(-1);
325 // Invalid daughter information for the parent of the first particle
326 // generated. This may or may not be the current primary according to
327 // whether decays have been recorded among the primaries
328 part = GetParticleMapEntry(fHgwmk+1);
329 fParticleMap.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
330 // Second pass, build map between old and new numbering
331 for(i=fHgwmk+1; i<fNtrack; i++) {
332 if(fParticleMap.At(i)->TestBit(kKeepBit)) {
333 // This particle has to be kept
334 fTrackLabelMap[i]=nkeep;
335 // If old and new are different, have to move the pointer
336 if(i!=nkeep) fParticleMap[nkeep]=fParticleMap.At(i);
337 part = GetParticleMapEntry(nkeep);
338 // as the parent is always *before*, it must be already
339 // in place. This is what we are checking anyway!
340 if((parent=part->GetFirstMother())>fHgwmk) {
341 if(fTrackLabelMap[parent]==-99) Fatal("PurifyKine","fTrackLabelMap[%d] = -99!\n",parent);
342 else part->SetFirstMother(fTrackLabelMap[parent]);}
347 // Fix daughters information
348 for (i=fHgwmk+1; i<nkeep; i++) {
349 part = GetParticleMapEntry(i);
350 parent = part->GetFirstMother();
352 father = GetParticleMapEntry(parent);
353 if(father->TestBit(kDaughtersBit)) {
355 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
356 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
358 // Initialise daughters info for first pass
359 father->SetFirstDaughter(i);
360 father->SetLastDaughter(i);
361 father->SetBit(kDaughtersBit);
366 // Now the output bit, from fHgwmk to nkeep we write everything and we erase
367 if(nkeep > fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
368 for (i=fHgwmk+1; i<nkeep; ++i) {
369 fParticleBuffer = GetParticleMapEntry(i);
370 fParticleFileMap[i]=static_cast<Int_t>(TreeK()->GetEntries());
372 fParticleMap[i]=fParticleBuffer=0;
375 for (i = nkeep; i < fNtrack; ++i) fParticleMap[i]=0;
377 Int_t toshrink = fNtrack-fHgwmk-1;
378 fLoadPoint-=toshrink;
380 for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles.RemoveAt(i);
387 Bool_t AliStack::ReorderKine()
390 // In some transport code children might not come in a continuous sequence.
391 // In this case the stack has to be reordered in order to establish the
392 // mother daughter relation using index ranges.
394 if(fHgwmk+1 == fNtrack) return kFALSE;
397 // Howmany secondaries have been produced ?
398 Int_t nNew = fNtrack - fHgwmk - 1;
404 // Copy pointers to temporary array
405 TParticle** tmp = new TParticle*[nNew];
407 for (i = 0; i < nNew; i++) {
408 if (fParticleMap.At(fHgwmk + 1 + i)) {
409 tmp[i] = GetParticleMapEntry(fHgwmk + 1 + i);
420 Int_t loadPoint = fHgwmk + 1;
422 // Re-Push particles into stack
423 // The outer loop is over parents, the inner over children.
424 // -1 refers to the primary particle
426 for (i = -1; i < nNew-1; i++) {
430 ipa = tmp[0]->GetFirstMother();
431 parP = GetParticleMapEntry(ipa);
433 ipa = (fHgwmk + 1 + i);
434 // Skip deleted particles
435 if (!tmp[i]) continue;
436 // Skip particles without children
437 if (tmp[i]->GetFirstDaughter() == -1) continue;
440 // Reset daughter information
442 Int_t idaumin = parP->GetFirstDaughter() - fHgwmk - 1;
443 Int_t idaumax = parP->GetLastDaughter() - fHgwmk - 1;
444 parP->SetFirstDaughter(-1);
445 parP->SetLastDaughter(-1);
446 for (j = idaumin; j <= idaumax; j++) {
447 // Skip deleted particles
448 if (!tmp[j]) continue;
449 // Skip particles already handled
450 if (map1[j] != -99) continue;
451 Int_t jpa = tmp[j]->GetFirstMother();
452 // Check if daughter of current parent
454 fParticleMap[loadPoint] = tmp[j];
455 // Re-establish daughter information
456 parP->SetLastDaughter(loadPoint);
457 if (parP->GetFirstDaughter() == -1) parP->SetFirstDaughter(loadPoint);
458 // Set Mother information
460 tmp[j]->SetFirstMother(map1[i]);
464 // Increase load point
473 // Build map for remapping of hits
475 fTrackLabelMap.Set(fNtrack);
476 for (i = 0; i < fNtrack; i ++) {
478 fTrackLabelMap[i] = i;
480 fTrackLabelMap[i] = map1[i - fHgwmk -1];
483 } // new particles poduced
488 Bool_t AliStack::KeepPhysics(const TParticle* part)
491 // Some particles have to kept on the stack for reasons motivated
492 // by physics analysis. Decision is put here.
494 Bool_t keep = kFALSE;
496 Int_t parent = part->GetFirstMother();
497 if (parent >= 0 && parent <= fHgwmk) {
498 TParticle* father = GetParticleMapEntry(parent);
500 // Keep first-generation daughter from primaries with heavy flavor
502 Int_t kf = father->GetPdgCode();
506 if (kfl > 10) kfl/=100;
508 if (kfl > 10) kfl/=10;
509 if (kfl > 10) kfl/=10;
514 // e+e- from pair production of primary gammas
516 if ((part->GetUniqueID()) == kPPair) keep = kTRUE;
521 //_____________________________________________________________________________
522 void AliStack::FinishEvent()
525 // Write out the kinematics that was not yet filled
528 // Update event header
531 // Fatal("FinishEvent", "No kinematics tree is defined.");
532 // Don't panic this is a probably a lego run
537 if(TreeK()->GetEntries() ==0) {
538 // set the fParticleFileMap size for the first time
539 fParticleFileMap.Set(fHgwmk+1);
542 Bool_t allFilled = kFALSE;
544 for(Int_t i=0; i<fHgwmk+1; ++i) {
545 if((part=GetParticleMapEntry(i))) {
546 fParticleBuffer = part;
547 fParticleFileMap[i]= static_cast<Int_t>(TreeK()->GetEntries());
550 fParticleMap.AddAt(0,i);
552 // When all primaries were filled no particle!=0
553 // should be left => to be removed later.
554 if (allFilled) AliWarning(Form("Why != 0 part # %d?\n",i));
558 // // printf("Why = 0 part # %d?\n",i); => We know.
560 // we don't break now in order to be sure there is no
561 // particle !=0 left.
562 // To be removed later and replaced with break.
563 if(!allFilled) allFilled = kTRUE;
567 //_____________________________________________________________________________
569 void AliStack::FlagTrack(Int_t track)
572 // Flags a track and all its family tree to be kept
579 particle = GetParticleMapEntry(curr);
581 // If the particle is flagged the three from here upward is saved already
582 if(particle->TestBit(kKeepBit)) return;
584 // Save this particle
585 particle->SetBit(kKeepBit);
587 // Move to father if any
588 if((curr=particle->GetFirstMother())==-1) return;
592 //_____________________________________________________________________________
593 void AliStack::KeepTrack(Int_t track)
596 // Flags a track to be kept
599 fParticleMap.At(track)->SetBit(kKeepBit);
602 //_____________________________________________________________________________
603 void AliStack::Clean(Int_t size)
606 // Reset stack data except for fTreeK
617 //_____________________________________________________________________________
618 void AliStack::Reset(Int_t size)
621 // Reset stack data including fTreeK
625 delete fParticleBuffer; fParticleBuffer = 0;
629 //_____________________________________________________________________________
630 void AliStack::ResetArrays(Int_t size)
633 // Resets stack arrays
636 fParticleMap.Clear();
637 if (size>0) fParticleMap.Expand(size);
640 //_____________________________________________________________________________
641 void AliStack::SetHighWaterMark(Int_t)
644 // Set high water mark for last track in event
648 fCurrentPrimary=fHgwmk;
649 // Set also number of primary tracks
650 fNprimary = fHgwmk+1;
653 //_____________________________________________________________________________
654 TParticle* AliStack::Particle(Int_t i)
657 // Return particle with specified ID
659 if(!fParticleMap.At(i)) {
660 Int_t nentries = fParticles.GetEntriesFast();
661 // algorithmic way of getting entry index
662 // (primary particles are filled after secondaries)
663 Int_t entry = TreeKEntry(i);
664 // check whether algorithmic way and
665 // and the fParticleFileMap[i] give the same;
666 // give the fatal error if not
667 if (entry != fParticleFileMap[i]) {
669 "!! The algorithmic way and map are different: !!\n entry: %d map: %d",
670 entry, fParticleFileMap[i]));
672 // Load particle at entry into fParticleBuffer
673 TreeK()->GetEntry(entry);
674 // Add to the TClonesarray
675 new (fParticles[nentries]) TParticle(*fParticleBuffer);
676 // Store a pointer in the TObjArray
677 fParticleMap.AddAt(fParticles[nentries],i);
679 return GetParticleMapEntry(i);
682 //_____________________________________________________________________________
683 TParticle* AliStack::ParticleFromTreeK(Int_t id) const
686 // return pointer to TParticle with label id
689 if ((entry = TreeKEntry(id)) < 0) return 0;
690 if (fTreeK->GetEntry(entry)<=0) return 0;
691 return fParticleBuffer;
694 //_____________________________________________________________________________
695 Int_t AliStack::TreeKEntry(Int_t id) const
698 // Return entry number in the TreeK for particle with label id
699 // Return negative number if label>fNtrack
701 // The order of particles in TreeK reflects the order of the transport of primaries and production of secondaries:
703 // Before transport there are fNprimary particles on the stack.
704 // They are transported one by one and secondaries (fNtrack - fNprimary) are produced.
705 // After the transport of each particles secondaries are written to the TreeK
706 // They occupy the entries 0 ... fNtrack - fNprimary - 1
707 // The primaries are written after they have been transported and occupy
708 // fNtrack - fNprimary .. fNtrack - 1
712 entry = id+fNtrack-fNprimary;
714 entry = id-fNprimary;
718 //_____________________________________________________________________________
719 Int_t AliStack::GetCurrentParentTrackNumber() const
722 // Return number of the parent of the current track
725 TParticle* current = GetParticleMapEntry(fCurrent);
728 return current->GetFirstMother();
730 AliWarning("Current track not found in the stack");
735 //_____________________________________________________________________________
736 Int_t AliStack::GetPrimary(Int_t id)
739 // Return number of primary that has generated track
747 parent=Particle(current)->GetFirstMother();
748 if(parent<0) return current;
752 //_____________________________________________________________________________
753 void AliStack::DumpPart (Int_t i) const
756 // Dumps particle i in the stack
758 GetParticleMapEntry(i)->Print();
761 //_____________________________________________________________________________
762 void AliStack::DumpPStack ()
765 // Dumps the particle stack
770 printf("\n\n=======================================================================\n");
771 for (i=0;i<fNtrack;i++)
773 TParticle* particle = Particle(i);
775 printf("-> %d ",i); particle->Print();
776 printf("--------------------------------------------------------------\n");
779 Warning("DumpPStack", "No particle with id %d.", i);
782 printf("\n=======================================================================\n\n");
784 // print particle file map
785 // printf("\nParticle file map: \n");
786 // for (i=0; i<fNtrack; i++)
787 // printf(" %d th entry: %d \n",i,fParticleFileMap[i]);
791 //_____________________________________________________________________________
792 void AliStack::DumpLoadedStack() const
795 // Dumps the particle in the stack
796 // that are loaded in memory.
800 "\n\n=======================================================================\n");
801 for (Int_t i=0;i<fNtrack;i++)
803 TParticle* particle = GetParticleMapEntry(i);
805 printf("-> %d ",i); particle->Print();
806 printf("--------------------------------------------------------------\n");
809 printf("-> %d Particle not loaded.\n",i);
810 printf("--------------------------------------------------------------\n");
814 "\n=======================================================================\n\n");
817 //_____________________________________________________________________________
818 void AliStack::SetCurrentTrack(Int_t track)
821 if (fCurrent < fNprimary) fCurrentTrack = Particle(track);
825 //_____________________________________________________________________________
830 //_____________________________________________________________________________
831 void AliStack::CleanParents()
834 // Clean particles stack
835 // Set parent/daughter relations
840 for(i=0; i<fHgwmk+1; i++) {
841 part = GetParticleMapEntry(i);
842 if(part) if(!part->TestBit(kDaughtersBit)) {
843 part->SetFirstDaughter(-1);
844 part->SetLastDaughter(-1);
849 //_____________________________________________________________________________
850 TParticle* AliStack::GetNextParticle()
853 // Return next particle from stack of particles
856 TParticle* particle = 0;
858 // search secondaries
859 //for(Int_t i=fNtrack-1; i>=0; i--) {
860 for(Int_t i=fNtrack-1; i>fHgwmk; i--) {
861 particle = GetParticleMapEntry(i);
862 if ((particle) && (!particle->TestBit(kDoneBit))) {
868 // take next primary if all secondaries were done
869 while (fCurrentPrimary>=0) {
870 fCurrent = fCurrentPrimary;
871 particle = GetParticleMapEntry(fCurrentPrimary--);
872 if ((particle) && (!particle->TestBit(kDoneBit))) {
877 // nothing to be tracked
883 //__________________________________________________________________________________________
885 void AliStack::ConnectTree(TTree* tree)
888 // Creates branch for writing particles
893 AliDebug(1, "Connecting TreeK");
898 AliFatal("Parameter is NULL");//we don't like such a jokes
901 return;//in this case TreeK() calls back this method (ConnectTree)
902 //tree after setting fTreeK, the rest was already executed
903 //it is safe to return now
906 // Create a branch for particles
908 AliDebug(2, Form("Tree name is %s",fTreeK->GetName()));
910 if (fTreeK->GetDirectory())
912 AliDebug(2, Form("and dir is %s",fTreeK->GetDirectory()->GetName()));
915 AliWarning("DIR IS NOT SET !!!");
917 TBranch *branch=fTreeK->GetBranch("Particles");
920 branch = fTreeK->Branch("Particles", "TParticle", &fParticleBuffer, 4000);
921 AliDebug(2, "Creating Branch in Tree");
925 AliDebug(2, "Branch Found in Tree");
926 branch->SetAddress(&fParticleBuffer);
928 if (branch->GetDirectory())
930 AliDebug(1, Form("Branch Dir Name is %s",branch->GetDirectory()->GetName()));
933 AliWarning("Branch Dir is NOT SET");
936 //_____________________________________________________________________________
938 Bool_t AliStack::GetEvent()
941 // Get new event from TreeK
943 // Reset/Create the particle stack
944 Int_t size = (Int_t)TreeK()->GetEntries();
948 //_____________________________________________________________________________
950 Bool_t AliStack::IsStable(Int_t pdg) const
953 // Decide whether particle (pdg) is stable
957 // All ions/nucleons are considered as stable
958 // Nuclear code is 10LZZZAAAI
959 if(pdg>1000000000)return kTRUE;
961 const Int_t kNstable = 18;
964 Int_t pdgStable[kNstable] = {
966 kElectron, // Electron
974 kLambda0, // Lambda_0
975 kSigmaMinus, // Sigma Minus
976 kSigmaPlus, // Sigma Plus
980 kNuE, // Electron Neutrino
981 kNuMu, // Muon Neutrino
982 kNuTau // Tau Neutrino
985 Bool_t isStable = kFALSE;
986 for (i = 0; i < kNstable; i++) {
987 if (pdg == TMath::Abs(pdgStable[i])) {
996 //_____________________________________________________________________________
997 Bool_t AliStack::IsPhysicalPrimary(Int_t index)
1000 // Test if a particle is a physical primary according to the following definition:
1001 // Particles produced in the collision including products of strong and
1002 // electromagnetic decay and excluding feed-down from weak decays of strange
1005 TParticle* p = Particle(index);
1006 Int_t ist = p->GetStatusCode();
1009 // Initial state particle
1010 if (ist > 1) return kFALSE;
1012 Int_t pdg = TMath::Abs(p->GetPdgCode());
1014 if (!IsStable(pdg)) return kFALSE;
1016 if (index < GetNprimary()) {
1018 // Particle produced by generator
1022 // Particle produced during transport
1025 Int_t imo = p->GetFirstMother();
1026 TParticle* pm = Particle(imo);
1027 Int_t mpdg = TMath::Abs(pm->GetPdgCode());
1028 // Check if it comes from a pi0 decay
1030 // What about the pi0 Dalitz ??
1031 // if ((mpdg == kPi0) && (imo < GetNprimary())) return kTRUE;
1033 // Check if this is a heavy flavor decay product
1034 Int_t mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1037 if (mfl < 4) return kFALSE;
1040 // Heavy flavor hadron produced by generator
1041 if (imo < GetNprimary()) {
1045 // To be sure that heavy flavor has not been produced in a secondary interaction
1046 // Loop back to the generated mother
1047 while (imo >= GetNprimary()) {
1048 imo = pm->GetFirstMother();
1051 mpdg = TMath::Abs(pm->GetPdgCode());
1052 mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1059 } // produced by generator ?