]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
TVirtualFitter in Cluster finding + derivatives analitically calculated. More robust...
authordibari <dibari@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 22 Aug 2007 12:13:26 +0000 (12:13 +0000)
committerdibari <dibari@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 22 Aug 2007 12:13:26 +0000 (12:13 +0000)
HMPID/AliHMPIDCluster.cxx
HMPID/AliHMPIDDigit.h
HMPID/AliHMPIDPreprocessor.cxx
HMPID/AliHMPIDRecon.cxx
HMPID/Hconfig.C

index b610d3ca0ab01e1af617ac4f1a0deab1c1a8471a..7ce3a5d68adb896d174b976612b3aa12bf077b2d 100644 (file)
 //  * provided "as is" without express or implied warranty.                  *
 //  **************************************************************************
 
-#include "AliHMPIDCluster.h"  //class header
+#include <TVirtualFitter.h>  //Solve()
 #include <TMinuit.h>         //Solve()
 #include <TClonesArray.h>    //Solve()
 #include <TMarker.h>         //Draw()
 
+#include "AliLog.h"          //FitFunc()
+
+#include "AliHMPIDCluster.h"  //class header
+
 Bool_t AliHMPIDCluster::fgDoCorrSin=kTRUE;
 
 ClassImp(AliHMPIDCluster)
