PreXtCalib Class Reference

#include <PreXtCalib.h>

Inheritance diagram for PreXtCalib:

CalibBase List of all members.

Public Member Functions

 PreXtCalib ()
 ~PreXtCalib ()
void init (TObjArray *hlist, MdcCosGeom *pGeom)
void mergeHist (TFile *fhist)
void calib (MdcCalibConst *calconst, TObjArray *newXtList, TObjArray *r2tList)

Private Member Functions

void renameHist ()

Static Private Member Functions

static Double_t xtfun (Double_t *x, Double_t *par)

Private Attributes

MdcCosGeomm_pGeom
TFolder * m_fdPreXt
TFolder * m_fdNhit
TGraph * m_grXt [NLAYER]
TH2F * m_haxis
TH1F * m_htrec [NLAYER]
TH1F * m_nhitBin [NLAYER]
TH1F * m_nhitTot
double m_tbin [NPreXtBin]

Static Private Attributes

static const int NPreXtBin = 40

Detailed Description

Definition at line 12 of file PreXtCalib.h.


Constructor & Destructor Documentation

PreXtCalib::PreXtCalib (  ) 

Definition at line 8 of file PreXtCalib.cpp.

00008                       {
00009      cout << "Calibration type: PreXtCalib" << endl;
00010 }

PreXtCalib::~PreXtCalib (  ) 

Definition at line 12 of file PreXtCalib.cpp.

00012                        {
00013 }


Member Function Documentation

void PreXtCalib::calib ( MdcCalibConst calconst,
TObjArray *  newXtList,
TObjArray *  r2tList 
) [virtual]

Implements CalibBase.

Definition at line 66 of file PreXtCalib.cpp.

References bin, check_raw_filter::dist, MdcCalibConst::getXtpar(), gFgCalib, m_fdPreXt, m_grXt, m_nhitBin, m_nhitTot, m_pGeom, m_tbin, NENTRXT, NLAYER, NLR, NPreXtBin, pi, renameHist(), and xtfun().

00066                                                                                        {
00067      double pi = 3.141592653;
00068      double dist[NPreXtBin];
00069      double xtpar[6];
00070      char hname[200];
00071 
00072      TF1* funXt = new TF1("funXt", xtfun, 0, 300, 6);
00073      funXt -> FixParameter(0, 0.0);
00074      funXt -> SetParameter(1, 0.03);
00075      funXt -> SetParameter(2, 0.0);
00076      funXt -> SetParameter(3, 0.0);
00077      funXt -> SetParameter(4, 0.0);
00078      funXt -> SetParameter(5, 0.0);
00079 
00080      ofstream fxtlog("preXtpar.dat");
00081      for(int lay=0; lay<NLAYER; lay++){
00082           sprintf(hname, "mgrPreXt%02d", lay);
00083           m_grXt[lay] = new TGraph();
00084           m_grXt[lay] -> SetName(hname);
00085           m_grXt[lay] -> SetMarkerStyle(20);
00086           m_fdPreXt -> Add(m_grXt[lay]);
00087 
00088           double layRad = m_pGeom -> getLayer(lay) -> getLayerRad();
00089           int ncel = m_pGeom -> getLayer(lay) -> getNcell();
00090           double dm = pi * layRad / (double)ncel;
00091           Double_t nTot = m_nhitTot->GetBinContent(lay+1);
00092           double tm = calconst->getXtpar(lay, 0, 0, 6);
00093 
00094           fxtlog << "layer " << lay << endl;
00095           for(int bin=0; bin<NPreXtBin; bin++){
00096                Double_t nhitBin = m_nhitBin[lay]->GetBinContent(bin+1);
00097                dist[bin] = dm * nhitBin / nTot;
00098                m_grXt[lay] -> SetPoint(bin, m_tbin[bin], dist[bin]);
00099                fxtlog << setw(4) << bin << setw(15) << m_tbin[bin]
00100                       << setw(15) << dist[bin] << setw(15) << dm
00101                       << setw(10) << nhitBin
00102                       << setw(10) << nTot << endl;
00103 
00104                if(m_tbin[bin] >= tm) break;
00105           }
00106 
00107           if(1 == gFgCalib[lay]){
00108                m_grXt[lay] -> Fit(funXt, "Q", "", 0.0, tm);
00109                for(int ord=0; ord<6; ord++) xtpar[ord] = funXt -> GetParameter(ord);
00110 
00111                for(int iEntr=0; iEntr<NENTRXT; iEntr++){
00112                     for(int iLR=0; iLR<NLR; iLR++){
00113                          for(int ord=0; ord<6; ord++){
00114                               calconst -> resetXtpar(lay, iEntr, iLR, ord, xtpar[ord]);
00115                          }
00116                     }
00117                }
00118           } else{
00119                for(int ord=0; ord<6; ord++) xtpar[ord] = calconst->getXtpar(lay, 0, 0, ord);
00120           }
00121 
00122           for(int ord=0; ord<6; ord++)  fxtlog << setw(14) << xtpar[ord];
00123           fxtlog << setw(10) << tm << "  0" << endl;
00124      } // end of layer loop
00125      fxtlog.close();
00126      cout << "preXt.dat was written." << endl;
00127 
00128      renameHist();
00129      delete funXt;
00130 }

void PreXtCalib::init ( TObjArray *  hlist,
MdcCosGeom pGeom 
) [virtual]

Implements CalibBase.

Definition at line 15 of file PreXtCalib.cpp.

References bin, m_fdNhit, m_fdPreXt, m_haxis, m_htrec, m_nhitBin, m_nhitTot, m_pGeom, m_tbin, NLAYER, and NPreXtBin.

