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 //----------------------------------------------------------------------------
19 // Implementation of the class to calculate the parton energy loss
20 // Based on the "BDMPS" quenching weights by C.A.Salgado and U.A.Wiedemann
23 // C.A.Salgado and U.A.Wiedemann, Phys.Rev.D68 (2003) 014008 [hep-ph/0302184]
24 // A.Dainese, Eur.Phys.J.C, in press, [nucl-ex/0312005]
27 // Origin: C. Loizides constantinos.loizides@cern.ch
28 // A. Dainese andrea.dainese@pd.infn.it
30 //=================== Added by C. Loizides 27/03/04 ===========================
32 // Added support for k-Quenching, where wc=I1*k and R=2I1^2/I0*k
33 // (see the AliFastGlauber class for definition of I0/I1)
34 //-----------------------------------------------------------------------------
36 #include <Riostream.h>
45 #include "AliQuenchingWeights.h"
47 ClassImp(AliQuenchingWeights)
49 // conversion from fm to GeV^-1: 1 fm = fmGeV GeV^-1
50 const Double_t AliQuenchingWeights::fgkConvFmToInvGeV = 1./0.197;
53 const Double_t AliQuenchingWeights::fgkRMax = 1.e6;
56 const Int_t AliQuenchingWeights::fgkBins = 1300;
57 const Double_t AliQuenchingWeights::fgkMaxBin = 1.3;
59 // counter for histogram labels
60 Int_t AliQuenchingWeights::fgCounter = 0;
63 AliQuenchingWeights::AliQuenchingWeights()
78 fInstanceNumber=fgCounter++;
80 sprintf(name,"hhistoqw_%d",fInstanceNumber);
81 fHisto = new TH1F(name,"",fgkBins,0.,fgkMaxBin);
82 for(Int_t bin=1;bin<=fgkBins;bin++)
83 fHisto->SetBinContent(bin,0.);
86 AliQuenchingWeights::AliQuenchingWeights(const AliQuenchingWeights& a)
94 fMultSoft=a.GetMultSoft();;
97 fQTransport=a.GetQTransport();
98 fECMethod=(kECMethod)a.GetECMethod();
99 fLengthMax=a.GetLengthMax();
100 fInstanceNumber=fgCounter++;
102 sprintf(name,"hhistoqw_%d",fInstanceNumber);
103 fHisto = new TH1F(name,"",fgkBins,0.,fgkMaxBin);
104 for(Int_t bin=1;bin<=fgkBins;bin++)
105 fHisto->SetBinContent(bin,0.);
107 //Missing in the class is the pathname
108 //to the tables, can be added if needed
111 AliQuenchingWeights::~AliQuenchingWeights()
117 void AliQuenchingWeights::Reset()
119 //reset tables if there were used
122 for(Int_t l=0;l<4*fLengthMaxOld;l++){
123 delete fHistos[0][l];
124 delete fHistos[1][l];
131 void AliQuenchingWeights::SetECMethod(kECMethod type)
133 //set energy constraint method
136 if(fECMethod==kDefault)
137 Info("SetECMethod","Energy Constraint Method set to DEFAULT:\nIf (sampled energy loss > parton energy) then sampled energy loss = parton energy.");
139 Info("SetECMethod","Energy Constraint Method set to REWEIGHT:\nRequire sampled energy loss <= parton energy.");
142 Int_t AliQuenchingWeights::InitMult(const Char_t *contall,const Char_t *discall)
144 // read in tables for multiple scattering approximation
145 // path to continuum and to discrete part
147 fTablesLoaded = kFALSE;
151 sprintf(fname,"%s",gSystem->ExpandPathName(contall));
152 //PH ifstream fincont(fname);
153 fstream fincont(fname,ios::in);
154 #if defined(__HP_aCC) || defined(__DECCXX)
155 if(!fincont.rdbuf()->is_open()) return -1;
157 if(!fincont.is_open()) return -1;
161 while(fincont>>fxx[nn]>>fcaq[0][nn]>>fcaq[1][nn]>>fcaq[2][nn]>>fcaq[3][nn]>>
162 fcaq[4][nn]>>fcaq[5][nn]>>fcaq[6][nn]>>fcaq[7][nn]>>fcaq[8][nn]>>
163 fcaq[9][nn]>>fcaq[10][nn]>>fcaq[11][nn]>>fcaq[12][nn]>>fcaq[13][nn]>>
164 fcaq[14][nn]>>fcaq[15][nn]>>fcaq[16][nn]>>fcaq[17][nn]>>fcaq[18][nn]>>
165 fcaq[19][nn]>>fcaq[20][nn]>>fcaq[21][nn]>>fcaq[22][nn]>>fcaq[23][nn]>>
166 fcaq[24][nn]>>fcaq[25][nn]>>fcaq[26][nn]>>fcaq[27][nn]>>fcaq[28][nn]>>
167 fcaq[29][nn]>>fcaq[30][nn]>>fcaq[31][nn]>>fcaq[32][nn]>>fcaq[33][nn])
174 while(fincont>>fxxg[nn]>>fcag[0][nn]>>fcag[1][nn]>>fcag[2][nn]>>fcag[3][nn]>>
175 fcag[4][nn]>>fcag[5][nn]>>fcag[6][nn]>>fcag[7][nn]>>fcag[8][nn]>>
176 fcag[9][nn]>>fcag[10][nn]>>fcag[11][nn]>>fcag[12][nn]>>fcag[13][nn]>>
177 fcag[14][nn]>>fcag[15][nn]>>fcag[16][nn]>>fcag[17][nn]>>fcag[18][nn]>>
178 fcag[19][nn]>>fcag[20][nn]>>fcag[21][nn]>>fcag[22][nn]>>fcag[23][nn]>>
179 fcag[24][nn]>>fcag[25][nn]>>fcag[26][nn]>>fcag[27][nn]>>fcag[28][nn]>>
180 fcag[29][nn]>>fcag[30][nn]>>fcag[31][nn]>>fcag[32][nn]>>fcag[33][nn])
187 sprintf(fname,"%s",gSystem->ExpandPathName(discall));
188 //PH ifstream findisc(fname);
189 fstream findisc(fname,ios::in);
190 #if defined(__HP_aCC) || defined(__DECCXX)
191 if(!findisc.rdbuf()->is_open()) return -1;
193 if(!findisc.is_open()) return -1;
197 while(findisc>>frrr[nn]>>fdaq[nn]) {
202 while(findisc>>frrrg[nn]>>fdag[nn]) {
207 fTablesLoaded = kTRUE;
212 C***************************************************************************
213 C Quenching Weights for Multiple Soft Scattering
218 C Carlos A. Salgado and Urs A. Wiedemann, hep-ph/0302184.
220 C Carlos A. Salgado and Urs A. Wiedemann Phys.Rev.Lett.89:092303,2002.
223 C This package contains quenching weights for gluon radiation in the
224 C multiple soft scattering approximation.
226 C swqmult returns the quenching weight for a quark (ipart=1) or
227 C a gluon (ipart=2) traversing a medium with transport coeficient q and
228 C length L. The input values are rrrr=0.5*q*L^3 and xxxx=w/wc, where
229 C wc=0.5*q*L^2 and w is the energy radiated. The output values are
230 C the continuous and discrete (prefactor of the delta function) parts
231 C of the quenching weights.
233 C In order to use this routine, the files cont.all and disc.all need to be
234 C in the working directory.
236 C An initialization of the tables is needed by doing call initmult before
239 C Please, send us any comment:
241 C urs.wiedemann@cern.ch
242 C carlos.salgado@cern.ch
245 C-------------------------------------------------------------------
247 SUBROUTINE swqmult(ipart,rrrr,xxxx,continuous,discrete)
249 REAL*8 xx(400), daq(34), caq(34,261), rrr(34)
250 COMMON /dataqua/ xx, daq, caq, rrr
252 REAL*8 xxg(400), dag(34), cag(34,261), rrrg(34)
253 COMMON /dataglu/ xxg, dag, cag, rrrg
255 REAL*8 rrrr,xxxx, continuous, discrete
257 INTEGER nrlow, nrhigh, nxlow, nxhigh
258 REAL*8 rrhigh, rrlow, rfraclow, rfrachigh
259 REAL*8 xfraclow, xfrachigh
270 if (rrin.lt.rrr(nr)) then
282 rfraclow = (rrhigh-rrin)/(rrhigh-rrlow)
283 rfrachigh = (rrin-rrlow)/(rrhigh-rrlow)
284 if (rrin.gt.10000d0) then
285 rfraclow = dlog(rrhigh/rrin)/dlog(rrhigh/rrlow)
286 rfrachigh = dlog(rrin/rrlow)/dlog(rrhigh/rrlow)
289 if (ipart.eq.1.and.rrin.ge.rrr(1)) then
296 if (ipart.ne.1.and.rrin.ge.rrrg(1)) then
303 if (xxxx.ge.xx(261)) go to 245
305 nxlow = int(xxin/0.01) + 1
307 xfraclow = (xx(nxhigh)-xxin)/0.01
308 xfrachigh = (xxin - xx(nxlow))/0.01
311 clow = xfraclow*caq(nrlow,nxlow)+xfrachigh*caq(nrlow,nxhigh)
312 chigh = xfraclow*caq(nrhigh,nxlow)+xfrachigh*caq(nrhigh,nxhigh)
314 clow = xfraclow*cag(nrlow,nxlow)+xfrachigh*cag(nrlow,nxhigh)
315 chigh = xfraclow*cag(nrhigh,nxlow)+xfrachigh*cag(nrhigh,nxhigh)
318 continuous = rfraclow*clow + rfrachigh*chigh
323 discrete = rfraclow*daq(nrlow) + rfrachigh*daq(nrhigh)
325 discrete = rfraclow*dag(nrlow) + rfrachigh*dag(nrhigh)
331 REAL*8 xxq(400), daq(34), caq(34,261), rrr(34)
332 COMMON /dataqua/ xxq, daq, caq, rrr
334 REAL*8 xxg(400), dag(34), cag(34,261), rrrg(34)
335 COMMON /dataglu/ xxg, dag, cag, rrrg
337 OPEN(UNIT=20,FILE='contnew.all',STATUS='OLD',ERR=90)
339 read (20,*) xxq(nn), caq(1,nn), caq(2,nn), caq(3,nn),
340 + caq(4,nn), caq(5,nn), caq(6,nn), caq(7,nn), caq(8,nn),
341 + caq(9,nn), caq(10,nn), caq(11,nn), caq(12,nn),
343 + caq(14,nn), caq(15,nn), caq(16,nn), caq(17,nn),
345 + caq(19,nn), caq(20,nn), caq(21,nn), caq(22,nn),
347 + caq(24,nn), caq(25,nn), caq(26,nn), caq(27,nn),
349 + caq(29,nn), caq(30,nn), caq(31,nn), caq(32,nn),
350 + caq(33,nn), caq(34,nn)
353 read (20,*) xxg(nn), cag(1,nn), cag(2,nn), cag(3,nn),
354 + cag(4,nn), cag(5,nn), cag(6,nn), cag(7,nn), cag(8,nn),
355 + cag(9,nn), cag(10,nn), cag(11,nn), cag(12,nn),
357 + cag(14,nn), cag(15,nn), cag(16,nn), cag(17,nn),
359 + cag(19,nn), cag(20,nn), cag(21,nn), cag(22,nn),
361 + cag(24,nn), cag(25,nn), cag(26,nn), cag(27,nn),
363 + cag(29,nn), cag(30,nn), cag(31,nn), cag(32,nn),
364 + cag(33,nn), cag(34,nn)
368 OPEN(UNIT=21,FILE='discnew.all',STATUS='OLD',ERR=91)
370 read (21,*) rrr(nn), daq(nn)
373 read (21,*) rrrg(nn), dag(nn)
378 90 PRINT*, 'input - output error'
379 91 PRINT*, 'input - output error #2'
385 =======================================================================
387 Adapted to ROOT macro by A. Dainese - 13/07/2003
388 Ported to class by C. Loizides - 12/02/2004
389 New version for extended R values added - 06/03/2004
392 Int_t AliQuenchingWeights::CalcMult(Int_t ipart, Double_t rrrr,Double_t xxxx,
393 Double_t &continuous,Double_t &discrete) const
395 // Calculate Multiple Scattering approx.
396 // weights for given parton type,
397 // rrrr=0.5*q*L^3 and xxxx=w/wc, wc=0.5*q*L^2
403 //read-in data before first call
405 Error("CalcMult","Tables are not loaded.");
409 Error("CalcMult","Tables are not loaded for Multiple Scattering.");
413 Double_t rrin = rrrr;
414 Double_t xxin = xxxx;
416 if(xxin>fxx[260]) return -1;
417 Int_t nxlow = (Int_t)(xxin/0.01) + 1;
418 Int_t nxhigh = nxlow + 1;
419 Double_t xfraclow = (fxx[nxhigh-1]-xxin)/0.01;
420 Double_t xfrachigh = (xxin - fxx[nxlow-1])/0.01;
423 if(rrin<=frrr[33]) rrin = 1.05*frrr[33]; // AD
424 if(rrin>=frrr[0]) rrin = 0.95*frrr[0]; // AD
426 Int_t nrlow=0,nrhigh=0;
427 Double_t rrhigh=0,rrlow=0;
428 for(Int_t nr=1; nr<=34; nr++) {
429 if(rrin<frrr[nr-1]) {
432 rrhigh = frrr[nr-1-1];
442 Double_t rfraclow = (rrhigh-rrin)/(rrhigh-rrlow);
443 Double_t rfrachigh = (rrin-rrlow)/(rrhigh-rrlow);
446 rfraclow = TMath::Log2(rrhigh/rrin)/TMath::Log2(rrhigh/rrlow);
447 rfrachigh = TMath::Log2(rrin/rrlow)/TMath::Log2(rrhigh/rrlow);
449 if((ipart==1) && (rrin>=frrr[0]))
456 if((ipart==2) && (rrin>=frrrg[0]))
464 //printf("R = %f,\nRlow = %f, Rhigh = %f,\nRfraclow = %f, Rfrachigh = %f\n",rrin,rrlow,rrhigh,rfraclow,rfrachigh); // AD
466 Double_t clow=0,chigh=0;
468 clow = xfraclow*fcaq[nrlow-1][nxlow-1]+xfrachigh*fcaq[nrlow-1][nxhigh-1];
469 chigh = xfraclow*fcaq[nrhigh-1][nxlow-1]+xfrachigh*fcaq[nrhigh-1][nxhigh-1];
471 clow = xfraclow*fcag[nrlow-1][nxlow-1]+xfrachigh*fcag[nrlow-1][nxhigh-1];
472 chigh = xfraclow*fcag[nrhigh-1][nxlow-1]+xfrachigh*fcag[nrhigh-1][nxhigh-1];
475 continuous = rfraclow*clow + rfrachigh*chigh;
476 //printf("rfraclow %f, clow %f, rfrachigh %f, chigh %f,\n continuous %f\n",
477 //rfraclow,clow,rfrachigh,chigh,continuous);
480 discrete = rfraclow*fdaq[nrlow-1] + rfrachigh*fdaq[nrhigh-1];
482 discrete = rfraclow*fdag[nrlow-1] + rfrachigh*fdag[nrhigh-1];
488 Int_t AliQuenchingWeights::InitSingleHard(const Char_t *contall,const Char_t *discall)
490 // read in tables for Single Hard Approx.
491 // path to continuum and to discrete part
493 fTablesLoaded = kFALSE;
497 sprintf(fname,"%s",gSystem->ExpandPathName(contall));
498 //PH ifstream fincont(fname);
499 fstream fincont(fname,ios::in);
500 #if defined(__HP_aCC) || defined(__DECCXX)
501 if(!fincont.rdbuf()->is_open()) return -1;
503 if(!fincont.is_open()) return -1;
507 while(fincont>>fxx[nn]>>fcaq[0][nn]>>fcaq[1][nn]>>fcaq[2][nn]>>fcaq[3][nn]>>
508 fcaq[4][nn]>>fcaq[5][nn]>>fcaq[6][nn]>>fcaq[7][nn]>>fcaq[8][nn]>>
509 fcaq[9][nn]>>fcaq[10][nn]>>fcaq[11][nn]>>fcaq[12][nn]>>
511 fcaq[14][nn]>>fcaq[15][nn]>>fcaq[16][nn]>>fcaq[17][nn]>>
513 fcaq[19][nn]>>fcaq[20][nn]>>fcaq[21][nn]>>fcaq[22][nn]>>
515 fcaq[24][nn]>>fcaq[25][nn]>>fcaq[26][nn]>>fcaq[27][nn]>>
524 while(fincont>>fxxg[nn]>>fcag[0][nn]>>fcag[1][nn]>>fcag[2][nn]>>fcag[3][nn]>>
525 fcag[4][nn]>>fcag[5][nn]>>fcag[6][nn]>>fcag[7][nn]>>fcag[8][nn]>>
526 fcag[9][nn]>>fcag[10][nn]>>fcag[11][nn]>>fcag[12][nn]>>
528 fcag[14][nn]>>fcag[15][nn]>>fcag[16][nn]>>fcag[17][nn]>>
530 fcag[19][nn]>>fcag[20][nn]>>fcag[21][nn]>>fcag[22][nn]>>
532 fcag[24][nn]>>fcag[25][nn]>>fcag[26][nn]>>fcag[27][nn]>>
540 sprintf(fname,"%s",gSystem->ExpandPathName(discall));
541 //PH ifstream findisc(fname);
542 fstream findisc(fname,ios::in);
543 #if defined(__HP_aCC) || defined(__DECCXX)
544 if(!findisc.rdbuf()->is_open()) return -1;
546 if(!findisc.is_open()) return -1;
550 while(findisc>>frrr[nn]>>fdaq[nn]) {
555 while(findisc>>frrrg[nn]>>fdag[nn]) {
561 fTablesLoaded = kTRUE;
566 C***************************************************************************
567 C Quenching Weights for Single Hard Scattering
572 C Carlos A. Salgado and Urs A. Wiedemann, hep-ph/0302184.
574 C Carlos A. Salgado and Urs A. Wiedemann Phys.Rev.Lett.89:092303,2002.
577 C This package contains quenching weights for gluon radiation in the
578 C single hard scattering approximation.
580 C swqlin returns the quenching weight for a quark (ipart=1) or
581 C a gluon (ipart=2) traversing a medium with Debye screening mass mu and
582 C length L. The input values are rrrr=0.5*mu^2*L^2 and xxxx=w/wc, where
583 C wc=0.5*mu^2*L and w is the energy radiated. The output values are
584 C the continuous and discrete (prefactor of the delta function) parts
585 C of the quenching weights.
587 C In order to use this routine, the files contlin.all and disclin.all
588 C need to be in the working directory.
590 C An initialization of the tables is needed by doing call initlin before
593 C Please, send us any comment:
595 C urs.wiedemann@cern.ch
596 C carlos.salgado@cern.ch
599 C-------------------------------------------------------------------
602 SUBROUTINE swqlin(ipart,rrrr,xxxx,continuous,discrete)
604 REAL*8 xx(400), dalq(30), calq(30,261), rrr(30)
605 COMMON /datalinqua/ xx, dalq, calq, rrr
607 REAL*8 xxlg(400), dalg(30), calg(30,261), rrrlg(30)
608 COMMON /datalinglu/ xxlg, dalg, calg, rrrlg
610 REAL*8 rrrr,xxxx, continuous, discrete
612 INTEGER nrlow, nrhigh, nxlow, nxhigh
613 REAL*8 rrhigh, rrlow, rfraclow, rfrachigh
614 REAL*8 xfraclow, xfrachigh
620 nxlow = int(xxin/0.038) + 1
622 xfraclow = (xx(nxhigh)-xxin)/0.038
623 xfrachigh = (xxin - xx(nxlow))/0.038
626 if (rrin.lt.rrr(nr)) then
638 rfraclow = (rrhigh-rrin)/(rrhigh-rrlow)
639 rfrachigh = (rrin-rrlow)/(rrhigh-rrlow)
642 clow = xfraclow*calq(nrlow,nxlow)+xfrachigh*calq(nrlow,nxhigh)
643 chigh = xfraclow*calq(nrhigh,nxlow)+xfrachigh*calq(nrhigh,nxhigh)
645 clow = xfraclow*calg(nrlow,nxlow)+xfrachigh*calg(nrlow,nxhigh)
646 chigh = xfraclow*calg(nrhigh,nxlow)+xfrachigh*calg(nrhigh,nxhigh)
649 continuous = rfraclow*clow + rfrachigh*chigh
652 discrete = rfraclow*dalq(nrlow) + rfrachigh*dalq(nrhigh)
654 discrete = rfraclow*dalg(nrlow) + rfrachigh*dalg(nrhigh)
660 REAL*8 xxlq(400), dalq(30), calq(30,261), rrr(30)
661 COMMON /datalinqua/ xxlq, dalq, calq, rrr
663 REAL*8 xxlg(400), dalg(30), calg(30,261), rrrlg(30)
664 COMMON /datalinglu/ xxlg, dalg, calg, rrrlg
666 OPEN(UNIT=20,FILE='contlin.all',STATUS='OLD',ERR=90)
668 read (20,*) xxlq(nn), calq(1,nn), calq(2,nn), calq(3,nn),
669 + calq(4,nn), calq(5,nn), calq(6,nn), calq(7,nn), calq(8,nn),
670 + calq(9,nn), calq(10,nn), calq(11,nn), calq(12,nn),
672 + calq(14,nn), calq(15,nn), calq(16,nn), calq(17,nn),
674 + calq(19,nn), calq(20,nn), calq(21,nn), calq(22,nn),
676 + calq(24,nn), calq(25,nn), calq(26,nn), calq(27,nn),
678 + calq(29,nn), calq(30,nn)
681 read (20,*) xxlg(nn), calg(1,nn), calg(2,nn), calg(3,nn),
682 + calg(4,nn), calg(5,nn), calg(6,nn), calg(7,nn), calg(8,nn),
683 + calg(9,nn), calg(10,nn), calg(11,nn), calg(12,nn),
685 + calg(14,nn), calg(15,nn), calg(16,nn), calg(17,nn),
687 + calg(19,nn), calg(20,nn), calg(21,nn), calg(22,nn),
689 + calg(24,nn), calg(25,nn), calg(26,nn), calg(27,nn),
691 + calg(29,nn), calg(30,nn)
695 OPEN(UNIT=21,FILE='disclin.all',STATUS='OLD',ERR=91)
697 read (21,*) rrr(nn), dalq(nn)
700 read (21,*) rrrlg(nn), dalg(nn)
705 90 PRINT*, 'input - output error'
706 91 PRINT*, 'input - output error #2'
711 =======================================================================
713 Ported to class by C. Loizides - 17/02/2004
717 Int_t AliQuenchingWeights::CalcSingleHard(Int_t ipart, Double_t rrrr,Double_t xxxx,
718 Double_t &continuous,Double_t &discrete) const
720 // calculate Single Hard approx.
721 // weights for given parton type,
722 // rrrr=0.5*mu^2*L^2 and xxxx=w/wc, wc=0.5*mu^2*L
724 // read-in data before first call
726 Error("CalcSingleHard","Tables are not loaded.");
730 Error("CalcSingleHard","Tables are not loaded for Single Hard Scattering.");
734 Double_t rrin = rrrr;
735 Double_t xxin = xxxx;
737 Int_t nxlow = (Int_t)(xxin/0.038) + 1;
738 Int_t nxhigh = nxlow + 1;
739 Double_t xfraclow = (fxx[nxhigh-1]-xxin)/0.038;
740 Double_t xfrachigh = (xxin - fxx[nxlow-1])/0.038;
743 if(rrin<=frrr[29]) rrin = 1.05*frrr[29]; // AD
744 if(rrin>=frrr[0]) rrin = 0.95*frrr[0]; // AD
746 Int_t nrlow=0,nrhigh=0;
747 Double_t rrhigh=0,rrlow=0;
748 for(Int_t nr=1; nr<=30; nr++) {
749 if(rrin<frrr[nr-1]) {
752 rrhigh = frrr[nr-1-1];
762 Double_t rfraclow = (rrhigh-rrin)/(rrhigh-rrlow);
763 Double_t rfrachigh = (rrin-rrlow)/(rrhigh-rrlow);
765 //printf("R = %f,\nRlow = %f, Rhigh = %f,\nRfraclow = %f, Rfrachigh = %f\n",rrin,rrlow,rrhigh,rfraclow,rfrachigh); // AD
767 Double_t clow=0,chigh=0;
769 clow = xfraclow*fcaq[nrlow-1][nxlow-1]+xfrachigh*fcaq[nrlow-1][nxhigh-1];
770 chigh = xfraclow*fcaq[nrhigh-1][nxlow-1]+xfrachigh*fcaq[nrhigh-1][nxhigh-1];
772 clow = xfraclow*fcag[nrlow-1][nxlow-1]+xfrachigh*fcag[nrlow-1][nxhigh-1];
773 chigh = xfraclow*fcag[nrhigh-1][nxlow-1]+xfrachigh*fcag[nrhigh-1][nxhigh-1];
776 continuous = rfraclow*clow + rfrachigh*chigh;
777 //printf("rfraclow %f, clow %f, rfrachigh %f, chigh %f,\n continuous %f\n",
778 // rfraclow,clow,rfrachigh,chigh,continuous);
781 discrete = rfraclow*fdaq[nrlow-1] + rfrachigh*fdaq[nrhigh-1];
783 discrete = rfraclow*fdag[nrlow-1] + rfrachigh*fdag[nrhigh-1];
789 Int_t AliQuenchingWeights::CalcMult(Int_t ipart,
790 Double_t w,Double_t qtransp,Double_t length,
791 Double_t &continuous,Double_t &discrete) const
793 // Calculate Multiple Scattering approx.
794 // weights for given parton type,
795 // rrrr=0.5*q*L^3 and xxxx=w/wc, wc=0.5*q*L^2
797 Double_t wc=CalcWC(qtransp,length);
798 Double_t rrrr=CalcR(wc,length);
800 return CalcMult(ipart,rrrr,xxxx,continuous,discrete);
803 Int_t AliQuenchingWeights::CalcSingleHard(Int_t ipart,
804 Double_t w,Double_t mu,Double_t length,
805 Double_t &continuous,Double_t &discrete) const
807 // calculate Single Hard approx.
808 // weights for given parton type,
809 // rrrr=0.5*mu^2*L^2 and xxxx=w/wc, wc=0.5*mu^2*L
811 Double_t wcbar=CalcWCbar(mu,length);
812 Double_t rrrr=CalcR(wcbar,length);
813 Double_t xxxx=w/wcbar;
814 return CalcSingleHard(ipart,rrrr,xxxx,continuous,discrete);
817 Double_t AliQuenchingWeights::CalcR(Double_t wc, Double_t l) const
819 //calculate R value and
820 //check if it is less then maximum
822 Double_t R = wc*l*fgkConvFmToInvGeV;
824 Warning("CalcR","Value of R = %.2f; should be less than %.2f",R,fgkRMax);
830 Double_t AliQuenchingWeights::CalcRk(Double_t k, Double_t I0, Double_t I1) const
832 //calculate R value and
833 //check if it is less then maximum
835 Double_t R = fgkRMax-1;
839 Warning("CalcRk","Value of R = %.2f; should be less than %.2f",R,fgkRMax);
845 Double_t AliQuenchingWeights::GetELossRandom(Int_t ipart, Double_t length, Double_t e) const
847 // return DeltaE for MS or SH scattering
848 // for given parton type, length and energy
849 // Dependant on ECM (energy constraint method)
850 // e is used to determine where to set bins to zero.
853 Fatal("GetELossRandom","Call SampleEnergyLoss method before!");
856 if((ipart<1) || (ipart>2)) {
857 Fatal("GetELossRandom","ipart =%d; but has to be 1 (quark) or 2 (gluon)",ipart);
861 Int_t l=GetIndex(length);
864 if(fECMethod==kReweight){
868 ret=fHistos[ipart-1][l-1]->GetRandom();
870 Warning("GetELossRandom",
871 "Aborted reweighting; maximum loss assigned after 1e6 trials.");
878 Double_t ret=fHistos[ipart-1][l-1]->GetRandom();
883 Double_t AliQuenchingWeights::CalcQuenchedEnergy(Int_t ipart, Double_t length, Double_t e) const
885 //return quenched parton energy
886 //for given parton type, length and energy
888 Double_t loss=GetELossRandom(ipart,length,e);
892 Double_t AliQuenchingWeights::GetELossRandom(Int_t ipart, TH1F *hell, Double_t e) const
894 // return DeltaE for MS or SH scattering
895 // for given parton type, length distribution and energy
896 // Dependant on ECM (energy constraint method)
897 // e is used to determine where to set bins to zero.
900 Warning("GetELossRandom","Pointer to length distribution is NULL.");
903 Double_t ell=hell->GetRandom();
904 return GetELossRandom(ipart,ell,e);
907 Double_t AliQuenchingWeights::CalcQuenchedEnergy(Int_t ipart, TH1F *hell, Double_t e) const
909 //return quenched parton energy
910 //for given parton type, length distribution and energy
912 Double_t loss=GetELossRandom(ipart,hell,e);
916 Double_t AliQuenchingWeights::GetELossRandomK(Int_t ipart, Double_t I0, Double_t I1, Double_t e)
918 // return DeltaE for new dynamic version
919 // for given parton type, I0 and I1 value and energy
920 // Dependant on ECM (energy constraint method)
921 // e is used to determine where to set bins to zero.
923 // read-in data before first call
925 Fatal("GetELossRandomK","Tables are not loaded.");
928 if((ipart<1) || (ipart>2)) {
929 Fatal("GetELossRandomK","ipart =%d; but has to be 1 (quark) or 2 (gluon)",ipart);
933 Double_t R=CalcRk(I0,I1);
935 Fatal("GetELossRandomK","R should not be negative");
938 Double_t wc=CalcWCk(I1);
940 Fatal("GetELossRandomK","wc should be greater than zero");
943 if(SampleEnergyLoss(ipart,R)!=0){
944 Fatal("GetELossRandomK","Could not sample energy loss");
948 if(fECMethod==kReweight){
952 ret=fHisto->GetRandom();
954 Warning("GetELossRandomK",
955 "Aborted reweighting; maximum loss assigned after 1e6 trials.");
963 Double_t ret=fHisto->GetRandom()*wc;
968 Double_t AliQuenchingWeights::CalcQuenchedEnergyK(Int_t ipart, Double_t I0, Double_t I1, Double_t e)
970 //return quenched parton energy
971 //for given parton type, I0 and I1 value and energy
973 Double_t loss=GetELossRandomK(ipart,I0,I1,e);
977 Double_t AliQuenchingWeights::GetELossRandomKFast(Int_t ipart, Double_t I0, Double_t I1, Double_t e)
979 // return DeltaE for new dynamic version
980 // for given parton type, I0 and I1 value and energy
981 // Dependant on ECM (energy constraint method)
982 // e is used to determine where to set bins to zero.
983 // method is optimized and should only be used if
984 // all parameters are well within the bounds.
985 // read-in data tables before first call
987 Double_t R=CalcRk(I0,I1);
992 Double_t wc=CalcWCk(I1);
997 return GetELossRandomKFastR(ipart,R,wc,e);
1000 Double_t AliQuenchingWeights::GetELossRandomKFastR(Int_t ipart, Double_t R, Double_t wc, Double_t e)
1002 // return DeltaE for new dynamic version
1003 // for given parton type, R and wc value and energy
1004 // Dependant on ECM (energy constraint method)
1005 // e is used to determine where to set bins to zero.
1006 // method is optimized and should only be used if
1007 // all parameters are well within the bounds.
1008 // read-in data tables before first call
1014 Double_t discrete=0.;
1015 Double_t continuous=0.;
1017 Double_t xxxx = fHisto->GetBinCenter(bin);
1019 CalcMult(ipart,R,xxxx,continuous,discrete);
1021 CalcSingleHard(ipart,R,xxxx,continuous,discrete);
1023 if(discrete>0.999) {
1024 return 0.; //no energy loss
1027 fHisto->SetBinContent(bin,continuous);
1028 const Int_t kbinmax=fHisto->FindBin(e/wc);
1029 if(kbinmax==1) return e; //maximum energy loss
1032 for(Int_t bin=2; bin<=kbinmax; bin++) {
1033 xxxx = fHisto->GetBinCenter(bin);
1034 CalcMult(ipart,R,xxxx,continuous,discrete);
1035 fHisto->SetBinContent(bin,continuous);
1038 for(Int_t bin=2; bin<=kbinmax; bin++) {
1039 xxxx = fHisto->GetBinCenter(bin);
1040 CalcSingleHard(ipart,R,xxxx,continuous,discrete);
1041 fHisto->SetBinContent(bin,continuous);
1045 const Double_t kdelta=fHisto->Integral(1,kbinmax);
1046 Double_t val=discrete*fgkBins/fgkMaxBin;
1047 fHisto->Fill(0.,val);
1049 if(fECMethod==kReweight)
1050 fHisto->SetBinContent(kbinmax+1,0);
1052 fHisto->SetBinContent(kbinmax+1,(1-discrete)*fgkBins/fgkMaxBin-kdelta);
1054 for(Int_t bin=kbinmax+2; bin<=fgkBins; bin++) {
1055 fHisto->SetBinContent(bin,0);
1058 Double_t ret=fHisto->GetRandom()*wc;
1063 Double_t AliQuenchingWeights::CalcQuenchedEnergyKFast(Int_t ipart, Double_t I0, Double_t I1, Double_t e)
1065 //return quenched parton energy (fast method)
1066 //for given parton type, I0 and I1 value and energy
1068 Double_t loss=GetELossRandomKFast(ipart,I0,I1,e);
1072 Int_t AliQuenchingWeights::SampleEnergyLoss()
1074 // Has to be called to fill the histograms
1076 // For stored values fQTransport loop over
1077 // particle type and length = 1 to fMaxLength (fm)
1078 // to fill energy loss histos
1080 // Take histogram of continuous weights
1081 // Take discrete_weight
1082 // If discrete_weight > 1, put all channels to 0, except channel 1
1083 // Fill channel 1 with discrete_weight/(1-discrete_weight)*integral
1085 // read-in data before first call
1087 Error("SampleEnergyLoss","Tables are not loaded.");
1092 Int_t lmax=CalcLengthMax(fQTransport);
1093 if(fLengthMax>lmax){
1094 Info("SampleEnergyLoss","Maximum length changed from %d to %d;\nin order to have R < %.f",fLengthMax,lmax,fgkRMax);
1098 Warning("SampleEnergyLoss","Maximum length not checked,\nbecause SingeHard is not yet tested.");
1102 fHistos=new TH1F**[2];
1103 fHistos[0]=new TH1F*[4*fLengthMax];
1104 fHistos[1]=new TH1F*[4*fLengthMax];
1105 fLengthMaxOld=fLengthMax; //remember old value in case
1106 //user wants to reset
1109 Char_t meddesc[100];
1111 medvalue=(Int_t)(fQTransport*1000.);
1112 sprintf(meddesc,"MS");
1114 medvalue=(Int_t)(fMu*1000.);
1115 sprintf(meddesc,"SH");
1118 for(Int_t ipart=1;ipart<=2;ipart++){
1119 for(Int_t l=1;l<=4*fLengthMax;l++){
1121 sprintf(hname,"hDisc-ContQW_%s_%d_%d_%d_%d",meddesc,fInstanceNumber,ipart,medvalue,l);
1123 Double_t wc = CalcWC(len);
1124 fHistos[ipart-1][l-1] = new TH1F(hname,hname,fgkBins,0.,fgkMaxBin*wc);
1125 fHistos[ipart-1][l-1]->SetXTitle("#Delta E [GeV]");
1126 fHistos[ipart-1][l-1]->SetYTitle("p(#Delta E)");
1127 fHistos[ipart-1][l-1]->SetLineColor(4);
1129 Double_t rrrr = CalcR(wc,len);
1130 Double_t discrete=0.;
1131 // loop on histogram channels
1132 for(Int_t bin=1; bin<=fgkBins; bin++) {
1133 Double_t xxxx = fHistos[ipart-1][l-1]->GetBinCenter(bin)/wc;
1134 Double_t continuous;
1136 CalcMult(ipart,rrrr,xxxx,continuous,discrete);
1138 CalcSingleHard(ipart,rrrr,xxxx,continuous,discrete);
1139 fHistos[ipart-1][l-1]->SetBinContent(bin,continuous);
1141 // add discrete part to distribution
1143 for(Int_t bin=2;bin<=fgkBins;bin++)
1144 fHistos[ipart-1][l-1]->SetBinContent(bin,0.);
1146 Double_t val=discrete/(1.-discrete)*fHistos[ipart-1][l-1]->Integral(1,fgkBins);
1147 fHistos[ipart-1][l-1]->Fill(0.,val);
1149 Double_t hint=fHistos[ipart-1][l-1]->Integral(1,fgkBins);
1150 fHistos[ipart-1][l-1]->Scale(1./hint);
1156 Int_t AliQuenchingWeights::SampleEnergyLoss(Int_t ipart, Double_t R)
1158 // Sample energy loss directly for one particle type
1159 // choses R (safe it and keep it until next call of function)
1161 // read-in data before first call
1163 Error("SampleEnergyLoss","Tables are not loaded.");
1167 Double_t discrete=0.;
1168 Double_t continuous=0;;
1170 Double_t xxxx = fHisto->GetBinCenter(bin);
1172 CalcMult(ipart,R,xxxx,continuous,discrete);
1174 CalcSingleHard(ipart,R,xxxx,continuous,discrete);
1177 fHisto->SetBinContent(1,1.);
1178 for(Int_t bin=2;bin<=fgkBins;bin++)
1179 fHisto->SetBinContent(bin,0.);
1183 fHisto->SetBinContent(bin,continuous);
1184 for(Int_t bin=2; bin<=fgkBins; bin++) {
1185 xxxx = fHisto->GetBinCenter(bin);
1187 CalcMult(ipart,R,xxxx,continuous,discrete);
1189 CalcSingleHard(ipart,R,xxxx,continuous,discrete);
1190 fHisto->SetBinContent(bin,continuous);
1193 Double_t val=discrete/(1.-discrete)*fHisto->Integral(1,fgkBins);
1194 fHisto->Fill(0.,val);
1195 Double_t hint=fHisto->Integral(1,fgkBins);
1197 fHisto->Scale(1./hint);
1199 //cout << discrete << " " << hint << " " << continuous << endl;
1205 const TH1F* AliQuenchingWeights::GetHisto(Int_t ipart,Double_t length) const
1207 //return quenching histograms
1208 //for ipart and length
1211 Fatal("GetELossRandom","Call SampleEnergyLoss method before!");
1214 if((ipart<1) || (ipart>2)) {
1215 Fatal("GetELossRandom","ipart =%d; but has to be 1 (quark) or 2 (gluon)",ipart);
1219 Int_t l=GetIndex(length);
1221 return fHistos[ipart-1][l-1];
1224 TH1F* AliQuenchingWeights::ComputeQWHisto(Int_t ipart,Double_t medval,Double_t length) const
1226 // ipart = 1 for quark, 2 for gluon
1227 // medval a) qtransp = transport coefficient (GeV^2/fm)
1228 // b) mu = Debye mass (GeV)
1229 // length = path length in medium (fm)
1230 // Get from SW tables:
1231 // - continuous weight, as a function of dE/wc
1234 Char_t meddesc[100];
1236 wc=CalcWC(medval,length);
1237 sprintf(meddesc,"MS");
1239 wc=CalcWCbar(medval,length);
1240 sprintf(meddesc,"SH");
1244 sprintf(hname,"hContQWHisto_%s_%d_%d_%d",meddesc,ipart,
1245 (Int_t)(medval*1000.),(Int_t)length);
1247 TH1F *hist = new TH1F("hist",hname,fgkBins,0.,fgkMaxBin*wc);
1248 hist->SetXTitle("#Delta E [GeV]");
1249 hist->SetYTitle("p(#Delta E)");
1250 hist->SetLineColor(4);
1252 Double_t rrrr = CalcR(wc,length);
1253 //loop on histogram channels
1254 for(Int_t bin=1; bin<=fgkBins; bin++) {
1255 Double_t xxxx = hist->GetBinCenter(bin)/wc;
1256 Double_t continuous,discrete;
1258 if(fMultSoft) ret=CalcMult(ipart,rrrr,xxxx,continuous,discrete);
1259 else ret=CalcSingleHard(ipart,rrrr,xxxx,continuous,discrete);
1264 hist->SetBinContent(bin,continuous);
1269 TH1F* AliQuenchingWeights::ComputeQWHistoX(Int_t ipart,Double_t medval,Double_t length) const
1271 // ipart = 1 for quark, 2 for gluon
1272 // medval a) qtransp = transport coefficient (GeV^2/fm)
1273 // b) mu = Debye mass (GeV)
1274 // length = path length in medium (fm)
1275 // Get from SW tables:
1276 // - continuous weight, as a function of dE/wc
1279 Char_t meddesc[100];
1281 wc=CalcWC(medval,length);
1282 sprintf(meddesc,"MS");
1284 wc=CalcWCbar(medval,length);
1285 sprintf(meddesc,"SH");
1289 sprintf(hname,"hContQWHistox_%s_%d_%d_%d",meddesc,ipart,
1290 (Int_t)(medval*1000.),(Int_t)length);
1292 TH1F *histx = new TH1F("histx",hname,fgkBins,0.,fgkMaxBin);
1293 histx->SetXTitle("x = #Delta E/#omega_{c}");
1295 histx->SetYTitle("p(#Delta E/#omega_{c})");
1297 histx->SetYTitle("p(#Delta E/#bar#omega_{c})");
1298 histx->SetLineColor(4);
1300 Double_t rrrr = CalcR(wc,length);
1301 //loop on histogram channels
1302 for(Int_t bin=1; bin<=fgkBins; bin++) {
1303 Double_t xxxx = histx->GetBinCenter(bin);
1304 Double_t continuous,discrete;
1306 if(fMultSoft) ret=CalcMult(ipart,rrrr,xxxx,continuous,discrete);
1307 else ret=CalcSingleHard(ipart,rrrr,xxxx,continuous,discrete);
1312 histx->SetBinContent(bin,continuous);
1317 TH1F* AliQuenchingWeights::ComputeQWHistoX(Int_t ipart,Double_t R) const
1319 // compute P(E) distribution for
1320 // given ipart = 1 for quark, 2 for gluon
1323 Char_t meddesc[100];
1325 sprintf(meddesc,"MS");
1327 sprintf(meddesc,"SH");
1331 sprintf(hname,"hQWHistox_%s_%d_%.2f",meddesc,ipart,R);
1332 TH1F *histx = new TH1F("histx",hname,fgkBins,0.,fgkMaxBin);
1333 histx->SetXTitle("x = #Delta E/#omega_{c}");
1335 histx->SetYTitle("p(#Delta E/#omega_{c})");
1337 histx->SetYTitle("p(#Delta E/#bar#omega_{c})");
1338 histx->SetLineColor(4);
1341 Double_t continuous=0.,discrete=0.;
1342 //loop on histogram channels
1343 for(Int_t bin=1; bin<=fgkBins; bin++) {
1344 Double_t xxxx = histx->GetBinCenter(bin);
1346 if(fMultSoft) ret=CalcMult(ipart,rrrr,xxxx,continuous,discrete);
1347 else ret=CalcSingleHard(ipart,rrrr,xxxx,continuous,discrete);
1352 histx->SetBinContent(bin,continuous);
1355 //add discrete part to distribution
1357 for(Int_t bin=2;bin<=fgkBins;bin++)
1358 histx->SetBinContent(bin,0.);
1360 Double_t val=discrete/(1.-discrete)*histx->Integral(1,fgkBins);
1361 histx->Fill(0.,val);
1363 Double_t hint=histx->Integral(1,fgkBins);
1364 if(hint!=0) histx->Scale(1./hint);
1369 TH1F* AliQuenchingWeights::ComputeELossHisto(Int_t ipart,Double_t medval,Double_t l,Double_t e) const
1371 // compute energy loss histogram for
1372 // parton type, medium value, length and energy
1374 AliQuenchingWeights *dummy=new AliQuenchingWeights(*this);
1376 dummy->SetQTransport(medval);
1379 dummy->SetMu(medval);
1380 dummy->InitSingleHard();
1382 dummy->SampleEnergyLoss();
1387 sprintf(name,"Energy Loss Distribution - Quarks;E_{loss} (GeV);#");
1388 sprintf(hname,"hLossQuarks");
1390 sprintf(name,"Energy Loss Distribution - Gluons;E_{loss} (GeV);#");
1391 sprintf(hname,"hLossGluons");
1394 TH1F *h = new TH1F(hname,name,250,0,250);
1395 for(Int_t i=0;i<100000;i++){
1396 //if(i % 1000 == 0) cout << "." << flush;
1397 Double_t loss=dummy->GetELossRandom(ipart,l,e);
1405 TH1F* AliQuenchingWeights::ComputeELossHisto(Int_t ipart,Double_t medval,TH1F *hEll,Double_t e) const
1407 // compute energy loss histogram for
1408 // parton type, medium value,
1409 // length distribution and energy
1411 AliQuenchingWeights *dummy=new AliQuenchingWeights(*this);
1413 dummy->SetQTransport(medval);
1416 dummy->SetMu(medval);
1417 dummy->InitSingleHard();
1419 dummy->SampleEnergyLoss();
1424 sprintf(name,"Energy Loss Distribution - Quarks;E_{loss} (GeV);#");
1425 sprintf(hname,"hLossQuarks");
1427 sprintf(name,"Energy Loss Distribution - Gluons;E_{loss} (GeV);#");
1428 sprintf(hname,"hLossGluons");
1431 TH1F *h = new TH1F(hname,name,250,0,250);
1432 for(Int_t i=0;i<100000;i++){
1433 //if(i % 1000 == 0) cout << "." << flush;
1434 Double_t loss=dummy->GetELossRandom(ipart,hEll,e);
1442 TH1F* AliQuenchingWeights::ComputeELossHisto(Int_t ipart,Double_t R) const
1444 // compute energy loss histogram for
1445 // parton type and given R
1447 TH1F *dummy = ComputeQWHistoX(ipart,R);
1448 if(!dummy) return 0;
1451 sprintf(hname,"hELossHistox_%d_%.2f",ipart,R);
1452 TH1F *histx = new TH1F("histxr",hname,fgkBins,0.,fgkMaxBin);
1453 for(Int_t i=0;i<100000;i++){
1454 //if(i % 1000 == 0) cout << "." << flush;
1455 Double_t loss=dummy->GetRandom();
1462 Double_t AliQuenchingWeights::GetMeanELoss(Int_t ipart,Double_t medval,Double_t l) const
1464 // compute average energy loss for
1465 // parton type, medium value, length and energy
1467 TH1F *dummy = ComputeELossHisto(ipart,medval,l);
1468 if(!dummy) return 0;
1469 Double_t ret=dummy->GetMean();
1474 Double_t AliQuenchingWeights::GetMeanELoss(Int_t ipart,Double_t medval,TH1F *hEll) const
1476 // compute average energy loss for
1477 // parton type, medium value,
1478 // length distribution and energy
1480 TH1F *dummy = ComputeELossHisto(ipart,medval,hEll);
1481 if(!dummy) return 0;
1482 Double_t ret=dummy->GetMean();
1487 Double_t AliQuenchingWeights::GetMeanELoss(Int_t ipart,Double_t R) const
1489 // compute average energy loss over wc
1490 // for parton type and given R
1492 TH1F *dummy = ComputeELossHisto(ipart,R);
1493 if(!dummy) return 0;
1494 Double_t ret=dummy->GetMean();
1499 void AliQuenchingWeights::PlotDiscreteWeights(Double_t len) const
1501 // plot discrete weights for given length
1505 c = new TCanvas("cdiscms","Discrete Weight for Multiple Scattering",0,0,500,400);
1507 c = new TCanvas("cdiscsh","Discrete Weight for Single Hard Scattering",0,0,500,400);
1510 TH2F *hframe = new TH2F("hdisc","",2,0,5.1,2,0,1);
1511 hframe->SetStats(0);
1513 hframe->SetXTitle("#hat{q} [GeV^{2}/fm]");
1515 hframe->SetXTitle("#mu [GeV]");
1516 hframe->SetYTitle("Probability #Delta E = 0 , p_{0}");
1519 TGraph *gq=new TGraph(20);
1522 for(Double_t q=0.05;q<=5.05;q+=0.25){
1524 CalcMult(1,1.0,q,len,cont,disc);
1525 gq->SetPoint(i,q,disc);i++;
1528 for(Double_t m=0.05;m<=5.05;m+=0.25){
1530 CalcSingleHard(1,1.0,m,len,cont, disc);
1531 gq->SetPoint(i,m,disc);i++;
1534 gq->SetMarkerStyle(20);
1537 TGraph *gg=new TGraph(20);
1540 for(Double_t q=0.05;q<=5.05;q+=0.25){
1542 CalcMult(2,1.0,q,len,cont,disc);
1543 gg->SetPoint(i,q,disc);i++;
1546 for(Double_t m=0.05;m<=5.05;m+=0.25){
1548 CalcSingleHard(2,1.0,m,len,cont,disc);
1549 gg->SetPoint(i,m,disc);i++;
1552 gg->SetMarkerStyle(24);
1555 TLegend *l1a = new TLegend(0.5,0.6,.95,0.8);
1556 l1a->SetFillStyle(0);
1557 l1a->SetBorderSize(0);
1559 sprintf(label,"L = %.1f fm",len);
1560 l1a->AddEntry(gq,label,"");
1561 l1a->AddEntry(gq,"quark","pl");
1562 l1a->AddEntry(gg,"gluon","pl");
1568 void AliQuenchingWeights::PlotContWeights(Int_t itype,Double_t ell) const
1570 // plot continous weights for
1571 // given parton type and length
1578 sprintf(title,"Cont. Weight for Multiple Scattering - Quarks");
1580 sprintf(title,"Cont. Weight for Multiple Scattering - Gluons");
1582 medvals[0]=4;medvals[1]=1;medvals[2]=0.5;
1583 sprintf(name,"ccont-ms-%d",itype);
1586 sprintf(title,"Cont. Weight for Single Hard Scattering - Quarks");
1588 sprintf(title,"Cont. Weight for Single Hard Scattering - Gluons");
1590 medvals[0]=2;medvals[1]=1;medvals[2]=0.5;
1591 sprintf(name,"ccont-ms-%d",itype);
1594 TCanvas *c = new TCanvas(name,title,0,0,500,400);
1596 TH1F *h1=ComputeQWHisto(itype,medvals[0],ell);
1598 h1->SetTitle(title);
1600 h1->SetLineColor(1);
1602 TH1F *h2=ComputeQWHisto(itype,medvals[1],ell);
1604 h2->SetLineColor(2);
1605 h2->DrawCopy("SAME");
1606 TH1F *h3=ComputeQWHisto(itype,medvals[2],ell);
1608 h3->SetLineColor(3);
1609 h3->DrawCopy("SAME");
1611 TLegend *l1a = new TLegend(0.5,0.6,.95,0.8);
1612 l1a->SetFillStyle(0);
1613 l1a->SetBorderSize(0);
1615 sprintf(label,"L = %.1f fm",ell);
1616 l1a->AddEntry(h1,label,"");
1618 sprintf(label,"#hat{q} = %.1f GeV^{2}/fm",medvals[0]);
1619 l1a->AddEntry(h1,label,"pl");
1620 sprintf(label,"#hat{q} = %.1f GeV^{2}/fm",medvals[1]);
1621 l1a->AddEntry(h2,label,"pl");
1622 sprintf(label,"#hat{q} = %.1f GeV^{2}/fm",medvals[2]);
1623 l1a->AddEntry(h3,label,"pl");
1625 sprintf(label,"#mu = %.1f GeV",medvals[0]);
1626 l1a->AddEntry(h1,label,"pl");
1627 sprintf(label,"#mu = %.1f GeV",medvals[1]);
1628 l1a->AddEntry(h2,label,"pl");
1629 sprintf(label,"#mu = %.1f GeV",medvals[2]);
1630 l1a->AddEntry(h3,label,"pl");
1637 void AliQuenchingWeights::PlotContWeightsVsL(Int_t itype,Double_t medval) const
1639 // plot continous weights for
1640 // given parton type and medium value
1646 sprintf(title,"Cont. Weight for Multiple Scattering - Quarks");
1648 sprintf(title,"Cont. Weight for Multiple Scattering - Gluons");
1650 sprintf(name,"ccont2-ms-%d",itype);
1653 sprintf(title,"Cont. Weight for Single Hard Scattering - Quarks");
1655 sprintf(title,"Cont. Weight for Single Hard Scattering - Gluons");
1657 sprintf(name,"ccont2-sh-%d",itype);
1659 TCanvas *c = new TCanvas(name,title,0,0,500,400);
1661 TH1F *h1=ComputeQWHisto(itype,medval,8);
1663 h1->SetTitle(title);
1665 h1->SetLineColor(1);
1667 TH1F *h2=ComputeQWHisto(itype,medval,5);
1669 h2->SetLineColor(2);
1670 h2->DrawCopy("SAME");
1671 TH1F *h3=ComputeQWHisto(itype,medval,2);
1673 h3->SetLineColor(3);
1674 h3->DrawCopy("SAME");
1676 TLegend *l1a = new TLegend(0.5,0.6,.95,0.8);
1677 l1a->SetFillStyle(0);
1678 l1a->SetBorderSize(0);
1681 sprintf(label,"#hat{q} = %.1f GeV^{2}/fm",medval);
1683 sprintf(label,"#mu = %.1f GeV",medval);
1685 l1a->AddEntry(h1,label,"");
1686 l1a->AddEntry(h1,"L = 8 fm","pl");
1687 l1a->AddEntry(h2,"L = 5 fm","pl");
1688 l1a->AddEntry(h3,"L = 2 fm","pl");
1694 void AliQuenchingWeights::PlotAvgELoss(Double_t len,Double_t e) const
1696 // plot average energy loss for given length
1697 // and parton energy
1700 Error("PlotAvgELoss","Tables are not loaded.");
1707 sprintf(title,"Average Energy Loss for Multiple Scattering");
1708 sprintf(name,"cavgelossms");
1710 sprintf(title,"Average Energy Loss for Single Hard Scattering");
1711 sprintf(name,"cavgelosssh");
1714 TCanvas *c = new TCanvas(name,title,0,0,500,400);
1716 TH2F *hframe = new TH2F("avgloss",title,2,0,5.1,2,0,100);
1717 hframe->SetStats(0);
1719 hframe->SetXTitle("#hat{q} [GeV^{2}/fm]");
1721 hframe->SetXTitle("#mu [GeV]");
1722 hframe->SetYTitle("<E_{loss}> [GeV]");
1725 TGraph *gq=new TGraph(20);
1727 for(Double_t v=0.05;v<=5.05;v+=0.25){
1728 TH1F *dummy=ComputeELossHisto(1,v,len,e);
1729 Double_t avgloss=dummy->GetMean();
1730 gq->SetPoint(i,v,avgloss);i++;
1733 gq->SetMarkerStyle(20);
1736 TGraph *gg=new TGraph(20);
1738 for(Double_t v=0.05;v<=5.05;v+=0.25){
1739 TH1F *dummy=ComputeELossHisto(2,v,len,e);
1740 Double_t avgloss=dummy->GetMean();
1741 gg->SetPoint(i,v,avgloss);i++;
1744 gg->SetMarkerStyle(24);
1747 TGraph *gratio=new TGraph(20);
1748 for(Int_t i=0;i<20;i++){
1750 gg->GetPoint(i,x,y);
1751 gq->GetPoint(i,x2,y2);
1753 gratio->SetPoint(i,x,y/y2*10/2.25);
1754 else gratio->SetPoint(i,x,0);
1756 gratio->SetLineStyle(4);
1758 TLegend *l1a = new TLegend(0.5,0.6,.95,0.8);
1759 l1a->SetFillStyle(0);
1760 l1a->SetBorderSize(0);
1762 sprintf(label,"L = %.1f fm",len);
1763 l1a->AddEntry(gq,label,"");
1764 l1a->AddEntry(gq,"quark","pl");
1765 l1a->AddEntry(gg,"gluon","pl");
1766 l1a->AddEntry(gratio,"gluon/quark/2.25*10","pl");
1772 void AliQuenchingWeights::PlotAvgELoss(TH1F *hEll,Double_t e) const
1774 // plot average energy loss for given
1775 // length distribution and parton energy
1778 Error("PlotAvgELossVs","Tables are not loaded.");
1785 sprintf(title,"Average Energy Loss for Multiple Scattering");
1786 sprintf(name,"cavgelossms2");
1788 sprintf(title,"Average Energy Loss for Single Hard Scattering");
1789 sprintf(name,"cavgelosssh2");
1792 TCanvas *c = new TCanvas(name,title,0,0,500,400);
1794 TH2F *hframe = new TH2F("havgloss",title,2,0,5.1,2,0,100);
1795 hframe->SetStats(0);
1797 hframe->SetXTitle("#hat{q} [GeV^{2}/fm]");
1799 hframe->SetXTitle("#mu [GeV]");
1800 hframe->SetYTitle("<E_{loss}> [GeV]");
1803 TGraph *gq=new TGraph(20);
1805 for(Double_t v=0.05;v<=5.05;v+=0.25){
1806 TH1F *dummy=ComputeELossHisto(1,v,hEll,e);
1807 Double_t avgloss=dummy->GetMean();
1808 gq->SetPoint(i,v,avgloss);i++;
1811 gq->SetMarkerStyle(20);
1814 TGraph *gg=new TGraph(20);
1816 for(Double_t v=0.05;v<=5.05;v+=0.25){
1817 TH1F *dummy=ComputeELossHisto(2,v,hEll,e);
1818 Double_t avgloss=dummy->GetMean();
1819 gg->SetPoint(i,v,avgloss);i++;
1822 gg->SetMarkerStyle(24);
1825 TGraph *gratio=new TGraph(20);
1826 for(Int_t i=0;i<20;i++){
1828 gg->GetPoint(i,x,y);
1829 gq->GetPoint(i,x2,y2);
1831 gratio->SetPoint(i,x,y/y2*10/2.25);
1832 else gratio->SetPoint(i,x,0);
1834 gratio->SetLineStyle(4);
1837 TLegend *l1a = new TLegend(0.5,0.6,.95,0.8);
1838 l1a->SetFillStyle(0);
1839 l1a->SetBorderSize(0);
1841 sprintf(label,"<L> = %.2f fm",hEll->GetMean());
1842 l1a->AddEntry(gq,label,"");
1843 l1a->AddEntry(gq,"quark","pl");
1844 l1a->AddEntry(gg,"gluon","pl");
1845 l1a->AddEntry(gratio,"gluon/quark/2.25*10","pl");
1851 void AliQuenchingWeights::PlotAvgELossVsL(Double_t e) const
1853 // plot average energy loss versus ell
1856 Error("PlotAvgELossVsEll","Tables are not loaded.");
1864 sprintf(title,"Average Energy Loss for Multiple Scattering");
1865 sprintf(name,"cavgelosslms");
1868 sprintf(title,"Average Energy Loss for Single Hard Scattering");
1869 sprintf(name,"cavgelosslsh");
1873 TCanvas *c = new TCanvas(name,title,0,0,600,400);
1875 TH2F *hframe = new TH2F("avglossell",title,2,0,fLengthMax,2,0,250);
1876 hframe->SetStats(0);
1877 hframe->SetXTitle("length [fm]");
1878 hframe->SetYTitle("<E_{loss}> [GeV]");
1881 TGraph *gq=new TGraph((Int_t)fLengthMax*4);
1883 for(Double_t v=0.25;v<=fLengthMax;v+=0.25){
1884 TH1F *dummy=ComputeELossHisto(1,medval,v,e);
1885 Double_t avgloss=dummy->GetMean();
1886 gq->SetPoint(i,v,avgloss);i++;
1889 gq->SetMarkerStyle(20);
1892 TGraph *gg=new TGraph((Int_t)fLengthMax*4);
1894 for(Double_t v=0.25;v<=fLengthMax;v+=0.25){
1895 TH1F *dummy=ComputeELossHisto(2,medval,v,e);
1896 Double_t avgloss=dummy->GetMean();
1897 gg->SetPoint(i,v,avgloss);i++;
1900 gg->SetMarkerStyle(24);
1903 TGraph *gratio=new TGraph((Int_t)fLengthMax*4);
1904 for(Int_t i=0;i<=(Int_t)fLengthMax*4;i++){
1906 gg->GetPoint(i,x,y);
1907 gq->GetPoint(i,x2,y2);
1909 gratio->SetPoint(i,x,y/y2*100/2.25);
1910 else gratio->SetPoint(i,x,0);
1912 gratio->SetLineStyle(1);
1913 gratio->SetLineWidth(2);
1915 TLegend *l1a = new TLegend(0.15,0.65,.60,0.85);
1916 l1a->SetFillStyle(0);
1917 l1a->SetBorderSize(0);
1920 sprintf(label,"#hat{q} = %.2f GeV^{2}/fm",medval);
1922 sprintf(label,"#mu = %.2f GeV",medval);
1923 l1a->AddEntry(gq,label,"");
1924 l1a->AddEntry(gq,"quark","pl");
1925 l1a->AddEntry(gg,"gluon","pl");
1926 l1a->AddEntry(gratio,"gluon/quark/2.25*100","pl");
1929 TF1 *f=new TF1("f","100",0,fLengthMax);
1936 void AliQuenchingWeights::PlotAvgELossVsPt(Double_t medval,Double_t len) const
1938 // plot relative energy loss for given
1939 // length and parton energy versus pt
1942 Error("PlotAvgELossVsPt","Tables are not loaded.");
1949 sprintf(title,"Relative Energy Loss for Multiple Scattering");
1950 sprintf(name,"cavgelossvsptms");
1952 sprintf(title,"Relative Energy Loss for Single Hard Scattering");
1953 sprintf(name,"cavgelossvsptsh");
1956 TCanvas *c = new TCanvas(name,title,0,0,500,400);
1958 TH2F *hframe = new TH2F("havglossvspt",title,2,0,100,2,0,1);
1959 hframe->SetStats(0);
1960 hframe->SetXTitle("p_{T} [GeV]");
1961 hframe->SetYTitle("<E_{loss}>/p_{T} [GeV]");
1964 TGraph *gq=new TGraph(40);
1966 for(Double_t pt=2.5;pt<=100.05;pt+=2.5){
1967 TH1F *dummy=ComputeELossHisto(1,medval,len,pt);
1968 Double_t avgloss=dummy->GetMean();
1969 gq->SetPoint(i,pt,avgloss/pt);i++;
1972 gq->SetMarkerStyle(20);
1975 TGraph *gg=new TGraph(40);
1977 for(Double_t pt=2.5;pt<=100.05;pt+=2.5){
1978 TH1F *dummy=ComputeELossHisto(2,medval,len,pt);
1979 Double_t avgloss=dummy->GetMean();
1980 gg->SetPoint(i,pt,avgloss/pt);i++;
1983 gg->SetMarkerStyle(24);
1986 TLegend *l1a = new TLegend(0.5,0.6,.95,0.8);
1987 l1a->SetFillStyle(0);
1988 l1a->SetBorderSize(0);
1990 sprintf(label,"L = %.1f fm",len);
1991 l1a->AddEntry(gq,label,"");
1992 l1a->AddEntry(gq,"quark","pl");
1993 l1a->AddEntry(gg,"gluon","pl");
1999 void AliQuenchingWeights::PlotAvgELossVsPt(Double_t medval,TH1F *hEll) const
2001 // plot relative energy loss for given
2002 // length distribution and parton energy versus pt
2005 Error("PlotAvgELossVsPt","Tables are not loaded.");
2012 sprintf(title,"Relative Energy Loss for Multiple Scattering");
2013 sprintf(name,"cavgelossvsptms2");
2015 sprintf(title,"Relative Energy Loss for Single Hard Scattering");
2016 sprintf(name,"cavgelossvsptsh2");
2018 TCanvas *c = new TCanvas(name,title,0,0,500,400);
2020 TH2F *hframe = new TH2F("havglossvspt",title,2,0,100,2,0,1);
2021 hframe->SetStats(0);
2022 hframe->SetXTitle("p_{T} [GeV]");
2023 hframe->SetYTitle("<E_{loss}>/p_{T} [GeV]");
2026 TGraph *gq=new TGraph(40);
2028 for(Double_t pt=2.5;pt<=100.05;pt+=2.5){
2029 TH1F *dummy=ComputeELossHisto(1,medval,hEll,pt);
2030 Double_t avgloss=dummy->GetMean();
2031 gq->SetPoint(i,pt,avgloss/pt);i++;
2034 gq->SetMarkerStyle(20);
2037 TGraph *gg=new TGraph(40);
2039 for(Double_t pt=2.5;pt<=100.05;pt+=2.5){
2040 TH1F *dummy=ComputeELossHisto(2,medval,hEll,pt);
2041 Double_t avgloss=dummy->GetMean();
2042 gg->SetPoint(i,pt,avgloss/pt);i++;
2045 gg->SetMarkerStyle(24);
2048 TLegend *l1a = new TLegend(0.5,0.6,.95,0.8);
2049 l1a->SetFillStyle(0);
2050 l1a->SetBorderSize(0);
2052 sprintf(label,"<L> = %.2f fm",hEll->GetMean());
2053 l1a->AddEntry(gq,label,"");
2054 l1a->AddEntry(gq,"quark","pl");
2055 l1a->AddEntry(gg,"gluon","pl");
2061 Int_t AliQuenchingWeights::GetIndex(Double_t len) const
2063 if(len>fLengthMax) len=fLengthMax;
2065 Int_t l=Int_t(len/0.25);
2066 if((len-l*0.25)>0.125) l++;