BOSS 7.0.4
BESIII Offline Software System
Loading...
Searching...
No Matches
PreXtCalib Class Reference

#include <PreXtCalib.h>

+ Inheritance diagram for PreXtCalib:

Public Member Functions

 PreXtCalib ()
 
 ~PreXtCalib ()
 
void init (TObjArray *hlist, MdcCosGeom *pGeom)
 
void mergeHist (TFile *fhist)
 
void calib (MdcCalibConst *calconst, TObjArray *newXtList, TObjArray *r2tList)
 
- Public Member Functions inherited from CalibBase
 CalibBase ()
 
virtual ~CalibBase ()
 
virtual void init (TObjArray *hlist, MdcCosGeom *pGeom)=0
 
virtual void mergeHist (TFile *fhist)=0
 
virtual void calib (MdcCalibConst *calconst, TObjArray *newXtList, TObjArray *r2tList)=0
 

Detailed Description

Definition at line 12 of file PreXtCalib.h.

Constructor & Destructor Documentation

◆ PreXtCalib()

PreXtCalib::PreXtCalib ( )

Definition at line 8 of file PreXtCalib.cpp.

8 {
9 cout << "Calibration type: PreXtCalib" << endl;
10}

◆ ~PreXtCalib()

PreXtCalib::~PreXtCalib ( )

Definition at line 12 of file PreXtCalib.cpp.

12 {
13}

Member Function Documentation

◆ calib()

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

Implements CalibBase.

Definition at line 66 of file PreXtCalib.cpp.

66 {
67 double pi = 3.141592653;
68 double dist[NPreXtBin];
69 double xtpar[6];
70 char hname[200];
71
72 TF1* funXt = new TF1("funXt", xtfun, 0, 300, 6);
73 funXt -> FixParameter(0, 0.0);
74 funXt -> SetParameter(1, 0.03);
75 funXt -> SetParameter(2, 0.0);
76 funXt -> SetParameter(3, 0.0);
77 funXt -> SetParameter(4, 0.0);
78 funXt -> SetParameter(5, 0.0);
79
80 ofstream fxtlog("preXtpar.dat");
81 for(int lay=0; lay<NLAYER; lay++){
82 sprintf(hname, "mgrPreXt%02d", lay);
83 m_grXt[lay] = new TGraph();
84 m_grXt[lay] -> SetName(hname);
85 m_grXt[lay] -> SetMarkerStyle(20);
86 m_fdPreXt -> Add(m_grXt[lay]);
87
88 double layRad = m_pGeom -> getLayer(lay) -> getLayerRad();
89 int ncel = m_pGeom -> getLayer(lay) -> getNcell();
90 double dm = pi * layRad / (double)ncel;
91 Double_t nTot = m_nhitTot->GetBinContent(lay+1);
92 double tm = calconst->getXtpar(lay, 0, 0, 6);
93
94 fxtlog << "layer " << lay << endl;
95 for(int bin=0; bin<NPreXtBin; bin++){
96 Double_t nhitBin = m_nhitBin[lay]->GetBinContent(bin+1);
97 dist[bin] = dm * nhitBin / nTot;
98 m_grXt[lay] -> SetPoint(bin, m_tbin[bin], dist[bin]);
99 fxtlog << setw(4) << bin << setw(15) << m_tbin[bin]
100 << setw(15) << dist[bin] << setw(15) << dm
101 << setw(10) << nhitBin
102 << setw(10) << nTot << endl;
103
104 if(m_tbin[bin] >= tm) break;
105 }
106
107 if(1 == gFgCalib[lay]){
108 m_grXt[lay] -> Fit(funXt, "Q", "", 0.0, tm);
109 for(int ord=0; ord<6; ord++) xtpar[ord] = funXt -> GetParameter(ord);
110
111 for(int iEntr=0; iEntr<NENTRXT; iEntr++){
112 for(int iLR=0; iLR<NLR; iLR++){
113 for(int ord=0; ord<6; ord++){
114 calconst -> resetXtpar(lay, iEntr, iLR, ord, xtpar[ord]);
115 }
116 }
117 }
118 } else{
119 for(int ord=0; ord<6; ord++) xtpar[ord] = calconst->getXtpar(lay, 0, 0, ord);
120 }
121
122 for(int ord=0; ord<6; ord++) fxtlog << setw(14) << xtpar[ord];
123 fxtlog << setw(10) << tm << " 0" << endl;
124 } // end of layer loop
125 fxtlog.close();
126 cout << "preXt.dat was written." << endl;
127
128 renameHist();
129 delete funXt;
130}
gr Fit("g1","R")
*******INTEGER m_nBinMax INTEGER m_NdiMax !No of bins in histogram for cell exploration division $ !Last vertex $ !Last active cell $ !Last cell in buffer $ !No of sampling when dividing cell $ !No of function total $ !Flag for random ceel for $ !Flag for type of for WtMax $ !Flag which decides whether vertices are included in the sampling $ entire domain is hyp !Maximum effective eevents per bin
Definition: FoamA.h:85
double getXtpar(int lay, int entr, int lr, int order)
sprintf(cut,"kal_costheta0_em>-0.93&&kal_costheta0_em<0.93&&kal_pxy0_em>=0.05+%d*0.1&&kal_pxy0_em<0.15+%d*0.1&&NGch>=2", j, j)
gr1 SetMarkerStyle(8)
mg Add(gr3)