+    
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 void AliHMPIDCluster::CoG()
 {
@@ -76,30 +81,71 @@ void AliHMPIDCluster::Draw(Option_t*)
   TMarker *pMark=new TMarker(X(),Y(),5); pMark->SetUniqueID(fSt);pMark->SetMarkerColor(kBlue); pMark->Draw();
 }
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
-void AliHMPIDCluster::FitFunc(Int_t &iNpars, Double_t* /*deriv*/, Double_t &chi2, Double_t *par, Int_t /* */)
+void AliHMPIDCluster::FitFunc(Int_t &iNpars, Double_t* deriv, Double_t &chi2, Double_t *par, Int_t /* */)
 {
 // Cluster fit function 
 // par[0]=x par[1]=y par[2]=q for the first Mathieson shape
 // par[3]=x par[4]=y par[5]=q for the second Mathieson shape and so on up to iNpars/3 Mathieson shapes
-// For each pad of the cluster calculates the difference between actual pad charge and the charge induced to this pad by all Mathieson distributions 
+// For each pad of the cluster calculates the difference between actual pad charge and the charge induced to this pad by all Mathieson distributions
 // Then the chi2 is calculated as the sum of this value squared for all pad in the cluster.  
 // Arguments: iNpars - number of parameters which is number of local maxima of cluster * 3
 //            chi2   - function result to be minimised 
 //            par   - parameters array of size iNpars            
 //   Returns: none  
-  AliHMPIDCluster *pClu=(AliHMPIDCluster*)gMinuit->GetObjectFit();
-  Int_t iNshape = iNpars/3;
-    
+  
+  AliHMPIDCluster *pClu=(AliHMPIDCluster*)TVirtualFitter::GetFitter()->GetObjectFit();
+
+  Int_t nPads = pClu->Size();  
+  Double_t **derivPart;
+  
+  derivPart = new Double_t*[iNpars];
   chi2 = 0;
-  for(Int_t i=0;i<pClu->Size();i++){                                                   //loop on all pads of the cluster
+  
+  Int_t iNshape = iNpars/3;
+  
+  for(Int_t j=0;j<iNpars;j++){                                                      
+    deriv[j] = 0;
+    derivPart[j] = new Double_t[nPads];
+    for(Int_t i=0;i<nPads;i++){                                                          
+      derivPart[j][i] = 0;
+    }
+  }
+  
+  if(iNshape>6) {Printf("HMPID Error!!: n. of clusters in FitFunc %i",iNshape);return;}
+  for(Int_t i=0;i<nPads;i++){                                                          //loop on all pads of the cluster
+    Double_t dQpadMath = 0;
+    for(Int_t j=0;j<iNshape;j++){                                                      //Mathiesons loop as all of them may contribute to this pad
+      Double_t fracMathi = pClu->Dig(i)->IntMathieson(par[3*j],par[3*j+1]);
+      dQpadMath+=par[3*j+2]*fracMathi;                                                 // par[3*j+2] is charge par[3*j] is x par[3*j+1] is y of current Mathieson
+      
+      derivPart[3*j  ][i] += par[3*j+2]*(pClu->Dig(i)->Mathieson(par[3*j]-pClu->Dig(i)->LorsX()-0.5*AliHMPIDParam::SizePadX())-
+                                         pClu->Dig(i)->Mathieson(par[3*j]-pClu->Dig(i)->LorsX()+0.5*AliHMPIDParam::SizePadX()))*
+                                         pClu->Dig(i)->IntPartMathi(par[3*j+1],2);
+      derivPart[3*j+1][i] += par[3*j+2]*(pClu->Dig(i)->Mathieson(par[3*j+1]-pClu->Dig(i)->LorsY()-0.5*AliHMPIDParam::SizePadY())-
+                                         pClu->Dig(i)->Mathieson(par[3*j+1]-pClu->Dig(i)->LorsY()+0.5*AliHMPIDParam::SizePadY()))*
+                                         pClu->Dig(i)->IntPartMathi(par[3*j],1);
+      derivPart[3*j+2][i] += fracMathi;
+    }
+    if(dQpadMath>0 && pClu->Dig(i)->Q()>0) {
+      chi2 +=TMath::Power((pClu->Dig(i)->Q()-dQpadMath),2)/pClu->Dig(i)->Q();          //chi2 function to be minimized
+    }
+  }
+                                                                                       //loop on all pads of the cluster     
+  for(Int_t i=0;i<nPads;i++){                                                          //loop on all pads of the cluster
     Double_t dQpadMath = 0;                                                            //pad charge collector  
     for(Int_t j=0;j<iNshape;j++){                                                      //Mathiesons loop as all of them may contribute to this pad
-      dQpadMath+=par[3*j+2]*pClu->Dig(i)->IntMathieson(par[3*j],par[3*j+1]);           // par[3*j+2] is charge par[3*j] is x par[3*j+1] is y of current Mathieson
+      Double_t fracMathi = pClu->Dig(i)->IntMathieson(par[3*j],par[3*j+1]);
+      dQpadMath+=par[3*j+2]*fracMathi;                                                 
+      if(dQpadMath>0 && pClu->Dig(i)->Q()>0) {
+        deriv[3*j]   += 2/pClu->Dig(i)->Q()*(pClu->Dig(i)->Q()-dQpadMath)*derivPart[3*j  ][i];
+        deriv[3*j+1] += 2/pClu->Dig(i)->Q()*(pClu->Dig(i)->Q()-dQpadMath)*derivPart[3*j+1][i];
+        deriv[3*j+2] += 2/pClu->Dig(i)->Q()*(pClu->Dig(i)->Q()-dQpadMath)*derivPart[3*j+2][i];
+      }
     }
-//    if(dQpadMath>0)chi2 +=TMath::Power((pClu->Dig(i)->Q()-dQpadMath),2)/dQpadMath;   //
-    if(dQpadMath>0 && pClu->Dig(i)->Q())
-      chi2 +=TMath::Power((pClu->Dig(i)->Q()-dQpadMath),2)/pClu->Dig(i)->Q(); //
-  }                                                                                    //loop on all pads of the cluster     
+  }
+  //delete array...
+  for(Int_t i=0;i<iNpars;i++) delete derivPart[i]; delete derivPart;
+  
 }//FitFunction()
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 void AliHMPIDCluster::Print(Option_t* opt)const
@@ -138,90 +184,116 @@ Int_t AliHMPIDCluster::Solve(TClonesArray *pCluLst,Bool_t isTryUnfold)
 //Arguments: pCluLst     - cluster list pointer where to add new cluster(s)
 //           isTryUnfold - flag to switch on/off unfolding   
 //  Returns: number of local maxima of original cluster
-  CoG();
+  const Int_t kMaxLocMax=6;                                                              //max allowed number of loc max for fitting
+//  
+  CoG();                                                                                 //First calculate CoG for the given cluster
   Int_t iCluCnt=pCluLst->GetEntriesFast();                                               //get current number of clusters already stored in the list by previous operations
   if(isTryUnfold==kFALSE || Size()==1) {                                                 //if cluster contains single pad there is no way to improve the knowledge 
     (isTryUnfold)?fSt=kSi1:fSt=kNot;
     new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this);  //add this raw cluster 
     return 1;
   } 
