]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - HMPID/AliHMPIDRecon.cxx
Routines needed for radiation scoring. (B. Pastircak)
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDRecon.cxx
index 0173f0a88ae79250a41c0f75212e786a8d232d1c..2260ae38dfc5b079f6b1f83b8c9950baeb05bd9d 100644 (file)
@@ -54,10 +54,12 @@ AliHMPIDRecon::AliHMPIDRecon():TTask("RichRec","RichPat"),
     fPhotWei [i] =  0;
   }
 //hidden algorithm
-  fMipX=fMipY=fThTrkFit=fPhTrkFit=fCkovFit=-999;
+  fMipX=fMipY=fThTrkFit=fPhTrkFit=fCkovFit=fMipQ=fRadX=fRadY=-999;
   fIdxMip=fNClu=0;
-  for (Int_t i=0; i<1000; i++) {
+  fCkovSig2=0;
+  for (Int_t i=0; i<100; i++) {
     fXClu[i] = fYClu[i] = 0;
+    fClCk[i] = kTRUE;
   }
 }
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
@@ -357,6 +359,9 @@ Double_t AliHMPIDRecon::Sigma2(Double_t ckovTh, Double_t ckovPh)const
   
   TVector3 v(-999,-999,-999);
   Double_t trkBeta = 1./(TMath::Cos(ckovTh)*fRadNmean);
+  
+  if(trkBeta > 1) trkBeta = 1;                 //protection against bad measured thetaCer  
+  if(trkBeta < 0) trkBeta = 0.0001;            //
 
   v.SetX(SigLoc (ckovTh,ckovPh,trkBeta));
   v.SetY(SigGeom(ckovTh,ckovPh,trkBeta));
@@ -373,20 +378,29 @@ Double_t AliHMPIDRecon::SigLoc(Double_t thetaC, Double_t phiC,Double_t betaM)con
 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
 //            MIP beta
 //   Returns: absolute error on Cerenkov angle, [radians]    
+  
   Double_t phiDelta = phiC - fTrkDir.Phi();
 
-  Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
-  Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);
+  Double_t sint     = TMath::Sin(fTrkDir.Theta());
+  Double_t cost     = TMath::Cos(fTrkDir.Theta());
+  Double_t sinf     = TMath::Sin(fTrkDir.Phi());
+  Double_t cosf     = TMath::Cos(fTrkDir.Phi());
+  Double_t sinfd    = TMath::Sin(phiDelta);
+  Double_t cosfd    = TMath::Cos(phiDelta);
+  Double_t tantheta = TMath::Tan(thetaC);
+  
+  Double_t alpha =cost-tantheta*cosfd*sint;                                                 // formula (11)
+  Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);                            // formula (after 8 in the text)
   if (k<0) return 1e10;
+  Double_t mu =sint*sinf+tantheta*(cost*cosfd*sinf+sinfd*cosf);                             // formula (10)
+  Double_t e  =sint*cosf+tantheta*(cost*cosfd*cosf-sinfd*sinf);                             // formula (9)
 
-  Double_t mu =TMath::Sin(fTrkDir.Theta())*TMath::Sin(fTrkDir.Phi())+TMath::Tan(thetaC)*(TMath::Cos(fTrkDir.Theta())*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Phi())+TMath::Sin(phiDelta)*TMath::Cos(fTrkDir.Phi()));
-  Double_t e  =TMath::Sin(fTrkDir.Theta())*TMath::Cos(fTrkDir.Phi())+TMath::Tan(thetaC)*(TMath::Cos(fTrkDir.Theta())*TMath::Cos(phiDelta)*TMath::Cos(fTrkDir.Phi())-TMath::Sin(phiDelta)*TMath::Sin(fTrkDir.Phi()));
+  Double_t kk = betaM*TMath::Sqrt(k)/(fgkGapThick*alpha);                                   // formula (6) and (7)
+  Double_t dtdxc = kk*(k*(cosfd*cosf-cost*sinfd*sinf)-(alpha*mu/(betaM*betaM))*sint*sinfd); // formula (6)           
+  Double_t dtdyc = kk*(k*(cosfd*sinf+cost*sinfd*cosf)+(alpha* e/(betaM*betaM))*sint*sinfd); // formula (7)            pag.4
 
