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 **************************************************************************/
16 //************************************************************************
18 // ESD track and V0 fast resolution parameterization
19 // Fast algorithm to make parameterization
20 // Track covariance used
23 // Origin: Marian Ivanov marian.ivanov@cern.ch
24 //-------------------------------------------------------------------------
32 gSystem->Load("libSTAT.so");
33 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
34 gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+");
35 AliXRDPROOFtoolkit tool;
36 TChain * tree = tool.MakeChain("esd.txt","esdTree",0,1000);
41 .L $ALICE_ROOT/PWG1/AliESDresolParams.cxx+
42 .L $ALICE_ROOT/PWG1/AliESDresolMakerFast.cxx+
43 AliESDresolParams params;
44 TCut cutDCA="Tracks[].fCchi2<100&&abs(Tracks[].fP[4])<8&&abs(Tracks[].fP[3])<1&&sqrt(Tracks[].fC[0])/(0.2+abs(Tracks[].fP[4]))<0.02&&abs(Tracks[].fX)<3&&Tracks[].fITSncls>4&&Tracks.fTPCncls>40"
46 TObjArray * array = AliESDresolMakerFast::MakeParamPrimFast(tree,cutDCA,0.95,2000);
54 #include "../STAT/TStatToolkit.h"
59 #include "AliESDresolParams.h"
60 #include "AliESDresolMakerFast.h"
63 ClassImp(AliESDresolMakerFast)
68 AliESDresolMakerFast::AliESDresolMakerFast() :
72 // Default constructor
77 TObjArray * AliESDresolMakerFast::MakeParamPrimFast(TTree * tree, TCut & cutDCA, Float_t fraction, Int_t entries){
79 // DCA resolution linear parameterization
80 // Only valid in ITS acceptance
82 // tree - esdTree or chain
83 // cutDCA - track selection criteria
84 // fraction - robust fitter fraction
85 // entries - total number of entries (tracks) used
89 TCut cutDCA="Tracks[].fCchi2<100&&abs(Tracks[].fP[4])<8&&abs(Tracks[].fP[3])<1&&sqrt(Tracks[].fC[0])/(0.2+abs(Tracks[].fP[4]))<0.02&&abs(Tracks[].fX)<3&&Tracks[].fITSncls>4&&Tracks.fTPCncls>40"
94 TObjArray *array = new TObjArray();
102 TString * dcayParam= TStatToolkit::FitPlane(tree,"sqrt(Tracks[].fC[0])/(0.2+abs(Tracks[].fP[4]))","(abs(Tracks[].fP[4]))++(abs(Tracks[].fP[4]))^2++(abs(Tracks[].fP[3]))++(abs(Tracks[].fP[3]))^2++(abs(Tracks[].fP[4]))*(abs(Tracks[].fP[3]))", cutDCA, chi2,npoints,fitParam,covMatrix,fraction, 0, entries);
103 tree->SetAlias("dcayParam",dcayParam->Data());
104 array->AddAt(new TVectorD(fitParam),0);
105 printf("Y resol\t%s\n",dcayParam->Data());
109 TString * dcazParam= TStatToolkit::FitPlane(tree, "sqrt(Tracks[].fC[2])/(0.2+abs(Tracks[].fP[4]))", "(abs(Tracks[].fP[4]))++(abs(Tracks[].fP[4]))^2++(abs(Tracks[].fP[3]))++(abs(Tracks[].fP[3]))^2++(abs(Tracks[].fP[4]))*(abs(Tracks[].fP[3]))", cutDCA, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
110 printf("Z resol\t%s\n",dcazParam->Data());
111 tree->SetAlias("dcazParam",dcazParam->Data());
112 array->AddAt(new TVectorD(fitParam),1);
116 TString * dcaphiParam= TStatToolkit::FitPlane(tree, "sqrt(Tracks[].fC[5])/(0.1+abs(Tracks[].fP[4]))", "(abs(Tracks[].fP[4]))++(abs(Tracks[].fP[4]))^2++(abs(Tracks[].fP[3]))++(abs(Tracks[].fP[3]))^2++(abs(Tracks[].fP[4]))*(abs(Tracks[].fP[3]))", cutDCA, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
117 printf("Phi resol\t%s\n",dcaphiParam->Data());
118 tree->SetAlias("dcaphiParam",dcaphiParam->Data());
119 array->AddAt(new TVectorD(fitParam),2);
123 TString * dcathParam= TStatToolkit::FitPlane(tree, "sqrt(Tracks[].fC[9])/(0.1+abs(Tracks[].fP[4]))", "(abs(Tracks[].fP[4]))++(abs(Tracks[].fP[4]))^2++(abs(Tracks[].fP[3]))++(abs(Tracks[].fP[3]))^2++(abs(Tracks[].fP[4]))*(abs(Tracks[].fP[3]))", cutDCA, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
124 printf("Theta resol\t%s\n",dcathParam->Data());
125 tree->SetAlias("dcathParam",dcathParam->Data());
126 array->AddAt(new TVectorD(fitParam),3);
130 TString * dca1ptParam= TStatToolkit::FitPlane(tree, "sqrt(Tracks[].fC[14])/(1+abs(Tracks[].fP[4]))^2", "(abs(Tracks[].fP[4]))++(abs(Tracks[].fP[4]))^2++(abs(Tracks[].fP[3]))++(abs(Tracks[].fP[3]))^2++(abs(Tracks[].fP[4]))*(abs(Tracks[].fP[3]))", cutDCA, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
131 printf("1pt resol\t%s\n",dca1ptParam->Data());
132 tree->SetAlias("dca1ptParam",dca1ptParam->Data());
133 array->AddAt(new TVectorD(fitParam),4);
139 TObjArray * AliESDresolMakerFast::MakeParamRFast(TTree * tree, TCut &cutV0, Float_t fraction, Int_t entries){
150 TCut cutV0 = "abs((V0s[].fParamP.fC[0]))<3&&abs(V0s[].fParamP.fP[3])<1&&abs(V0s[].fParamP.fP[4])<8";//
152 TObjArray *array = new TObjArray;
156 TString * v0sigmaY= TStatToolkit::FitPlane(tree,"sqrt(sqrt((V0s[].fParamP.fC[0])))","abs(V0s[].fParamP.fP[4])++V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])^2++V0s[].fParamP.fX^2++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX^2", cutV0, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
157 tree->SetAlias("v0sigmaY",v0sigmaY->Data());
158 array->AddLast(new TVectorD(fitParam));
162 TString * v0sigmaZ= TStatToolkit::FitPlane(tree,"sqrt(sqrt((V0s[].fParamP.fC[2])))","abs(V0s[].fParamP.fP[4])++V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])^2++V0s[].fParamP.fX^2++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX^2", cutV0, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
163 tree->SetAlias("v0sigmaZ",v0sigmaZ->Data());
164 array->AddLast(new TVectorD(fitParam));
168 TString * v0sigmaPhi= TStatToolkit::FitPlane(tree,"sqrt(sqrt((V0s[].fParamP.fC[5])))","abs(V0s[].fParamP.fP[4])++V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])^2++V0s[].fParamP.fX^2++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX^2", cutV0, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
169 tree->SetAlias("v0sigmaPhi",v0sigmaPhi->Data());
170 array->AddLast(new TVectorD(fitParam));
174 TString * v0sigmaTh= TStatToolkit::FitPlane(tree,"sqrt(sqrt((V0s[].fParamP.fC[9])))","abs(V0s[].fParamP.fP[4])++V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])^2++V0s[].fParamP.fX^2++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX^2", cutV0, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
175 tree->SetAlias("v0sigmaTh",v0sigmaTh->Data());
176 array->AddLast(new TVectorD(fitParam));
180 TString * v0sigma1pt= TStatToolkit::FitPlane(tree,"sqrt(sqrt((V0s[].fParamP.fC[14])))","abs(V0s[].fParamP.fP[4])++V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])^2++V0s[].fParamP.fX^2++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX++abs(V0s[].fParamP.fP[4])*V0s[].fParamP.fX^2", cutV0, chi2,npoints,fitParam,covMatrix,fraction,0,entries);
181 tree->SetAlias("v0sigma1pt",v0sigma1pt->Data());
182 array->AddLast(new TVectorD(fitParam));