-//Phase 0. Initialise TMinuit  
-  const Int_t kMaxLocMax=6;                                                              //max allowed number of loc max for fitting
-  if(!gMinuit) gMinuit = new TMinuit(100);                                               //init MINUIT with this number of parameters (3 params per mathieson)
-  gMinuit->mncler();                                                                     // reset Minuit list of paramters
-  gMinuit->SetObjectFit((TObject*)this);  gMinuit->SetFCN(AliHMPIDCluster::FitFunc);     //set fit function
-  Double_t aArg=-1;                                                                    //tmp vars for TMinuit
-  Int_t iErrFlg;                   
-  gMinuit->mnexcm("SET PRI",&aArg,1,iErrFlg);                                          //suspend all printout from TMinuit 
-  gMinuit->mnexcm("SET NOW",&aArg,0,iErrFlg);                                          //suspend all warning printout from TMinuit
+  
+//Phase 0. Initialise Fitter  
+  Double_t arglist[10];
+  Int_t ierflg = 0;
+  TVirtualFitter *fitter = TVirtualFitter::Fitter(this,3*6);                            //initialize Fitter
+
+  arglist[0] = -1;
+  ierflg = fitter->ExecuteCommand("SET PRI", arglist, 1);                               // no printout
+  ierflg = fitter->ExecuteCommand("SET NOW", arglist, 0);                               //no warning messages
+  arglist[0] =  1;
+  ierflg = fitter->ExecuteCommand("SET GRA", arglist, 1);                               //force Fitter to use my gradient
+
+  fitter->SetFCN(AliHMPIDCluster::FitFunc);
+
+//  arglist[0] = 1;
+//  ierflg = fitter->ExecuteCommand("SET ERR", arglist ,1);
+  
+// Set starting values and step sizes for parameters
+    
 //Phase 1. Find number of local maxima. Strategy is to check if the current pad has QDC more then all neigbours. Also find the box contaning the cluster   
   fNlocMax=0;
 