00015                                                         {
00016      m_pGeom = pGeom;
00017 
00018      double twid = 10.0;        // ns
00019      for(int bin=0; bin<NPreXtBin; bin++)  m_tbin[bin] = (double)(bin+1) * twid;
00020 
00021      m_fdPreXt = new TFolder("mPreXt", "PreXt");
00022      hlist->Add(m_fdPreXt);
00023 
00024      m_fdNhit = new TFolder("mXtNhit", "XtNhit");
00025      hlist->Add(m_fdNhit);
00026 
00027      m_haxis = new TH2F("maxis", "", 50, 0, 300, 50, 0, 9);
00028      m_haxis -> SetStats(0);
00029      m_fdPreXt -> Add(m_haxis);
00030 
00031      m_nhitTot = new TH1F("mnhitTot", "", 43, -0.5, 42.5);
00032      m_fdNhit -> Add(m_nhitTot);
00033 
00034      char hname[200];
00035      for(int lay=0; lay<NLAYER; lay++){
00036           sprintf(hname, "mtrec%02d", lay);
00037           m_htrec[lay] = new TH1F(hname, "", 310, -20, 600);
00038           m_fdPreXt -> Add(m_htrec[lay]);
00039 
00040           sprintf(hname, "mnhitBin%02d", lay);
00041           m_nhitBin[lay] = new TH1F(hname, "", 40, 5.0, 405.0);
00042           m_fdNhit -> Add(m_nhitBin[lay]);
00043      }
00044 }

void PreXtCalib::mergeHist ( TFile *  fhist  )  [virtual]

Implements CalibBase.

Definition at line 46 of file PreXtCalib.cpp.

References m_htrec, m_nhitBin, m_nhitTot, and NLAYER.

00046                                       {
00047      TFolder* fdPreXt = (TFolder*)fhist->Get("PreXt");
00048      TFolder* fdNhit = (TFolder*)fhist->Get("XtNhit");
00049 
00050      char hname[200];
00051      TH1F* hist;
00052      for(int lay=0; lay<NLAYER; lay++){
00053           sprintf(hname, "trec%02d", lay);
00054           hist = (TH1F*)fdPreXt->FindObjectAny(hname);
00055           m_htrec[lay]->Add(hist);
00056 
00057           sprintf(hname, "nhitBin%02d", lay);
00058           hist = (TH1F*)fdNhit->FindObjectAny(hname);
00059           m_nhitBin[lay]->Add(hist);
00060      }
00061 
00062      hist = (TH1F*)fdNhit->FindObjectAny("nhitTot");
00063      m_nhitTot->Add(hist);
00064 }

void PreXtCalib::renameHist (  )  [private]

Reimplemented from CalibBase.

Definition at line 140 of file PreXtCalib.cpp.

References m_fdNhit, m_fdPreXt, m_grXt, m_haxis, m_htrec, m_nhitBin, m_nhitTot, and NLAYER.

Referenced by calib().

00140                            {
00141      m_fdPreXt->SetName("PreXt");
00142      m_fdNhit->SetName("XtNhit");
00143      m_haxis->SetName("axis");
00144      m_nhitTot->SetName("nhitTot");
00145 
00146      char hname[200];
00147      for(int lay=0; lay<NLAYER; lay++){
00148           sprintf(hname, "trec%02d", lay);
00149           m_htrec[lay]->SetName(hname);
00150 
00151           sprintf(hname, "nhitBin%02d", lay);
00152           m_nhitBin[lay]->SetName(hname);
00153 
00154           sprintf(hname, "grPreXt%02d", lay);
00155           m_grXt[lay] -> SetName(hname);
00156      }
00157 }

Double_t PreXtCalib::xtfun ( Double_t *  x,
Double_t *  par 
) [static, private]

Definition at line 132 of file PreXtCalib.cpp.

Referenced by calib().

00132                                                     {
00133      Double_t val = 0.0;
00134      for(Int_t ord=0; ord<6; ord++){
00135           val += par[ord] * pow(x[0], ord);
00136      }
00137      return val;
00138 }


Member Data Documentation

TFolder* PreXtCalib::m_fdNhit [private]

Definition at line 28 of file PreXtCalib.h.

Referenced by init(), and renameHist().

TFolder* PreXtCalib::m_fdPreXt [private]

Definition at line 27 of file PreXtCalib.h.

Referenced by calib(), init(), and renameHist().

TGraph* PreXtCalib::m_grXt[NLAYER] [private]

Definition at line 29 of file PreXtCalib.h.

Referenced by calib(), and renameHist().

TH2F* PreXtCalib::m_haxis [private]

Definition at line 30 of file PreXtCalib.h.

Referenced by init(), and renameHist().

TH1F* PreXtCalib::m_htrec[NLAYER] [private]

Definition at line 31 of file PreXtCalib.h.

Referenced by init(), mergeHist(), and renameHist().

TH1F* PreXtCalib::m_nhitBin[NLAYER] [private]

Definition at line 32 of file PreXtCalib.h.

Referenced by calib(), init(), mergeHist(), and renameHist().

TH1F* PreXtCalib::m_nhitTot [private]

Definition at line 33 of file PreXtCalib.h.

Referenced by calib(), init(), mergeHist(), and renameHist().

MdcCosGeom* PreXtCalib::m_pGeom [private]

Definition at line 26 of file PreXtCalib.h.

Referenced by calib(), and init().

double PreXtCalib::m_tbin[NPreXtBin] [private]

Definition at line 35 of file PreXtCalib.h.

Referenced by calib(), and init().

const int PreXtCalib::NPreXtBin = 40 [static, private]

Definition at line 24 of file PreXtCalib.h.

Referenced by calib(), and init().


Generated on Tue Nov 29 23:20:43 2016 for BOSS_7.0.2 by  doxygen 1.4.7