]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Modified plots and made jet finder use SDigits
authormhorner <mhorner@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 18 May 2004 16:33:27 +0000 (16:33 +0000)
committermhorner <mhorner@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 18 May 2004 16:33:27 +0000 (16:33 +0000)
EMCAL/AliEMCALJetFinderInputSimPrep.cxx
EMCAL/AliEMCALJetFinderInputSimPrep.h
EMCAL/AliEMCALJetFinderPlots.cxx
EMCAL/AliEMCALJetFinderPlots.h

index cfcfc52a54bb89d044f4901b1df3604a0ec80134..4940c519de460e56ea21c7b4e90adb2ceea18559 100755 (executable)
@@ -110,17 +110,32 @@ if (fDebug > 1) Info("Reset","Beginning Reset");
 
 }
 
-Int_t AliEMCALJetFinderInputSimPrep::FillFromFile(TString *filename, AliEMCALJetFinderFileType_t filetype,Int_t EventNumber)
+Int_t AliEMCALJetFinderInputSimPrep::FillFromFile(TString *filename, AliEMCALJetFinderFileType_t filetype,Int_t EventNumber,TString data="H")
 {
 if (fDebug > 1) Info("FillFromFile","Beginning FillFromFile");
    fFileType = filetype;       
    AliEMCALGetter *gime = AliEMCALGetter::Instance(filename->Data());
    if (fDebug > 1) Info("FillFromFile","Instantiated Getter with %s as the file",filename->Data());
-   gime->Event(EventNumber,"XH") ;
+   if (data.Contains("S"))
+   {
+          gime->Event(EventNumber,"XS") ;
+   }else
+   {
+          gime->Event(EventNumber,"XH") ;
+   }
    if (fDebug > 1) Info("FillFromFile","Got event %i with option \"XH\"",EventNumber);
 
    if (        fEMCALType == kHits ||
-       fEMCALType == kTimeCut )  FillHits();
+       fEMCALType == kTimeCut )  
+   {
+          if (data.Contains("S"))
+          {
+                  FillSDigits();
+          }else
+          {
+                   FillHits();
+          }
+   }
    if ( fTrackType != kNoTracks  )  FillTracks();      
    if ( fFileType  != kData){
           FillPartons();
@@ -574,6 +589,18 @@ void AliEMCALJetFinderInputSimPrep::FillDigits()
 {
        // Fill digits to input object
 
+}
+void AliEMCALJetFinderInputSimPrep::FillSDigits()  
+{
+Info("FillSDigits","Beginning FillSDigits");
+       AliEMCALGetter *gime = AliEMCALGetter::Instance();
+        
+       // Fill digits to input object
+       for (Int_t towerid=0; towerid < gime->SDigits()->GetEntries(); towerid++)
+       {
+               fInputObject.AddEnergyToDigit(gime->SDigit(towerid)->GetId(), gime->SDigit(towerid)->GetAmp());
+       }
+
 }
 
 void AliEMCALJetFinderInputSimPrep::Smear(TParticle *particle)
index a120c6acd1bb35865a0be26b91730dbe6dac7d82..41ddfabd32b89def6e8f3f82ed21be1f0414d5ec 100755 (executable)
@@ -26,14 +26,15 @@ class AliEMCALJetFinderInputSimPrep : public AliEMCALJetFinderInputPrep
        void Reset(AliEMCALJetFinderResetType_t resettype);
         void SetEMCALType(AliEMCALJetFinderEMCALType_t emcaltype )  {fEMCALType = emcaltype;}
        //void SetDebug(Int_t debug = 0)  {fDebug = debug;}
-        void SetSmearingType(AliEMCALJetFinderSmearingType_t smeartype )  {fSmearType = smeartype;}
+        void SetSmearingType(AliEMCALJetFinderSmearingType_t smeartype ) {fSmearType = smeartype;}
        void SetTrackType(AliEMCALJetFinderTrackType_t tracktype){fTrackType = tracktype;}  
-       void SetEfficiency(Float_t efficiency)  {fEfficiency = efficiency; }
+       void SetEfficiency(Float_t efficiency)  {fEfficiency = efficiency;}
        void SetTimeCut(Float_t timecut)  {fTimeCut = timecut; fEMCALType = kTimeCut;}
-       Int_t FillFromFile(TString *filename, AliEMCALJetFinderFileType_t filetype,Int_t EventNumber);
+       Int_t FillFromFile(TString *filename, AliEMCALJetFinderFileType_t filetype,Int_t EventNumber,TString data);
        AliEMCALJetFinderInput* GetJetFinderInput(){return &fInputObject;}
        private:
        void FillHits();                // Fill from the hits to input object from simulation
+       void FillSDigits();             // Fill from the hits to input object from simulation
        void FillTracks();              // Fill from particles simulating a TPC to input object from simulation
        void Smear(TParticle *particle);
        Bool_t Efficiency();
index edefd4504f535bf72398b4e7f98126f67b021942..ea445127429fd8a7218f37e3b246859a545d81b6 100755 (executable)
@@ -77,6 +77,7 @@ fhInputOutput=0;
 //     TH2F                            *fhInputOutput;  //("hJetEtRatio2","Ratio of Second Highest to Highest",100,0,1);
      
 fhRecoBinPt=0;        // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
+fhRecoBinPtNoBg=0;            // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
 fhRecoBinPartonPt=0;    // ("fhRecoBinPartonPt","Input Pt Distribution",100,0,1);
 fhRecoBinJetEt=0;       // ("fhRecoJetEt","E_{T}^{reco}",250,0.,250.);
 fhRecoBinInputJetEt=0;  // ("fhRecoInputJetEt","E_{T}^{reco}",250,0.,250.);
@@ -221,6 +222,8 @@ fhInputOutput=      new TH2F("hInputOutput","Input and Reconstruction Correlations;In
 
 fhRecoBinPt =new TH1F("fhRecoBinPt","Reconstructed Pt Distribution",200,0,200);
 fhRecoBinPt->Sumw2();
+fhRecoBinPtNoBg =new TH1F("fhRecoBinPtNoBg","Reconstructed Pt Distribution Background Subtracted",200,0,200);
+fhRecoBinPtNoBg->Sumw2();
 fhRecoBinPartonPt = new TH1F("fhRecoBinPartonPt","Input Pt Distribution",200,0,200);
 fhRecoBinPartonPt->Sumw2();
 fhRecoBinFragmFcn =new TH1F("fhRecoBinFragmFcn","Reconstructed Frag. Fcn",200,0,2);
@@ -288,6 +291,7 @@ delete        fhEtaPhiSpread;
        delete                          fhEtaPhiDist2;  //("hEtaPhiDist2","Angular Distance Between First and Second",100,0,3);
 
        delete fhRecoBinPt;          // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
+       delete fhRecoBinPtNoBg;          // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
        delete fhRecoBinPartonPt;    // ("fhRecoBinPartonPt","Input Pt Distribution",100,0,1);
        delete fhRecoBinJetEt;       // ("fhRecoJetEt","E_{T}^{reco}",250,0.,250.);
        delete fhRecoBinInputJetEt;  // ("fhRecoInputJetEt","E_{T}^{reco}",250,0.,250.);
@@ -310,6 +314,15 @@ void AliEMCALJetFinderPlots::FillFromOutput(AliEMCALJetFinderOutput* output)
 if (!fInitialised) InitPlots();        
   fOutput = output;
   if (!fOutput) return;
+// Make some temporary histograms to make sure we subtract 
+  // background properly
+/*
+tempFragmFcnNoBg =new TH1F("tempFragmFcnNoBg","Reconstructed Frag. Fcn With Background Removed",200,0,2);
+tempPtNoBg =new TH1F("tempPtNoBg","Reconstructed Frag. Fcn With Background Removed",200,0,200);
+tempFragmFcnNoBg->Fill(count/(jethighest->Energy()*fScaleFactor),-fhBackHisto->GetBinContent(count));
+tempPtNoBg->AddBinContent(count,-fhBackHisto->GetBinContent(count));
+*/
+  
   fhNJets->Fill(fOutput->GetNJets());
   Bool_t doesJetMeetBinCriteria = 0;
   AliEMCALJet* jethighest=0; 
@@ -519,6 +532,7 @@ if (numappjet > 1)
          if (doesJetMeetBinCriteria)
          {
                  fhRecoBinPt->Fill(correctp*sin(tt));          // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
+                 fhRecoBinPtNoBg->Fill(correctp*sin(tt));          // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
                  fhRecoBinFragmFcn->Fill(  correctp*sin(tt)/(jethighest->Energy()*fScaleFactor)  ); // This is the jet fragmentation function
                  fhRecoBinFragmFcnNoBg->Fill(  correctp*sin(tt)/(jethighest->Energy()*fScaleFactor)  ); // This is the jet fragmentation function
          }
@@ -636,7 +650,7 @@ if (numappjet > 1)
          Double_t correctp = pt[iT]/stt;
          fhJetPL->Fill( correctp*cos(alpha));      
          if ( (jet->Eta()-eta[iT])*(jet->Eta()-eta[iT]) +  
-                         (jet->Phi()-phi[iT])*(jet->Phi()-phi[iT]) < 0.2*0.2 )   
+              (jet->Phi()-phi[iT])*(jet->Phi()-phi[iT]) < 0.2*0.2 )   
                  fhJetJT->Fill( correctp*sin(alpha));
          fhJetPT->Fill(correctp*sin(tt));
          if (fNominalEnergy==0.0){
@@ -648,17 +662,19 @@ if (numappjet > 1)
          if (doesJetMeetBinCriteria)
          {
                  fhRecoBinPt->Fill(correctp*sin(tt));          // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
+                 fhRecoBinPtNoBg->Fill(correctp*sin(tt));          // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
                  fhRecoBinFragmFcn->Fill(  correctp*sin(tt)/(jet->Energy()*fScaleFactor)  ); // This is the jet fragmentation function
                  fhRecoBinFragmFcnNoBg->Fill(  correctp*sin(tt)/(jet->Energy()*fScaleFactor)  ); // This is the jet fragmentation function
          }
   }// loop over tracks
   }
 