- for(Int_t iDig1=0;iDig1<Size();iDig1++) {                                               //first digits loop
+  for(Int_t iDig1=0;iDig1<Size();iDig1++) {                                               //first digits loop
+    
     AliHMPIDDigit *pDig1 = Dig(iDig1);                                                   //take next digit    
     Int_t iCnt = 0;                                                                      //counts how many neighbouring pads has QDC more then current one
+    
     for(Int_t iDig2=0;iDig2<Size();iDig2++) {                                            //loop on all digits again
+      
       if(iDig1==iDig2) continue;                                                         //the same digit, no need to compare 
       AliHMPIDDigit *pDig2 = Dig(iDig2);                                                 //take second digit to compare with the first one
       Int_t dist = TMath::Sign(Int_t(pDig1->PadChX()-pDig2->PadChX()),1)+TMath::Sign(Int_t(pDig1->PadChY()-pDig2->PadChY()),1);//distance between pads
       if(dist==1)                                                                        //means dig2 is a neighbour of dig1
          if(pDig2->Q()>=pDig1->Q()) iCnt++;                                              //count number of pads with Q more then Q of current pad
+      
     }//second digits loop
+    
     if(iCnt==0&&fNlocMax<kMaxLocMax){                                                    //this pad has Q more then any neighbour so it's local maximum
+      
       Double_t xStart=pDig1->LorsX();Double_t yStart=pDig1->LorsY();
       Double_t xMin=xStart-AliHMPIDParam::SizePadX();
       Double_t xMax=xStart+AliHMPIDParam::SizePadX();
       Double_t yMin=yStart-AliHMPIDParam::SizePadY();
       Double_t yMax=yStart+AliHMPIDParam::SizePadY();
-      gMinuit->mnparm(3*fNlocMax  ,Form("x%i",fNlocMax),xStart,0.1,xMin,xMax,iErrFlg);   // X,Y,Q initial values of the loc max pad
-      gMinuit->mnparm(3*fNlocMax+1,Form("y%i",fNlocMax),yStart,0.1,yMin,yMax,iErrFlg);   // X, Y constrained to be near the loc max
-      gMinuit->mnparm(3*fNlocMax+2,Form("q%i",fNlocMax),pDig1->Q(),0.1,0,100000,iErrFlg);// Q constrained to be positive
+      
+      ierflg = fitter->SetParameter(3*fNlocMax  ,Form("x%i",fNlocMax),xStart,0.1,xMin,xMax);    // X,Y,Q initial values of the loc max pad
+      ierflg = fitter->SetParameter(3*fNlocMax+1,Form("y%i",fNlocMax),yStart,0.1,yMin,yMax);    // X, Y constrained to be near the loc max
+      ierflg = fitter->SetParameter(3*fNlocMax+2,Form("q%i",fNlocMax),pDig1->Q(),0.1,0,100000); // Q constrained to be positive
+      
       fNlocMax++;
+      
     }//if this pad is local maximum
   }//first digits loop
   
 //Phase 2. Fit loc max number of Mathiesons or add this current cluster to the list
 // case 1 -> no loc max found
- if ( fNlocMax == 0) {                                                                   // case of no local maxima found: pads with same charge...
-   gMinuit->mnparm(3*fNlocMax  ,Form("x%i",fNlocMax),fX,0.1,0,0,iErrFlg);                // Init values taken from CoG() -> fX,fY,fQRaw
-   gMinuit->mnparm(3*fNlocMax+1,Form("y%i",fNlocMax),fY,0.1,0,0,iErrFlg);                //
-   gMinuit->mnparm(3*fNlocMax+2,Form("q%i",fNlocMax),fQRaw,0.1,0,100000,iErrFlg);        //
+ if ( fNlocMax == 0) {                                                                       // case of no local maxima found: pads with same charge...
+   
+   ierflg = fitter->SetParameter(3*fNlocMax  ,Form("x%i",fNlocMax),fX,0.1,0,0);              // Init values taken from CoG() -> fX,fY,fQRaw
+   ierflg = fitter->SetParameter(3*fNlocMax+1,Form("y%i",fNlocMax),fY,0.1,0,0);              //
+   ierflg = fitter->SetParameter(3*fNlocMax+2,Form("q%i",fNlocMax),fQRaw,0.1,0,100000);      //
+   
    fNlocMax = 1;
    fSt=kNoLoc;
  }
 
 // case 2 -> loc max found. Check # of loc maxima 
- if ( fNlocMax >= kMaxLocMax)                                                            // if # of local maxima exceeds kMaxLocMax... 
-   {
-     fSt = kMax;   new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this);                   //...add this raw cluster  
-   }                                                                                     //or...
- else{                                                                                   //...resonable number of local maxima to fit and user requested it
-   Double_t arglist[10];arglist[0] = 10000;arglist[1] = 1.;                              //number of steps and sigma on pads charges  
-   gMinuit->mnexcm("SIMPLEX" ,arglist,2,iErrFlg);                                        //start fitting with Simplex
-   if (!iErrFlg)
-     gMinuit->mnexcm("MIGRAD" ,arglist,2,iErrFlg);                                       //fitting improved by Migrad
-   if(iErrFlg) {
+ if ( fNlocMax >= kMaxLocMax)  {                                                          // if # of local maxima exceeds kMaxLocMax...
+   fSt = kMax;   new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this);                      //...add this raw cluster  
+   } else {                                                                               //or resonable number of local maxima to fit and user requested it
+  // Now ready for minimization step
+   arglist[0] = 500;                                                                      //number of steps and sigma on pads charges
+   arglist[1] = 1.;                                                                       //
+
+   ierflg = fitter->ExecuteCommand("SIMPLEX",arglist,2);                                  //start fitting with Simplex
+   if (!ierflg)
+     fitter->ExecuteCommand("MIGRAD" ,arglist,2);                                         //fitting improved by Migrad
+   if(ierflg) {
      Double_t strategy=2;
-     gMinuit->mnexcm("SET STR",&strategy,1,iErrFlg);                                     //change level of strategy 
-     if(!iErrFlg) {
-       gMinuit->mnexcm("SIMPLEX" ,arglist,2,iErrFlg);
-       if (!iErrFlg)
-        gMinuit->mnexcm("MIGRAD" ,arglist,2,iErrFlg);                                   //fitting improved by Migrad
-//       Printf("Try to improve fit --> err %d",iErrFlg);
+     ierflg = fitter->ExecuteCommand("SET STR",&strategy,1);                              //change level of strategy 
+     if(!ierflg) {
+       ierflg = fitter->ExecuteCommand("SIMPLEX",arglist,2);                              //start fitting with Simplex
+       if (!ierflg)
+         fitter->ExecuteCommand("MIGRAD" ,arglist,2);                                     //fitting improved by Migrad
      }
    }        
-   if(iErrFlg) fSt=kAbn;                                                                 //no convergence of the fit...
-   Double_t dummy; TString sName;                                                        //vars to get results from Minuit
+   if(ierflg) fSt=kAbn;                                                                   //no convergence of the fit...
+   Double_t dummy; char sName[80];                                                        //vars to get results from Minuit
+   Double_t edm, errdef;
+   Int_t nvpar, nparx;
+  
    for(Int_t i=0;i<fNlocMax;i++){                                                        //store the local maxima parameters
-      gMinuit->mnpout(3*i   ,sName,  fX, fErrX , dummy, dummy, iErrFlg);                 // X 
-      gMinuit->mnpout(3*i+1 ,sName,  fY, fErrY , dummy, dummy, iErrFlg);                 // Y
-      gMinuit->mnpout(3*i+2 ,sName,  fQ, fErrQ , dummy, dummy, iErrFlg);                 // Q
-      gMinuit->mnstat(fChi2,dummy,dummy,iErrFlg,iErrFlg,iErrFlg);                        // Chi2 of the fit
+     fitter->GetParameter(3*i   ,sName,  fX, fErrX , dummy, dummy);                      // X
+     fitter->GetParameter(3*i+1 ,sName,  fY, fErrY , dummy, dummy);                      // Y
+     fitter->GetParameter(3*i+2 ,sName,  fQ, fErrQ , dummy, dummy);                      // Q
+     fitter->GetStats(fChi2, edm, errdef, nvpar, nparx);                                 //get fit infos
       if(fSt!=kAbn) {         
-        if(fNlocMax!=1)fSt=kUnf;                                                           // if unfolded
-        if(fNlocMax==1&&fSt!=kNoLoc) fSt=kLo1;                                             // if only 1 loc max
-        if ( !IsInPc()) fSt = kEdg;                                                        // if Out of Pc
-        if(fSt==kNoLoc) fNlocMax=0;                                                        // if with no loc max (pads with same charge..)
+        if(fNlocMax!=1)fSt=kUnf;                                                         // if unfolded
+        if(fNlocMax==1&&fSt!=kNoLoc) fSt=kLo1;                                           // if only 1 loc max
+        if ( !IsInPc()) fSt = kEdg;                                                      // if Out of Pc
+        if(fSt==kNoLoc) fNlocMax=0;                                                      // if with no loc max (pads with same charge..)
       }
       new ((*pCluLst)[iCluCnt++]) AliHMPIDCluster(*this);                               //add new unfolded cluster
    }
index 9a6fe7ecce30178fce69d0582db0c0b5fcbea363..116ffab4a1115c61bc4f06a6f26b3a48182955c7 100644 (file)
@@ -40,7 +40,9 @@ public:
 
          Float_t LorsY       (                               )const{return AliHMPIDParam::LorsY(AliHMPIDParam::A2P(fPad),AliHMPIDParam::A2Y(fPad));                               } //center of the pad y, [cm]
 //  
-  inline Float_t IntMathieson(Float_t x,Float_t y            )const;                                                                   //Mathieson distribution 
+  inline Float_t Mathieson   (Float_t x                      )const;                                                                   //Mathieson distribution 
+  inline Float_t IntPartMathi(Float_t z, Int_t axis          )const;                                                                   //integral in 1-dim of Mathieson
+  inline Float_t IntMathieson(Float_t x,Float_t y            )const;                                                                   //integral in 2-dim of Mathieson  
          Int_t   PadPcX      (                               )const{return AliHMPIDParam::A2X(fPad);}                                                 //pad pc x # 0..79
          Int_t   PadPcY      (                               )const{return AliHMPIDParam::A2Y(fPad);}                                                 //pad pc y # 0..47
          Int_t   PadChX      (                               )const{return (Pc()%2)*AliHMPIDParam::kPadPcX+PadPcX();}                                 //pad ch x # 0..159
@@ -81,19 +83,55 @@ Int_t AliHMPIDDigit::Compare(const TObject *pObj) const
 }
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 
+Float_t AliHMPIDDigit::Mathieson(Float_t x)const
+{
+// Mathieson function.
+// This is the answer to electrostatic problem of charge distrubution in MWPC described elsewhere. (NIM A370(1988)602-603)
+// Arguments: x- position of the center of Mathieson distribution
+//  Returns: value of the Mathieson function
+  Float_t  kK1=0.28278795,kK2=0.96242952, kSqrtK3 =0.77459667, kD=0.445;
+  Float_t lambda = x/kD;
+  Float_t a=1-TMath::TanH(kK2*lambda)*TMath::TanH(kK2*lambda);
+  Float_t b=1+kSqrtK3*kSqrtK3*TMath::TanH(kK2*lambda)*TMath::TanH(kK2*lambda);
+  Float_t mathi = kK1*a/b;
+  return mathi;
+}
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+
+Float_t AliHMPIDDigit::IntPartMathi(Float_t z, Int_t axis)const
+{
+// Integration of Mathieson.
+// This is the answer to electrostatic problem of charge distrubution in MWPC described elsewhere. (NIM A370(1988)602-603)
+// Arguments: x,y- position of the center of Mathieson distribution
+//  Returns: a charge fraction [0-1] imposed into the pad
+  Float_t shift1,shift2;
+  if(axis==1) {
+    shift1 = -LorsX()+0.5*AliHMPIDParam::SizePadX();
+    shift2 = -LorsX()-0.5*AliHMPIDParam::SizePadX();
+  } else {
+    shift1 = -LorsY()+0.5*AliHMPIDParam::SizePadY();
+    shift2 = -LorsY()-0.5*AliHMPIDParam::SizePadY();
+  }
+    
+  Float_t  kK2=0.96242952, kSqrtK3 =0.77459667,  kK4=0.37932926, kD=0.445;
+
+  Float_t ux1=kSqrtK3*TMath::TanH(kK2*(z+shift1)/kD);
+  Float_t ux2=kSqrtK3*TMath::TanH(kK2*(z+shift2)/kD);
+  
+  return kK4*(TMath::ATan(ux2)-TMath::ATan(ux1));
+}
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+
 Float_t AliHMPIDDigit::IntMathieson(Float_t x,Float_t y)const
 {
 // Integration of Mathieson.
 // This is the answer to electrostatic problem of charge distrubution in MWPC described elsewhere. (NIM A370(1988)602-603)
 // Arguments: x,y- position of the center of Mathieson distribution
 //  Returns: a charge fraction [0-1] imposed into the pad
-  Float_t  kK2=0.96242952, kSqrtK3 =0.77459667, kK4=0.37932926;
 
-  Float_t ux1=kSqrtK3*TMath::TanH(kK2*(x-LorsX()+0.5*AliHMPIDParam::SizePadX())/0.445);
-  Float_t ux2=kSqrtK3*TMath::TanH(kK2*(x-LorsX()-0.5*AliHMPIDParam::SizePadX())/0.445);
-  Float_t uy1=kSqrtK3*TMath::TanH(kK2*(y-LorsY()+0.5*AliHMPIDParam::SizePadY())/0.445);
-  Float_t uy2=kSqrtK3*TMath::TanH(kK2*(y-LorsY()-0.5*AliHMPIDParam::SizePadY())/0.445);
-  return 4*kK4*(TMath::ATan(ux2)-TMath::ATan(ux1))*kK4*(TMath::ATan(uy2)-TMath::ATan(uy1));
+  Float_t xm = IntPartMathi(x,1);
+  Float_t ym = IntPartMathi(y,2);
+  return 4*xm*ym;
 }
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
 void AliHMPIDDigit::Raw(UInt_t &w32,Int_t &ddl,Int_t &r,Int_t &d,Int_t &a)const
index 493fe51c37425ddb24e4614d8795264df5d3fda4..4062d1fbfb069f274deb8b72c260705b79043b9a 100644 (file)
@@ -118,9 +118,9 @@ Bool_t AliHMPIDPreprocessor::ProcDcs(TMap* pMap)
   metaData.SetResponsible("AliHMPIDPreprocessor"); 
   metaData.SetComment("SIMULATED");
 
-  stDcsStore =   Store("Calib","Qthre",&arQthre,&metaData,0,kTRUE) && 
-                 Store("Calib","Nmean",&arNmean,&metaData,0,kTRUE) &&
-                 Store("Calib","UserCut",&arUserCut,&metaData,0,kTRUE);
+  stDcsStore =   Store("Calib","Qthre",&arQthre,&metaData,0,kTRUE) &&    // from DCS 
+                 Store("Calib","Nmean",&arNmean,&metaData,0,kTRUE) &&    // from DCS
+                 Store("Calib","UserCut",&arUserCut,&metaData,0,kTRUE);  //really not from DCS...a method ProcManual maybe needed
   if(!stDcsStore) {
     Log("HMPID - failure to store DCS data results in OCDB");    
   }
index 2d7310337828890f2e64bdb263c9ffe01e3c2edb..e5b0c2385586b7d7f3ae46dc16360bc781ea65e5 100644 (file)
@@ -94,20 +94,25 @@ void AliHMPIDRecon::CkovAngle(AliESDtrack *pTrk,TClonesArray *pCluLst,Double_t n
       if(FindPhotCkov(pClu->X(),pClu->Y(),thetaCer,phiCer)){                                  //find ckov angle for this  photon candidate
         fPhotCkov[fPhotCnt]=thetaCer;                                                         //actual theta Cerenkov (in TRS)
         fPhotPhi [fPhotCnt]=phiCer;                                                           //actual phi   Cerenkov (in TRS): -pi to come back to "unusual" ref system (X,Y,-Z)
+        Printf("photon n. %i reconstructed theta = %f",fPhotCnt,fPhotCkov[fPhotCnt]);
         fPhotCnt++;                                                                           //increment counter of photon candidates
       }
     }
   }//clusters loop
+  if(fPhotCnt<=3) pTrk->SetHMPIDsignal(kNoPhotAccept);                                        //no reconstruction with <=3 photon candidates
   Int_t iNacc=FlagPhot(HoughResponse());                                                      //flag photons according to individual theta ckov with respect to most probable
   pTrk->SetHMPIDmip(mipX,mipY,mipQ,iNacc);                                                    //store mip info 
 
   if(mipId==-1)              {pTrk->SetHMPIDsignal(kMipQdcCut);  return;}                     //no clusters with QDC more the threshold at all
   if(dMin>pParam->DistCut()) {pTrk->SetHMPIDsignal(kMipDistCut); return;}                     //closest cluster with enough charge is still too far from intersection
   pTrk->SetHMPIDcluIdx(chId,mipId);                                                           //set index of cluster
-  if(iNacc<1)    pTrk->SetHMPIDsignal(kNoPhotAccept);                                         //no photon candidates is accepted
-  else           pTrk->SetHMPIDsignal(FindRingCkov(pCluLst->GetEntries()));                   //find best Theta ckov for ring i.e. track
-
-  pTrk->SetHMPIDchi2(fCkovSigma2);                                                            //errors squared 
+  if(iNacc<1){
+    pTrk->SetHMPIDsignal(kNoPhotAccept);                                                      //no photon candidates is accepted
+  }
+  else {
+    pTrk->SetHMPIDsignal(FindRingCkov(pCluLst->GetEntries()));                                //find best Theta ckov for ring i.e. track
+    pTrk->SetHMPIDchi2(fCkovSigma2);                                                          //errors squared
+  }
 
 }//CkovAngle()
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
@@ -298,7 +303,7 @@ void AliHMPIDRecon::Refract(TVector3 &dir,Double_t n1,Double_t n2)const
 //   Returns: none
 //   On exit: dir is new direction
   Double_t sinref=(n1/n2)*TMath::Sin(dir.Theta());
-  if(sinref>1.)    dir.SetXYZ(-999,-999,-999);
+  if(TMath::Abs(sinref)>1.) dir.SetXYZ(-999,-999,-999);
   else             dir.SetTheta(TMath::ASin(sinref));
 }//Refract()
 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
index 4289dd21cfca1c36a7f01ac40494d3cf411cbb9c..d37497c70fb236f475eb25bc4d3333fb717dec2c 100644 (file)
@@ -672,7 +672,6 @@ void HmpConfig::WriteBatch()
                                                     fprintf(fp,"  cout<<\"!!!!!!!!!!!!Info in <my/Batch.C>: Stop  time: \";time.Set();  time.Print();\n");
                                                     fprintf(fp,"  gBenchmark->Show(\"ALICE\");\n");
   
-                                                    fprintf(fp,"  gSystem->Exec(\"play -c 2 ~/my/end.wav\");\n");
                                                     fprintf(fp,"  gSystem->Exec(\"touch ZZZ______finished_______ZZZ\");\n}\n");
   fclose(fp);  
 }//WriteBatch()