]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibLaser.cxx
Adding fits of the distortions (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibLaser.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /*
17   //
18   // FUNCTIONALITY:
19   //
20   // 1. The laser track is associated with the mirror
21   //    see function FindMirror
22   //
23   // 2. The laser track is accepted for the analysis under certain condition
24   //    (see function Accpet laser)
25   // 
26   // 3. The drift velocity and jitter is calculated event by event
27   //    (see function drift velocity) 
28   //
29   //
30   //
31   gSystem->Load("libANALYSIS");
32   gSystem->Load("libTPCcalib");
33   TFile fcalib("CalibObjects.root");
34   TObjArray * array = (TObjArray*)fcalib.Get("TPCCalib");
35   AliTPCcalibLaser * laser = ( AliTPCcalibLaser *)array->FindObject("laserTPC");
36   laser->DumpMeanInfo(-0.4)
37   TFile fmean("laserMean.root")
38
39   //  laser track clasification;
40
41   TCut cutT("cutT","abs(Tr.fP[3])<0.06");
42   TCut cutPt("cutPt","abs(Tr.fP[4])<0.1");
43   TCut cutN("cutN","fTPCncls>70");
44   TCut cutP("cutP","abs(atan2(x1,x0)-atan2(lx1,lx0))<0.03")
45   TCut cutA = cutT+cutPt+cutP;
46
47   TFile f("laserTPCDebug.root");
48   TTree * treeT = (TTree*)f.Get("Track");
49
50
51   // LASER scan 
52   gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
53   gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
54   AliXRDPROOFtoolkit tool; 
55   TChain * chain = tool.MakeChain("laserScan.txt","Mean",0,10200);
56   chain->Lookup();
57
58
59   treeT->Draw("(atan2(x1,x0)-atan2(lx1,lx0))*250.:fBundle","fSide==1&&fRod==0"+cutA,"prof") 
60
61   gSystem->Load("libSTAT.so")
62   TStatToolkit toolkit;
63   Double_t chi2;
64   TVectorD fitParam;
65   TMatrixD covMatrix;
66   Int_t npoints;
67   TString *strq0 = toolkit.FitPlane(treeT,"Tr.fP[1]-LTr.fP[1]","lx1++lx2", "fSide==1"+cutA, chi2,npoints,fitParam,covMatrix);
68
69 */
70
71
72
73 #include "TLinearFitter.h"
74 #include "AliTPCcalibLaser.h"
75 #include "AliExternalTrackParam.h"
76 #include "AliESDEvent.h"
77 #include "AliESDfriend.h"
78 #include "AliESDtrack.h"
79 #include "AliTPCTracklet.h"
80 #include "TH1D.h"
81 #include "TVectorD.h"
82 #include "TTreeStream.h"
83 #include "TFile.h"
84 #include "TF1.h"
85 #include "TGraphErrors.h"
86 #include "AliTPCclusterMI.h"
87 #include "AliTPCseed.h"
88 #include "AliTracker.h"
89 #include "TClonesArray.h"
90
91
92 #include "TTreeStream.h"
93 #include <iostream>
94 #include <sstream>
95 #include "AliTPCLaserTrack.h"
96
97 using namespace std;
98
99 ClassImp(AliTPCcalibLaser)
100
101 AliTPCcalibLaser::AliTPCcalibLaser():
102   AliTPCcalibBase(),
103   fESD(0),
104   fESDfriend(),
105   fTracksMirror(336),
106   fTracksEsd(336),
107   fTracksEsdParam(336),
108   fTracksTPC(336),
109   fDeltaZ(336),          // array of histograms of delta z for each track
110   fDeltaPhi(336),          // array of histograms of delta z for each track
111   fDeltaPhiP(336),          // array of histograms of delta z for each track
112   fSignals(336),           // array of dedx signals
113   fFitAside(new TVectorD(3)),        // drift fit - A side
114   fFitCside(new TVectorD(3)),        // drift fit - C- side
115   fRun(0)
116 {
117   //
118   // Constructor
119   //
120   fTracksEsdParam.SetOwner(kTRUE);
121 }
122
123 AliTPCcalibLaser::AliTPCcalibLaser(const Text_t *name, const Text_t *title):
124   AliTPCcalibBase(),
125   fESD(0),
126   fESDfriend(0),
127   fTracksMirror(336),
128   fTracksEsd(336),
129   fTracksEsdParam(336),
130   fTracksTPC(336),
131   fDeltaZ(336),          // array of histograms of delta z for each track
132   fDeltaPhi(336),          // array of histograms of delta z for each track
133   fDeltaPhiP(336),          // array of histograms of delta z for each track
134   fSignals(336),           // array of dedx signals
135   fFitAside(new TVectorD(3)),        // drift fit - A side
136   fFitCside(new TVectorD(3)),        // drift fit - C- side
137   fRun(0)
138 {
139   SetName(name);
140   SetTitle(title);
141   //
142   // Constructor
143   //
144   fTracksEsdParam.SetOwner(kTRUE);  
145 }
146
147 AliTPCcalibLaser::~AliTPCcalibLaser() {
148   //
149   // destructor
150   //
151 }
152
153
154
155 void AliTPCcalibLaser::Process(AliESDEvent * event) {
156   //
157   // 
158   // Loop over tracks and call  Process function
159   //
160   fESD = event;
161   if (!fESD) {
162     return;
163   }
164   fESDfriend=static_cast<AliESDfriend*>(fESD->FindListObject("AliESDfriend"));
165   if (!fESDfriend) {
166     return;
167   }
168   fTracksTPC.Clear();
169   fTracksEsd.Clear();
170   fTracksEsdParam.Delete();
171   //
172   Int_t n=fESD->GetNumberOfTracks();
173   Int_t run = fESD->GetRunNumber();
174   fRun = run;
175   for (Int_t i=0;i<n;++i) {
176     AliESDfriendTrack *friendTrack=fESDfriend->GetTrack(i);
177     AliESDtrack *track=fESD->GetTrack(i);
178     TObject *calibObject=0;
179     AliTPCseed *seed=0;
180     for (Int_t j=0;(calibObject=friendTrack->GetCalibObject(j));++j)
181       if ((seed=dynamic_cast<AliTPCseed*>(calibObject)))
182         break;
183     if (track&&seed) FindMirror(track,seed);
184     //
185   }
186   
187   FitDriftV();
188   MakeDistHisto();
189   //
190   for (Int_t id=0; id<336; id++){
191     //
192     //
193     if (!fTracksEsdParam.At(id)) continue;
194     DumpLaser(id);
195     RefitLaser(id);    
196     
197   }
198 }
199
200 void AliTPCcalibLaser::MakeDistHisto(){
201   //
202   //
203   //
204   for (Int_t id=0; id<336; id++){
205     //
206     //
207     if (!fTracksEsdParam.At(id)) continue;  
208     if (!AcceptLaser(id)) continue;
209     //
210     //
211     TH1F * hisdz = (TH1F*)fDeltaZ.At(id);
212     TH1F * hisdphi = (TH1F*)fDeltaPhi.At(id);
213     TH1F * hisdphiP = (TH1F*)fDeltaPhiP.At(id);
214     TH1F * hisSignal = (TH1F*)fSignals.At(id);
215
216     if (!hisdz){      
217       hisdz = new TH1F(Form("hisdz%d",id),Form("hisdz%d",id),500,-10,10);
218       fDeltaZ.AddAt(hisdz,id);
219       //
220       hisdphi = new TH1F(Form("hisdphi%d",id),Form("hisdphi%d",id),500,-1,1);
221       fDeltaPhi.AddAt(hisdphi,id);
222       //
223       hisdphiP = new TH1F(Form("hisdphiP%d",id),Form("hisdphiP%d",id),500,-0.01,0.01);
224       fDeltaPhiP.AddAt(hisdphiP,id);
225       hisSignal = new TH1F(Form("hisSignal%d",id),Form("hisSignal%d",id),500,0,1000);
226       fSignals.AddAt(hisSignal,id);
227     }
228
229     AliExternalTrackParam *param=(AliExternalTrackParam*)fTracksEsdParam.At(id);
230     AliTPCLaserTrack *ltrp = ( AliTPCLaserTrack*)fTracksMirror.At(id);
231     AliESDtrack   *track    = (AliESDtrack*)fTracksEsd.At(id);
232     Double_t xyz[3];
233     Double_t pxyz[3];
234     Double_t lxyz[3];
235     Double_t lpxyz[3];
236     param->GetXYZ(xyz);
237     param->GetPxPyPz(pxyz);
238     ltrp->GetXYZ(lxyz);
239     ltrp->GetPxPyPz(lpxyz);
240     //
241     Float_t dz   = param->GetZ()-ltrp->GetZ();
242     Float_t dphi = (TMath::ATan2(xyz[1],xyz[0])- TMath::ATan2(lxyz[1],lxyz[0]))*254.;
243     Float_t dphiP = param->GetParameter()[2]-ltrp->GetParameter()[2];
244     hisdz->Fill(dz);
245     hisdphi->Fill(dphi);
246     hisdphiP->Fill(dphiP); 
247     hisSignal->Fill(track->GetTPCsignal());
248   }
249 }
250
251 void AliTPCcalibLaser::FitDriftV(){
252   //
253   // Fit drift velocity - linear approximation in the z and global y
254   //
255   static TLinearFitter fdriftA(3,"hyp2");
256   static TLinearFitter fdriftC(3,"hyp2");
257   fdriftA.ClearPoints();
258   fdriftC.ClearPoints();
259   //
260   for (Int_t id=0; id<336; id++){
261     if (!fTracksEsdParam.At(id)) continue;  
262     if (!AcceptLaser(id)) continue;
263     AliExternalTrackParam *param=(AliExternalTrackParam*)fTracksEsdParam.At(id);
264     AliTPCLaserTrack *ltrp = ( AliTPCLaserTrack*)fTracksMirror.At(id);
265     Double_t xyz[3];
266     Double_t pxyz[3];
267     Double_t lxyz[3];
268     Double_t lpxyz[3];
269     param->GetXYZ(xyz);
270     param->GetPxPyPz(pxyz);
271     ltrp->GetXYZ(lxyz);
272     ltrp->GetPxPyPz(lpxyz);
273     Double_t xxx[2] = {lxyz[2],lxyz[1]};
274     if (ltrp->GetSide()==0){
275       fdriftA.AddPoint(xxx,xyz[2],1);
276     }else{
277       fdriftC.AddPoint(xxx,xyz[2],1);
278     }
279   }
280   Float_t chi2A = 0;
281   Float_t chi2C = 0;
282   Int_t npointsA=0;
283   Int_t npointsC=0;
284   //
285   if (fdriftA.GetNpoints()>10){
286     fdriftA.Eval();
287     fdriftA.EvalRobust(0.8);
288     fdriftA.GetParameters(*fFitAside);
289     npointsA= fdriftA.GetNpoints();
290     chi2A = fdriftA.GetChisquare()/fdriftA.GetNpoints();
291   }
292   if (fdriftC.GetNpoints()>10){
293     fdriftC.Eval();
294     fdriftC.EvalRobust(0.8);
295     fdriftC.GetParameters(*fFitCside);
296     npointsC= fdriftC.GetNpoints();
297     chi2C = fdriftC.GetChisquare()/fdriftC.GetNpoints();
298   }
299   
300   if (fStreamLevel>0){
301     TTreeSRedirector *cstream = GetDebugStreamer();
302     Int_t time = fESD->GetTimeStamp();
303     if (cstream){
304       (*cstream)<<"driftv"<<
305         "driftA.="<<fFitAside<<
306         "driftC.="<<fFitCside<<
307         "chi2A="<<chi2A<<
308         "chi2C="<<chi2C<<
309         "nA="<<npointsA<<
310         "nC="<<npointsC<<
311         "time="<<time<<
312         "\n";
313     }
314   }
315   //
316 }
317
318
319 Bool_t  AliTPCcalibLaser::AcceptLaser(Int_t id){
320   //
321   //
322   //
323   /*
324   TCut cutT("cutT","abs(Tr.fP[3])<0.06");
325   TCut cutPt("cutPt","abs(Tr.fP[4])<0.1");
326   TCut cutN("cutN","fTPCncls>70");
327   TCut cutP("cutP","abs(atan2(x1,x0)-atan2(lx1,lx0))<0.03")
328   TCut cutA = cutT+cutPt+cutP;
329   */
330   AliExternalTrackParam *param =(AliExternalTrackParam*)fTracksEsdParam.At(id);
331   AliTPCLaserTrack *ltrp  = ( AliTPCLaserTrack*)fTracksMirror.At(id);
332   AliESDtrack   *track    = (AliESDtrack*)fTracksEsd.At(id);
333
334   if (TMath::Abs(param->GetParameter()[4])>0.03) return kFALSE;
335   if (TMath::Abs(param->GetParameter()[3])>0.06) return kFALSE;  
336   if (TMath::Abs(param->GetParameter()[2]-ltrp->GetParameter()[2])>0.06) return kFALSE;
337   if (TMath::Abs(param->GetParameter()[1]-ltrp->GetParameter()[1])>10) return kFALSE;
338   //
339   // dedx cut
340   //
341   if (TMath::Abs(track->GetTPCsignal())<20) return kFALSE;
342   if (TMath::Abs(track->GetTPCsignal())>800) return kFALSE;
343   //
344   return kTRUE;  
345 }
346
347 Int_t  AliTPCcalibLaser::FindMirror(AliESDtrack *track, AliTPCseed *seed){
348   //
349   // Find corresponding mirror
350   // add the corresponding tracks
351   //
352   Float_t kRadius0 = 252;
353   Float_t kRadius  = 253.4;
354   if (!track->GetOuterParam()) return -1;
355   AliExternalTrackParam param(*(track->GetOuterParam()));
356   AliTracker::PropagateTrackTo(&param,kRadius0,0.10566,3,kTRUE);
357   AliTracker::PropagateTrackTo(&param,kRadius,0.10566,0.1,kTRUE);
358   AliTPCLaserTrack ltr;
359   AliTPCLaserTrack *ltrp=0x0;
360   //
361   Int_t id = AliTPCLaserTrack::IdentifyTrack(&param);
362   if (id!=-1 && (AliTPCLaserTrack::GetTracks()->UncheckedAt(id))) 
363     ltrp=(AliTPCLaserTrack*)AliTPCLaserTrack::GetTracks()->UncheckedAt(id);
364   else 
365     ltrp=&ltr;
366   
367   if (id>=0){
368     //
369     //
370     Float_t radius=TMath::Abs(ltrp->GetX());
371     AliTracker::PropagateTrackTo(&param,radius,0.10566,0.01,kTRUE);
372     //
373     if (!fTracksMirror.At(id)) fTracksMirror.AddAt(ltrp,id);
374     fTracksEsdParam.AddAt(param.Clone(),id);
375     fTracksEsd.AddAt(track,id);
376     fTracksTPC.AddAt(seed,id);
377     //
378   }
379
380   return id;
381 }
382
383
384
385 void AliTPCcalibLaser::DumpLaser(Int_t id) {
386   //
387   //  Dump Laser info to the tree
388   //
389   AliESDtrack   *track    = (AliESDtrack*)fTracksEsd.At(id);
390   AliExternalTrackParam *param=(AliExternalTrackParam*)fTracksEsdParam.At(id);
391   AliTPCLaserTrack *ltrp = ( AliTPCLaserTrack*)fTracksMirror.At(id);
392   //
393   // Fast laser ID
394   //
395   Double_t xyz[3];
396   Double_t pxyz[3];
397   Double_t lxyz[3];
398   Double_t lpxyz[3];
399   param->GetXYZ(xyz);
400   param->GetPxPyPz(pxyz);
401   ltrp->GetXYZ(lxyz);
402   ltrp->GetPxPyPz(lpxyz);
403
404   if (fStreamLevel>0){
405     TTreeSRedirector *cstream = GetDebugStreamer();
406     Int_t time = fESD->GetTimeStamp();
407     Bool_t accept = AcceptLaser(id);
408     if (cstream){
409       (*cstream)<<"Track"<<
410         "run="<<fRun<<
411         "id="<<id<<
412         "accept="<<accept<<
413         "driftA.="<<fFitAside<<
414         "driftC.="<<fFitCside<<
415         "time="<<time<<
416         //
417         "LTr.="<<ltrp<<
418         "Esd.="<<track<<
419         "Tr.="<<param<<
420         "x0="<<xyz[0]<<
421         "x1="<<xyz[1]<<
422         "x2="<<xyz[2]<<
423         "px0="<<pxyz[0]<<
424         "px1="<<pxyz[1]<<
425         "px2="<<pxyz[2]<<
426         //
427         "lx0="<<lxyz[0]<<
428         "lx1="<<lxyz[1]<<
429         "lx2="<<lxyz[2]<<
430         "lpx0="<<lpxyz[0]<<
431         "lpx1="<<lpxyz[1]<<
432         "lpx2="<<lpxyz[2]<<
433         "\n";
434     }
435   }
436 }
437
438
439 void AliTPCcalibLaser::RefitLaser(Int_t id){
440   //
441   // Refit the track store residuals
442   //
443
444   AliTPCseed *track    = (AliTPCseed*)fTracksTPC.At(id);
445   AliExternalTrackParam *param=(AliExternalTrackParam*)fTracksEsdParam.At(id);
446   AliTPCLaserTrack *ltrp = (AliTPCLaserTrack*)fTracksMirror.At(id);
447                              
448   //
449   static TLinearFitter fy2(3,"hyp2");
450   static TLinearFitter fz2(3,"hyp2");
451   static TLinearFitter fy1(2,"hyp1");
452   static TLinearFitter fz1(2,"hyp1");
453   static TVectorD vecy2,vecz2,vecy1,vecz1;
454
455   const Int_t kMinClusters=20;
456   Int_t nclusters[72]; 
457   //
458   for (Int_t i=0;i<72;++i) nclusters[i]=0;
459
460   for (Int_t i=0;i<160;++i) {    
461     AliTPCclusterMI *c=track->GetClusterPointer(i);
462     if (c) nclusters[c->GetDetector()]++;
463   }
464    
465   for (Int_t isec=0; isec<72;isec++){
466     if (nclusters[isec]<kMinClusters) continue;
467     fy2.ClearPoints();
468     fz2.ClearPoints();
469     fy1.ClearPoints();
470     fz1.ClearPoints();
471     //
472     for (Int_t irow=0;irow<160;++irow) {      
473       AliTPCclusterMI *c=track->GetClusterPointer(irow);
474       //if (c && RejectCluster(c)) continue;
475       if (c&&c->GetDetector()==isec) {
476         Double_t xd = c->GetX()-120;;
477         Double_t x[2]={xd,xd*xd};
478         fy2.AddPoint(x,c->GetY());
479         fz2.AddPoint(x,c->GetZ());
480         //
481         fy1.AddPoint(x,c->GetY());
482         fz1.AddPoint(x,c->GetZ());
483       }
484     }
485     fy2.Eval();
486     fz2.Eval();
487     fy1.Eval();
488     fz1.Eval();
489     fy1.GetParameters(vecy1);
490     fy2.GetParameters(vecy2);
491     fz1.GetParameters(vecz1);
492     fz2.GetParameters(vecz2);
493     
494     if (fStreamLevel>0){
495       TTreeSRedirector *cstream = GetDebugStreamer();
496       if (cstream){
497         Float_t dedx = track->GetdEdx();
498         (*cstream)<<"Tracklet"<<
499           "LTr.="<<ltrp<<
500           "Tr.="<<param<<
501           "sec="<<isec<<
502           "ncl="<<nclusters[isec]<<
503           "dedx="<<dedx<<
504           "dedx="<<dedx<<
505           "vy1.="<<&vecy1<<
506           "vy2.="<<&vecy2<<
507           "vz1.="<<&vecz1<<
508           "vz2.="<<&vecz2<<
509           "\n";
510       }
511     }
512   }
513   //
514   //
515   //
516   //   for (Int_t irow=0;irow<160;++irow) {      
517   //       AliTPCclusterMI *c=track->GetClusterPointer(irow);
518   //       if (c && RejectCluster(c)) continue;
519   //       if (c&&c->GetDetector()==isec) {
520   //    Double_t x[2]={c->GetX(),c->GetX()*c->GetX()};
521   //    fy2.AddPoint(&x,c->GetY());
522   //    fz2.AddPoint(&x,c->GetZ());
523   //    //
524   //    fy1.AddPoint(&x,c->GetY());
525   //    fz1.AddPoint(&x,c->GetZ());
526   //       }
527   //     }    
528   
529 }
530
531
532 void AliTPCcalibLaser::DumpMeanInfo(Float_t bfield){
533   //
534   //
535   //
536   AliTPCcalibLaser *laser = this;
537   TTreeSRedirector *pcstream = new TTreeSRedirector("laserMean.root");
538   TF1 fg("fg","gaus");
539   //
540   //
541   for (Int_t id=0; id<336; id++){
542     TH1F * hisphi  = (TH1F*)laser->fDeltaPhi.At(id);
543     TH1F * hisphiP = (TH1F*)laser->fDeltaPhiP.At(id);
544     TH1F * hisZ    = (TH1F*)laser->fDeltaZ.At(id);
545     //    TH1F * hisS    = laser->fDeltaZ->At(id);
546     if (!hisphi) continue;;
547     Double_t entries = hisphi->GetEntries();
548     if (entries<30) continue;
549     //
550     AliTPCLaserTrack *ltrp = (AliTPCLaserTrack*)fTracksMirror.At(id);
551     if (!ltrp) {
552      AliTPCLaserTrack::LoadTracks();
553       ltrp =(AliTPCLaserTrack*)AliTPCLaserTrack::GetTracks()->UncheckedAt(id);
554     }
555     Float_t meanphi = hisphi->GetMean();
556     Float_t rmsphi = hisphi->GetRMS();
557     Float_t meanphiP = hisphiP->GetMean();
558     Float_t rmsphiP = hisphiP->GetRMS();
559     Float_t meanZ = hisZ->GetMean();
560     Float_t rmsZ = hisZ->GetRMS();
561     hisphi->Fit(&fg,"","",hisphi->GetMean()-4*hisphi->GetRMS(),hisphi->GetMean()+4*hisphi->GetRMS());
562     Double_t gphi1 = fg.GetParameter(1); 
563     Double_t gphi2 = fg.GetParameter(2); 
564     hisphiP->Fit(&fg,"","",hisphiP->GetMean()-4*hisphiP->GetRMS(),hisphiP->GetMean()+4*hisphiP->GetRMS());
565     Double_t gphiP1 = fg.GetParameter(1); 
566     Double_t gphiP2 = fg.GetParameter(2); 
567     hisZ->Fit(&fg,"","",hisZ->GetMean()-4*hisZ->GetRMS(),hisZ->GetMean()+4*hisZ->GetRMS());
568     Double_t gz1 = fg.GetParameter(1); 
569     Double_t gz2 = fg.GetParameter(2); 
570
571     //
572     Double_t lxyz[3];
573     Double_t lpxyz[3];
574     ltrp->GetXYZ(lxyz);
575     ltrp->GetPxPyPz(lpxyz);
576     //
577     (*pcstream)<<"Mean"<<
578       "entries="<<entries<<      // number of entries
579       "bz="<<bfield<<            // bfield
580       "LTr.="<<ltrp<<             // refernece track
581       "lx0="<<lxyz[0]<<          // reference x
582       "lx1="<<lxyz[1]<<          // reference y
583       "lx2="<<lxyz[2]<<          // refernece z      
584       "lpx0="<<lpxyz[0]<<          // reference x
585       "lpx1="<<lpxyz[1]<<          // reference y
586       "lpx2="<<lpxyz[2]<<          // refernece z            
587       "mphi="<<meanphi<<         //
588       "rmsphi="<<rmsphi<<        //
589       "gphi1="<<gphi1<<
590       "gphi2="<<gphi2<<
591       "mphiP="<<meanphiP<<       //
592       "rmsphiP="<<rmsphiP<<      //
593       "gphiP1="<<gphiP1<<
594       "gphiP2="<<gphiP2<<
595       "meanZ="<<meanZ<<
596       "rmsZ="<<rmsZ<<
597       "gz1="<<gz1<<
598       "gz2="<<gz2<<
599
600       "\n";
601   }
602   delete pcstream;
603 }
604
605 void AliTPCcalibLaser::Analyze(){
606   //
607   //
608   //
609 }
610
611
612
613