-  Double_t kk = betaM*TMath::Sqrt(k)/(8*alpha);
-  Double_t dtdxc = kk*(k*(TMath::Cos(phiDelta)*TMath::Cos(fTrkDir.Phi())-TMath::Cos(fTrkDir.Theta())*TMath::Sin(phiDelta)*TMath::Sin(fTrkDir.Phi()))-(alpha*mu/(betaM*betaM))*TMath::Sin(fTrkDir.Theta())*TMath::Sin(phiDelta));
-  Double_t dtdyc = kk*(k*(TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Phi())+TMath::Cos(fTrkDir.Theta())*TMath::Sin(phiDelta)*TMath::Cos(fTrkDir.Phi()))+(alpha* e/(betaM*betaM))*TMath::Sin(fTrkDir.Theta())*TMath::Sin(phiDelta));
-
-  return  TMath::Sqrt(0.2*0.2*dtdxc*dtdxc + 0.25*0.25*dtdyc*dtdyc);
+  Double_t errX = 0.2,errY=0.25;                                                            //end of page 7
+  return  TMath::Sqrt(errX*errX*dtdxc*dtdxc + errY*errY*dtdyc*dtdyc);
 }
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 Double_t AliHMPIDRecon::SigCrom(Double_t thetaC, Double_t phiC,Double_t betaM)const
@@ -397,12 +411,19 @@ Double_t AliHMPIDRecon::SigCrom(Double_t thetaC, Double_t phiC,Double_t betaM)co
 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
 //            MIP beta
 //   Returns: absolute error on Cerenkov angle, [radians]    
+  
   Double_t phiDelta = phiC - fTrkDir.Phi();
-  Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
 
-  Double_t dtdn = TMath::Cos(fTrkDir.Theta())*fRadNmean*betaM*betaM/(alpha*TMath::Tan(thetaC));
+  Double_t sint     = TMath::Sin(fTrkDir.Theta());
+  Double_t cost     = TMath::Cos(fTrkDir.Theta());
+  Double_t cosfd    = TMath::Cos(phiDelta);
+  Double_t tantheta = TMath::Tan(thetaC);
+  
+  Double_t alpha =cost-tantheta*cosfd*sint;                                                 // formula (11)
+  Double_t dtdn = cost*fRadNmean*betaM*betaM/(alpha*tantheta);                              // formula (12)
             
-  Double_t f = 0.00928*(7.75-5.635)/TMath::Sqrt(12.);
+//  Double_t f = 0.00928*(7.75-5.635)/TMath::Sqrt(12.);
+  Double_t f = 0.0172*(7.75-5.635)/TMath::Sqrt(24.);
 
   return f*dtdn;
 }//SigCrom()
@@ -417,29 +438,37 @@ Double_t AliHMPIDRecon::SigGeom(Double_t thetaC, Double_t phiC,Double_t betaM)co
 //   Returns: absolute error on Cerenkov angle, [radians]    
 
   Double_t phiDelta = phiC - fTrkDir.Phi();
-  Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
 
-  Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);
+  Double_t sint     = TMath::Sin(fTrkDir.Theta());
+  Double_t cost     = TMath::Cos(fTrkDir.Theta());
+  Double_t sinf     = TMath::Sin(fTrkDir.Phi());
+  Double_t cosfd    = TMath::Cos(phiDelta);
+  Double_t costheta = TMath::Cos(thetaC);
+  Double_t tantheta = TMath::Tan(thetaC);
+  
+  Double_t alpha =cost-tantheta*cosfd*sint;                                                  // formula (11)
+  
+  Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);                             // formula (after 8 in the text)
   if (k<0) return 1e10;
 
-  Double_t eTr = 0.5*1.5*betaM*TMath::Sqrt(k)/(8*alpha);
-  Double_t lambda = 1.-TMath::Sin(fTrkDir.Theta())*TMath::Sin(fTrkDir.Theta())*TMath::Sin(phiC)*TMath::Sin(phiC);
-
-  Double_t c = 1./(1.+ eTr*k/(alpha*alpha*TMath::Cos(thetaC)*TMath::Cos(thetaC)));
-  Double_t i = betaM*TMath::Tan(thetaC)*lambda*TMath::Power(k,1.5);
-  Double_t ii = 1.+eTr*betaM*i;
+  Double_t eTr = 0.5*fgkRadThick*betaM*TMath::Sqrt(k)/(fgkGapThick*alpha);                   // formula (14)
+  Double_t lambda = 1.-sint*sint*sinf*sinf;                                                  // formula (15)
 
