Changeset 1827 for trunk/MagicSoft/Mars


Ignore:
Timestamp:
03/18/03 18:13:32 (22 years ago)
Author:
moralejo
Message:
*** empty log message ***
Location:
trunk/MagicSoft/Mars
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/MagicSoft/Mars/Changelog

    r1826 r1827  
    1616    * mhist/MHMcCT1CollectionArea.[h,cc]
    1717      - Added arguments in constructor: number of bins and ranges of the
    18         x-axis (energy) of the 2-d histograms.
     18        x-axis (energy) of the 2-d histograms. Changed type of binning:
     19        now the x-axis is log10(energy) and bins have equal width.
    1920
    2021    * macros/CT1collarea.C
    2122      - The MHMcCT1CollectionArea object is now created and added to the
    22         parlist so that  we can choose the binning.
    23 
     23        parlist so that  we can choose the binning. Changed the way
     24        histograms are written to the output file.
    2425
    2526
  • trunk/MagicSoft/Mars/mhist/MHMcCT1CollectionArea.cc

    r1824 r1827  
    5858
    5959  fName  = name  ? name  : "MHMcCT1CollectionArea";
    60   fTitle = title ? title : "Collection Area vs. Energy";
     60  fTitle = title ? title : "Collection Area vs. log10 Energy";
    6161
    6262  fHistAll = new TH2D;
     
    7171
    7272  fHistCol->SetTitle(fTitle);
    73   fHistAll->SetTitle("All showers - Theta vs Energy distribution");
    74   fHistSel->SetTitle("Selected showers - Theta vs Energy distribution");
     73  fHistAll->SetTitle("All showers - Theta vs log10 Energy distribution");
     74  fHistSel->SetTitle("Selected showers - Theta vs log10 Energy distribution");
    7575
    7676  fHistAll->SetDirectory(NULL);
     
    7878  fHistCol->SetDirectory(NULL);
    7979
    80   fHistAll->SetXTitle("E [GeV]");
     80  fHistAll->SetXTitle("log10 E [GeV]");
    8181  fHistAll->SetYTitle("theta [deg]");
    8282  fHistAll->SetZTitle("N");
    8383
    84   fHistSel->SetXTitle("E [GeV]");
     84  fHistSel->SetXTitle("log10 E [GeV]");
    8585  fHistSel->SetYTitle("theta [deg]");
    86   fHistSel->SetYTitle("N");
    87 
    88   fHistCol->SetXTitle("E [GeV]");
     86  fHistSel->SetZTitle("N");
     87
     88  fHistCol->SetXTitle("log10 E [GeV]");
    8989  fHistCol->SetYTitle("theta [deg]");
    9090  fHistCol->SetZTitle("A [m^{2}]");
     
    9999{
    100100  MBinning binsx;
    101   binsx.SetEdgesLog(nbins, minEnergy, maxEnergy);
     101  binsx.SetEdges(nbins, minEnergy, maxEnergy);
    102102
    103103  MBinning binsy;
     
    129129void MHMcCT1CollectionArea::FillSel(Double_t energy, Double_t theta)
    130130{
    131   fHistSel->Fill(energy, theta);
     131  fHistSel->Fill(log10(energy), theta);
    132132}
    133133
     
    143143  fHistAll->Draw(option);
    144144
    145   gPad->SetLogx();
    146 
    147145  gPad->Modified();
    148146  gPad->Update();
     
    159157
    160158  fHistSel->Draw(option);
    161 
    162   gPad->SetLogx();
    163159
    164160  gPad->Modified();
     
    183179  fHistCol->DrawCopy(option);
    184180
    185   gPad->SetLogx();
    186 
    187181  c->Modified();
    188182  c->Update();
     
    197191
    198192  fHistCol->Draw(option);
    199 
    200   gPad->SetLogx();
    201193
    202194  gPad->Modified();
     
    306298      for (Int_t i=1; i <= fHistAll->GetNbinsX(); i++)
    307299        {
    308           const Float_t e1 = fHistAll->GetXaxis()->GetBinLowEdge(i);
    309           const Float_t e2 = fHistAll->GetXaxis()->GetBinLowEdge(i+1);
     300          const Float_t e1 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i));
     301          const Float_t e2 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i+1));
    310302
    311303          if (e1 < emin1 || e2 > emax2)
     
    336328
    337329          if (Na <= 0)
    338             continue;
     330            {
     331              //
     332              // If energy is large, this case means that no or very few events
     333              // were generated at this energy bin. In this case we assign it
     334              // the effective area of the bin below it in energy. If energy is
     335              // below 1E4, it means that no events triggered -> eff area = 0
     336              //
     337
     338              if (fHistSel->GetXaxis()->GetBinLowEdge(ix) > 4.)
     339                {
     340                  fHistCol->SetBinContent(ix, thetabin, fHistCol->GetBinContent(ix-1, thetabin));
     341                  fHistCol->SetBinError(ix, thetabin, fHistCol->GetBinError(ix-1, thetabin));
     342                }
     343              continue;
     344            }
    339345
    340346          const Float_t Ns = fHistSel->GetBinContent(ix,thetabin);
     
    350356          const Double_t err = sqrt((1.-eff)*Ns)/Na;
    351357
     358
    352359          const Float_t area = dr * cos(theta*TMath::Pi()/180.);
    353360
    354361          fHistCol->SetBinContent(ix, thetabin, eff*area);
    355362          fHistCol->SetBinError(ix, thetabin, err*area);
     363
    356364        }
    357365    }
  • trunk/MagicSoft/Mars/mhist/MHMcCT1CollectionArea.h

    r1824 r1827  
    1818
    1919public:
    20     MHMcCT1CollectionArea(const char *name=NULL, const char *title=NULL, Int_t nbins=45, Axis_t minEnergy=100., Axis_t maxEnergy=30000.);
     20    MHMcCT1CollectionArea(const char *name=NULL, const char *title=NULL, Int_t nbins=30, Axis_t minEnergy=2., Axis_t maxEnergy=5.);
    2121    ~MHMcCT1CollectionArea();
    2222
     
    2727
    2828    const TH2D *GetHist() const { return fHistCol; }
     29    const TH2D *GetHAll() const { return fHistAll; }
     30    const TH2D *GetHSel() const { return fHistSel; }
    2931
    3032    void Draw(Option_t *option="");
Note: See TracChangeset for help on using the changeset viewer.