-if (numappjet>=1 && fhBackHisto != 0)
+if (numappjet>=1 && fhBackHisto != 0 && doesJetMeetBinCriteria)
   {
     for (Int_t count=1;count<=100;count++)
     {
-       fhRecoBinFragmFcnNoBg->AddBinContent(count,-fhBackHisto->GetBinContent(count));
+       fhRecoBinFragmFcnNoBg->Fill( ((Float_t)count)/(jethighest->Energy()*fScaleFactor),-fhBackHisto->GetBinContent(count));
+       fhRecoBinPtNoBg->AddBinContent(count,-fhBackHisto->GetBinContent(count));
     } 
   }
 
index 305980e5418b1599602e57712859b15cbc700dd2..8837d197020d30f8594e63c5e41fef1d03893a91 100755 (executable)
@@ -93,6 +93,7 @@ class AliEMCALJetFinderPlots : public TObject
        //============================== Reconstruction Bin Comparison  ============================================
        
        TH1F* GetRecoBinPt(){return fhRecoBinPt;}                  // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
+       TH1F* GetRecoBinPtNoBg(){return fhRecoBinPtNoBg;}                  // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
        TH1F* GetRecoBinPartonPt(){return fhRecoBinPartonPt;}      // ("fhRecoBinPartonPt","Input Pt Distribution",100,0,1);
        TH1F* GetRecoBinJetEt(){return fhRecoBinJetEt;}            // ("fhRecoJetEt","E_{T}^{reco}",250,0.,250.);
        TH1F* GetRecoBinInputJetEt(){return fhRecoBinInputJetEt;}  // ("fhRecoInputJetEt","E_{T}^{reco}",250,0.,250.);
@@ -152,6 +153,7 @@ class AliEMCALJetFinderPlots : public TObject
        //============================== Reconstruction Bin Comparison  ============================================
        
        TH1F                            *fhRecoBinPt;          // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
+       TH1F                            *fhRecoBinPtNoBg;              // ("fhRecoBinPt","Reconstructed Pt Distribution",100,0,1);
        TH1F                            *fhRecoBinPartonPt;    // ("fhRecoBinPartonPt","Input Pt Distribution",100,0,1);
        TH1F                            *fhRecoBinJetEt;       // ("fhRecoJetEt","E_{T}^{reco}",250,0.,250.);
        TH1F                            *fhRecoBinInputJetEt;  // ("fhRecoInputJetEt","E_{T}^{reco}",250,0.,250.);