-  Double_t err = c * (i/(alpha*alpha*8) +  ii*(1.-lambda) / ( alpha*alpha*8*betaM*(1.+eTr)) );
-  Double_t trErr = 1.5/(TMath::Sqrt(12.)*TMath::Cos(fTrkDir.Theta()));
+  Double_t c1 = 1./(1.+ eTr*k/(alpha*alpha*costheta*costheta));                              // formula (13.a)
+  Double_t c2 = betaM*TMath::Power(k,1.5)*tantheta*lambda/(fgkGapThick*alpha*alpha);         // formula (13.b)
+  Double_t c3 = (1.+eTr*k*betaM*betaM)/((1+eTr)*alpha*alpha);                                // formula (13.c)
+  Double_t c4 = TMath::Sqrt(k)*tantheta*(1-lambda)/(fgkGapThick*betaM);                      // formula (13.d)
+  Double_t dtdT = c1 * (c2+c3*c4);
+  Double_t trErr = fgkRadThick/(TMath::Sqrt(12.)*cost);
 
-  return trErr*err;
+  return trErr*dtdT;
 }//SigGeom()
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 //
 // From here HTA....
 //
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
-Int_t AliHMPIDRecon::CkovHiddenTrk(AliESDtrack *pTrk,TClonesArray *pCluLst,Double_t nmean)
+Bool_t AliHMPIDRecon::CkovHiddenTrk(AliESDtrack *pTrk,TClonesArray *pCluLst,Double_t nmean)
 {
 // Pattern recognition method without any infos from tracking:HTA (Hidden Track Algorithm)...
 // The method finds in the chmber the cluster with the highest charge
@@ -453,32 +482,43 @@ Int_t AliHMPIDRecon::CkovHiddenTrk(AliESDtrack *pTrk,TClonesArray *pCluLst,Doubl
     
   fRadNmean=nmean;
 
+  if(pCluLst->GetEntriesFast()>100) return kFALSE;                                            //boundary check for CluX,CluY...
   Float_t mipX=-1,mipY=-1;Int_t mipId=-1,mipQ=-1;                                                                           
-  fPhotCnt=0;                                                      
   Double_t qRef = 0;
-  fNClu = pCluLst->GetEntriesFast();
-  for (Int_t iClu=0;iClu<fNClu;iClu++){                                                       //clusters loop
+  Int_t nCh=0;
+  for (Int_t iClu=0;iClu<pCluLst->GetEntriesFast();iClu++){                                   //clusters loop
     AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluLst->UncheckedAt(iClu);                       //get pointer to current cluster    
+    nCh = pClu->Ch();
     fXClu[iClu] = pClu->X();fYClu[iClu] = pClu->Y();                                          //store x,y for fitting procedure
+    fClCk[iClu] = kTRUE;                                                                      //all cluster are accepted at this stage to be reconstructed
     if(pClu->Q()>qRef){                                                                       //searching the highest charge to select a MIP      
       qRef = pClu->Q();
       mipId=iClu; mipX=pClu->X();mipY=pClu->Y();mipQ=(Int_t)pClu->Q();
     }                                                                                    
   }//clusters loop
 
-  if(qRef>pParam->QCut()){                                                                       //charge compartible with MIP clusters      
+  fNClu = pCluLst->GetEntriesFast();
+  if(qRef>pParam->QCut()){                                                                     //charge compartible with MIP clusters
     fIdxMip = mipId;
+    fClCk[mipId] = kFALSE;
     fMipX = mipX; fMipY=mipY; fMipQ = qRef;
-    if(!DoRecHiddenTrk()) return 1;                                                               //Do track and ring reconstruction,if problems returns 1
-    pTrk->SetHMPIDtrk(fRadX,fRadY,fThTrkFit,fPhTrkFit);                                          //store track intersection info
-    pTrk->SetHMPIDmip(fMipX,fMipY,(Int_t)fMipQ,fNClu);                                           //store mip info 
-    pTrk->SetHMPIDcluIdx(pCluLst->GetUniqueID(),fIdxMip);                                        //set cham number and index of cluster
-    pTrk->SetHMPIDsignal(fCkovFit);                                                              //find best Theta ckov for ring i.e. track
+    if(!DoRecHiddenTrk(pCluLst)) {
+      pTrk->SetHMPIDsignal(kNoPhotAccept);
+      return kFALSE;
+    }                                                                           //Do track and ring reconstruction,if problems returns 1
+    pTrk->SetHMPIDtrk(fRadX,fRadY,fThTrkFit,fPhTrkFit);                                        //store track intersection info
+    pTrk->SetHMPIDmip(fMipX,fMipY,(Int_t)fMipQ,fNClu);                                         //store mip info 
+    pTrk->SetHMPIDcluIdx(nCh,fIdxMip);                                                         //set cham number and index of cluster
+    pTrk->SetHMPIDsignal(fCkovFit);                                                            //find best Theta ckov for ring i.e. track
+    pTrk->SetHMPIDchi2(fCkovSig2);                                                             //errors squared
+//    Printf(" n clusters tot %i accepted %i",pCluLst->GetEntriesFast(),fNClu);
+    return kTRUE;
   }
-  return 0;
+  
+  return kFALSE;
 }//CkovHiddenTrk()
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
-Bool_t AliHMPIDRecon::DoRecHiddenTrk()
+Bool_t AliHMPIDRecon::DoRecHiddenTrk(TClonesArray *pCluLst)
 {
 // Pattern recognition method without any infos from tracking...
 // First a preclustering filter to avoid part of the noise
@@ -489,17 +529,27 @@ Bool_t AliHMPIDRecon::DoRecHiddenTrk()
 // Arguments:   none
 //   Returns:   none
   Double_t phiRec;
-  CluPreFilter();
+  if(!CluPreFilter(pCluLst)) {return kFALSE;}
   if(!FitEllipse(phiRec)) {return kFALSE;}
-  return FitFree(phiRec);
-}
+  Int_t nClTmp1 = pCluLst->GetEntriesFast()-1;  //minus MIP...
+  Int_t nClTmp2 = 0;
+  while(nClTmp1 != nClTmp2){
+    SetNClu(pCluLst->GetEntriesFast());
+    if(!FitFree(phiRec)) {return kFALSE;}
+    nClTmp2 = NClu();
+    if(nClTmp2!=nClTmp1) {nClTmp1=nClTmp2;nClTmp2=0;}
+  }
+  fNClu = nClTmp2;
+  return kTRUE;
+}//DoRecHiddenTrk()
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
-void AliHMPIDRecon::CluPreFilter()
+Bool_t AliHMPIDRecon::CluPreFilter(TClonesArray *pCluLst)
 {
 // Filter of bkg clusters
 // based on elliptical-shapes...
 //
-    ;
+  if(pCluLst->GetEntriesFast()>50||pCluLst->GetEntriesFast()<4) return kFALSE; 
+  else return kTRUE;
 }
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 Bool_t AliHMPIDRecon::FitEllipse(Double_t &phiRec)
@@ -525,7 +575,6 @@ Bool_t AliHMPIDRecon::FitEllipse(Double_t &phiRec)
   //       HG - AF
   //  b = --------
   //       AB - H^2
-  
   Double_t cA,cB,cF,cG,cH;
   Double_t aArg=-1;      Int_t iErrFlg;                                                //tmp vars for TMinuit
 
@@ -544,7 +593,7 @@ Bool_t AliHMPIDRecon::FitEllipse(Double_t &phiRec)
   gMinuit->mnparm(3," G ",1,0.01,0,0,iErrFlg);
   gMinuit->mnparm(4," F ",1,0.01,0,0,iErrFlg);
 
-  gMinuit->mnexcm("SIMPLEX" ,&aArg,0,iErrFlg);
+  gMinuit->mnexcm("SIMPLEX",&aArg,0,iErrFlg);
   gMinuit->mnexcm("MIGRAD" ,&aArg,0,iErrFlg);   
   gMinuit->mnpout(0,sName,cA,d1,d2,d3,iErrFlg);
   gMinuit->mnpout(1,sName,cB,d1,d2,d3,iErrFlg);
@@ -554,24 +603,26 @@ Bool_t AliHMPIDRecon::FitEllipse(Double_t &phiRec)
   delete gMinuit;
 
   Double_t i2 = cA*cB-cH*cH;                                       //quartic invariant : i2 > 0  ellipse, i2 < 0 hyperbola
+  if(i2<=0) return kFALSE;
   Double_t aX = (cH*cF-cB*cG)/i2;                                  //x centre of the canonical section 
   Double_t bY = (cH*cG-cA*cF)/i2;                                  //y centre of the canonical section 
   Double_t alfa1 = TMath::ATan(2*cH/(cA-cB));                      //alpha = angle of rotation of the conical section
   if(alfa1<0) alfa1+=TMath::Pi(); 
   alfa1*=0.5;
-  Double_t alfa2 = alfa1+TMath::Pi();
-  Double_t phiref = TMath::ATan2(bY-fMipY,aX-fMipX);               //evaluate in a unique way the angle of rotation comapring it
-  if(phiref<0) phiref+=TMath::TwoPi();                             //with the vector that poinst to the centre from the mip 
+//  Double_t alfa2 = alfa1+TMath::Pi();
+  Double_t phiref = TMath::ATan2(bY-fMipY,aX-fMipX);               //evaluate in a unique way the angle of rotation comparing it
+  if(phiref<0) phiref+=TMath::TwoPi();                             //with the vector that points to the centre from the mip 
   if(i2<0) phiref+=TMath::Pi();
   if(phiref>TMath::TwoPi()) phiref-=TMath::TwoPi();
 
 //  Printf(" alfa1 %f",alfa1*TMath::RadToDeg());
 //  Printf(" alfa2 %f",alfa2*TMath::RadToDeg());
 //  Printf(" firef %f",phiref*TMath::RadToDeg());
-  if(TMath::Abs(alfa1-phiref)<TMath::Abs(alfa2-phiref)) phiRec = alfa1; else phiRec = alfa2;  
+//  if(TMath::Abs(alfa1-phiref)<TMath::Abs(alfa2-phiref)) phiRec = alfa1; else phiRec = alfa2;  
   
-//  cout << " phi reconstructed " << phiRec*TMath::RadToDeg() << endl;
-  return (i2>0);
+//  Printf("FitEllipse: phi reconstructed %f",phiRec*TMath::RadToDeg());
+  phiRec=phiref;
+  return kTRUE;
 //
 }
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
@@ -594,9 +645,6 @@ Bool_t AliHMPIDRecon::FitFree(Double_t phiRec)
   TString sName;
   Double_t th,ph;
   
-  gMinuit->mnexcm("SET PRI",&aArg,1,iErrFlg);                                          //suspend all printout from TMinuit
-  gMinuit->mnexcm("SET NOW",&aArg,0,iErrFlg);
-
   gMinuit->mnparm(0," theta ",  0.01,0.01,0,TMath::PiOver2(),iErrFlg);
   gMinuit->mnparm(1," phi   ",phiRec,0.01,0,TMath::TwoPi()  ,iErrFlg);
   
@@ -646,21 +694,48 @@ void AliHMPIDRecon::FunMinPhot(Int_t &/* */,Double_t* /* */,Double_t &f,Double_t
   pRec->SetRadXY(xrad,yrad);
   pRec->SetTrack(xrad,yrad,thTrk,phTrk);
 
-  Double_t meanCkov=0;
+  Double_t meanCkov =0;
   Double_t meanCkov2=0;
   Double_t thetaCer,phiCer;
-  for(Int_t i=0;i<pRec->NClu();i++) {
+  Int_t nClAcc = 0;
+  Int_t nClTot=pRec->NClu();
+    
+  for(Int_t i=0;i<nClTot;i++) {
+    if(!(pRec->ClCk(i))) continue;
     pRec->FindPhotCkov(pRec->XClu(i),pRec->YClu(i),thetaCer,phiCer);  
-
     meanCkov  += thetaCer;
     meanCkov2 += thetaCer*thetaCer;
+    nClAcc++;
   }
-  meanCkov/=pRec->NClu();
-  Double_t rms = TMath::Sqrt(meanCkov2/pRec->NClu() - meanCkov*meanCkov);
-  f = rms/TMath::Sqrt(pRec->NClu());
-  Printf(" mean %f rms/sqrt(n) %f",meanCkov,f);  
-  if(iflag==3) pRec->SetCkovFit(meanCkov);
+  if(nClAcc==0) {f=999;return;}
+  meanCkov/=nClAcc;
+  Double_t rms = (meanCkov2 - meanCkov*meanCkov*nClAcc)/nClAcc;
+  if(rms<0) Printf(" rms2 = %f, strange!!!",rms);
+  rms = TMath::Sqrt(rms);
+  f = rms/TMath::Sqrt(nClAcc);
   
+  
+  if(iflag==3) {
+    Printf("FunMinPhot before: photons candidates %i used %i",nClTot,nClAcc);
+    nClAcc = 0;
+    Double_t meanCkov1=0;
+    Double_t meanCkov2=0;
+    for(Int_t i=0;i<nClTot;i++) {
+      if(!(pRec->ClCk(i))) continue;
+      pRec->FindPhotCkov(pRec->XClu(i),pRec->YClu(i),thetaCer,phiCer);  
+      if(TMath::Abs(thetaCer-meanCkov)<2*rms) {
+        meanCkov1 += thetaCer;
+        meanCkov2 += thetaCer*thetaCer;
+        nClAcc++;
+      } else pRec->SetClCk(i,kFALSE);
+    }
+    meanCkov1/=nClAcc;
+    Double_t rms2 = (meanCkov2 - meanCkov*meanCkov*nClAcc)/nClAcc;
+    Printf("FunMinPhot after: photons candidates %i used %i thetaCer %f",nClTot,nClAcc,meanCkov1);
+    pRec->SetCkovFit(meanCkov1);
+    pRec->SetCkovSig2(rms2);
+    pRec->SetNClu(nClAcc);
+  }
 }//FunMinPhot()
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 //