scripts/CheckAlign.C

Go to the documentation of this file.
00001 //____________________________________________________________________
00002 //
00003 // $Id: CheckAlign.C,v 1.1 2006/03/23 17:10:48 cholm Exp $
00004 //
00005 // Script that contains a class to compare the raw data written to the
00006 // digits it's created from.
00007 //
00008 // Use the script `Compile.C' to compile this class using ACLic. 
00009 //
00010 #include <AliLog.h>
00011 #include <AliFMDHit.h>
00012 #include <AliFMDDigit.h>
00013 #include <AliFMDInput.h>
00014 #include <AliFMDUShortMap.h>
00015 #include <AliFMDParameters.h>
00016 #include <AliFMDGeometry.h>
00017 #include <AliFMDRing.h>
00018 #include <AliFMDDetector.h>
00019 #include <iostream>
00020 #include <TH3D.h>
00021 #include <TStyle.h>
00022 #include <TArrayF.h>
00023 #include <TParticle.h>
00024 #include <TCanvas.h>
00025 
00036 class CheckAlign : public AliFMDInput
00037 {
00038 public:
00039   CheckAlign()
00040   {
00041     AddLoad(kHits);
00042     AddLoad(kDigits);
00043     AddLoad(kGeometry);
00044     AliFMDGeometry* geom = AliFMDGeometry::Instance();
00045     geom->Init();
00046     // geom->InitTransformations();
00047     Double_t iR  = geom->GetRing('I')->GetHighR()+5;
00048     Double_t oR  = geom->GetRing('O')->GetHighR()+5;
00049     Double_t z1l = geom->GetDetector(1)->GetInnerZ()-5;
00050     Double_t z1h = z1l+10;
00051     Double_t z2l = geom->GetDetector(2)->GetOuterZ()-5;
00052     Double_t z2h = geom->GetDetector(2)->GetInnerZ()+5;
00053     Double_t z3l = geom->GetDetector(3)->GetOuterZ()-5;
00054     Double_t z3h = geom->GetDetector(3)->GetInnerZ()+5;
00055     
00056     f1Hits   = new TH3D("hits1", "FMD1 hits", 
00057                         300,-iR,iR,300,-iR,iR,100,z1l,z1h);
00058     f1Hits->SetMarkerColor(2);
00059     f1Hits->SetMarkerStyle(3);
00060       
00061     f2Hits   = new TH3D("hits2", "FMD2 hits", 
00062                         300,-oR,oR,300,-oR,oR,100,z2l,z2h);
00063     f2Hits->SetMarkerColor(2);
00064     f2Hits->SetMarkerStyle(3);
00065 
00066     f3Hits   = new TH3D("hits3", "FMD3 hits", 
00067                         300,-oR,oR,300,-oR,oR,100,z3l,z3h);
00068     f3Hits->SetMarkerColor(2);
00069     f3Hits->SetMarkerStyle(3);
00070 
00071     f1Digits = new TH3D("digits1", "FMD1 digits", 
00072                         300,-iR,iR,300,-iR,iR,100,z1l,z1h); 
00073     f1Digits->SetMarkerColor(4);
00074     f1Digits->SetMarkerStyle(2);
00075 
00076     f2Digits = new TH3D("digits2", "FMD2 digits", 
00077                         300,-oR,oR,300,-oR,oR,100,z2l,z2h);
00078     f2Digits->SetMarkerColor(4);
00079     f2Digits->SetMarkerStyle(2);
00080 
00081     f3Digits = new TH3D("digits3", "FMD3 hits", 
00082                         300,-oR,oR,300,-oR,oR,100,z3l,z3h);
00083     f3Digits->SetMarkerColor(4);
00084     f3Digits->SetMarkerStyle(2);
00085   }
00086   Bool_t Init() 
00087   {
00088     Bool_t ret = AliFMDInput::Init();
00089     AliFMDGeometry* geom = AliFMDGeometry::Instance();
00090     geom->Init();
00091     geom->InitTransformations();
00092     AliFMDParameters* param = AliFMDParameters::Instance();
00093     param->Init();
00094     return ret;
00095   }
00096   
00097   Bool_t ProcessHit(AliFMDHit* hit, TParticle*)
00098   {
00099     // Cache the energy loss 
00100     if (!hit) return kFALSE;
00101     UShort_t det = hit->Detector();
00102     UShort_t str = hit->Strip();
00103     if (str > 511) {
00104       AliWarning(Form("Bad strip number %d in hit", str));
00105       return kTRUE;
00106     }
00107     switch (det) {
00108     case 1: f1Hits->Fill(hit->X(), hit->Y(), hit->Z()); break;
00109     case 2: f2Hits->Fill(hit->X(), hit->Y(), hit->Z()); break;
00110     case 3: f3Hits->Fill(hit->X(), hit->Y(), hit->Z()); break;
00111     }
00112     return kTRUE;
00113   }
00114   Bool_t ProcessDigit(AliFMDDigit* digit)
00115   {
00116     // Cache the energy loss 
00117     if (!digit) return kFALSE;
00118     UShort_t det = digit->Detector();
00119     Char_t   rng = digit->Ring();
00120     UShort_t sec = digit->Sector();
00121     UShort_t str = digit->Strip();
00122     if (str > 511) {
00123       AliWarning(Form("Bad strip number %d in digit", str));
00124       return kTRUE;
00125     }
00126     AliFMDParameters* param = AliFMDParameters::Instance();
00127     if (digit->Counts() < (param->GetPedestal(det, rng, sec, str) + 4 *
00128                            param->GetPedestalWidth(det, rng, sec, str)))
00129       return kTRUE;
00130     AliFMDGeometry* geom = AliFMDGeometry::Instance();
00131     Double_t x, y, z;
00132     geom->Detector2XYZ(det, rng, sec, str, x, y, z);
00133     switch (det) {
00134     case 1: f1Digits->Fill(x, y , z); break;
00135     case 2: f2Digits->Fill(x, y , z); break;
00136     case 3: f3Digits->Fill(x, y , z); break;
00137     }
00138     return kTRUE;
00139   }
00140   //__________________________________________________________________
00141   Bool_t Finish()
00142   {
00143     gStyle->SetPalette(1);
00144     gStyle->SetOptTitle(0);
00145     gStyle->SetCanvasColor(0);
00146     gStyle->SetCanvasBorderSize(0);
00147     gStyle->SetPadColor(0);
00148     gStyle->SetPadBorderSize(0);
00149 
00150     TCanvas* c1 = new TCanvas("FMD1","FMD1");
00151     c1->cd();
00152     f1Hits->Draw();
00153     f1Digits->Draw("same");
00154 
00155     TCanvas* c2 = new TCanvas("FMD2","FMD2");
00156     c2->cd();
00157     f2Hits->Draw();
00158     f2Digits->Draw("same");
00159 
00160     TCanvas* c3 = new TCanvas("FMD3","FMD3");
00161     c3->cd();
00162     f3Hits->Draw();
00163     f3Digits->Draw("same");
00164     return kTRUE;
00165   }
00166 protected:
00167   TH3D* f1Hits;
00168   TH3D* f2Hits;
00169   TH3D* f3Hits;
00170   TH3D* f1Digits;
00171   TH3D* f2Digits;
00172   TH3D* f3Digits;
00173 };
00174 
00175 
00176 //____________________________________________________________________
00177 //
00178 // EOF
00179 //
00180 
00181 
00182   
00183   

Generated on Fri Mar 24 17:11:21 2006 for ALICE FMD Off-line by  doxygen 1.4.6