◆ init()

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

Implements CalibBase.

Definition at line 15 of file PreXtCalib.cpp.

15 {
16 m_pGeom = pGeom;
17
18 double twid = 10.0; // ns
19 for(int bin=0; bin<NPreXtBin; bin++) m_tbin[bin] = (double)(bin+1) * twid;
20
21 m_fdPreXt = new TFolder("mPreXt", "PreXt");
22 hlist->Add(m_fdPreXt);
23
24 m_fdNhit = new TFolder("mXtNhit", "XtNhit");
25 hlist->Add(m_fdNhit);
26
27 m_haxis = new TH2F("maxis", "", 50, 0, 300, 50, 0, 9);
28 m_haxis -> SetStats(0);
29 m_fdPreXt -> Add(m_haxis);
30
31 m_nhitTot = new TH1F("mnhitTot", "", 43, -0.5, 42.5);
32 m_fdNhit -> Add(m_nhitTot);
33
34 char hname[200];
35 for(int lay=0; lay<NLAYER; lay++){
36 sprintf(hname, "mtrec%02d", lay);
37 m_htrec[lay] = new TH1F(hname, "", 310, -20, 600);
38 m_fdPreXt -> Add(m_htrec[lay]);
39
40 sprintf(hname, "mnhitBin%02d", lay);
41 m_nhitBin[lay] = new TH1F(hname, "", 40, 5.0, 405.0);
42 m_fdNhit -> Add(m_nhitBin[lay]);
43 }
44}

◆ mergeHist()

void PreXtCalib::mergeHist ( TFile *  fhist)
virtual

Implements CalibBase.

Definition at line 46 of file PreXtCalib.cpp.

46 {
47 TFolder* fdPreXt = (TFolder*)fhist->Get("PreXt");
48 TFolder* fdNhit = (TFolder*)fhist->Get("XtNhit");
49
50 char hname[200];
51 TH1F* hist;
52 for(int lay=0; lay<NLAYER; lay++){
53 sprintf(hname, "trec%02d", lay);
54 hist = (TH1F*)fdPreXt->FindObjectAny(hname);
55 m_htrec[lay]->Add(hist);
56
57 sprintf(hname, "nhitBin%02d", lay);
58 hist = (TH1F*)fdNhit->FindObjectAny(hname);
59 m_nhitBin[lay]->Add(hist);
60 }
61
62 hist = (TH1F*)fdNhit->FindObjectAny("nhitTot");
63 m_nhitTot->Add(hist);
64}

The documentation for this class was generated from the following files: