BOSS 7.0.9
BESIII Offline Software System
Loading...
Searching...
No Matches
PreXtMdcCalib.cxx
Go to the documentation of this file.
2
3#include "GaudiKernel/MsgStream.h"
4#include "GaudiKernel/IMessageSvc.h"
5#include "GaudiKernel/StatusCode.h"
6#include "GaudiKernel/ISvcLocator.h"
7#include "GaudiKernel/Bootstrap.h"
8
9#include "GaudiKernel/SmartDataPtr.h"
10#include "GaudiKernel/IDataProviderSvc.h"
11#include "GaudiKernel/PropertyMgr.h"
12
14#include "EventModel/Event.h"
15#include "MdcRawEvent/MdcDigi.h"
17#include "Identifier/MdcID.h"
18#include "TStyle.h"
19
20#include <iostream>
21#include <fstream>
22#include <iomanip>
23#include <string>
24#include <cstring>
25#include <cmath>
26
27#include "TF1.h"
28
29using namespace Event;
30
32 m_fgIniTm = false;
33}
34
36}
37
39 for(int lay=0; lay<MdcCalNLayer; lay++){
40 delete m_grXt[lay];
41 delete m_htrec[lay];
42 }
43 delete m_haxis;
44 delete m_fdPreXt;
45}
46
47void PreXtMdcCalib::initialize(TObjArray* hlist, IMdcGeomSvc* mdcGeomSvc,
48 IMdcCalibFunSvc* mdcFunSvc, IMdcUtilitySvc* mdcUtilitySvc) {
49 IMessageSvc* msgSvc;
50 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
51 MsgStream log(msgSvc, "PreXtMdcCalib");
52 log << MSG::INFO << "PreXtMdcCalib::initialize()" << endreq;
53
54 m_hlist = hlist;
55 m_mdcGeomSvc = mdcGeomSvc;
56 m_mdcFunSvc = mdcFunSvc;
57 m_mdcUtilitySvc = mdcUtilitySvc;
58
59 int lay;
60
61 m_nWire = m_mdcGeomSvc -> getWireSize();
62 m_nLayer = m_mdcGeomSvc -> getLayerSize();
63
64 m_fdPreXt = new TFolder("PreXt", "PreXt");
65 m_hlist->Add(m_fdPreXt);
66
67 m_fdNhit = new TFolder("XtNhit", "XtNhit");
68 m_hlist->Add(m_fdNhit);
69
70 m_haxis = new TH2F("axis", "", 50, 0, 300, 50, 0, 9);
71 m_haxis -> SetStats(0);
72 m_fdPreXt -> Add(m_haxis);
73
74 char hname[200];
75 for(lay=0; lay<MdcCalNLayer; lay++){
76 sprintf(hname, "trec%02d", lay);
77 m_htrec[lay] = new TH1F(hname, "", 310, -20, 600);
78 m_fdPreXt -> Add(m_htrec[lay]);
79 }
80
81 m_nhitTot = new TH1F("nhitTot", "", 43, -0.5, 42.5);
82 m_fdNhit -> Add(m_nhitTot);
83
84 for(lay=0; lay<MdcCalNLayer; lay++){
85 sprintf(hname, "nhitBin%02d", lay);
86 m_nhitBin[lay] = new TH1F(hname, "", 40, 5.0, 405.0);
87 m_fdNhit -> Add(m_nhitBin[lay]);
88 }
89
90 /* integrated time Spectrum for determination X-T */
91 int bin;
92 m_nXtBin = 40;
93 double twid = 10.0; // ns
94 for(bin=0; bin<m_nXtBin; bin++) m_tbin[bin] = (double)(bin+1) * twid;
95
96 for(lay=0; lay<MdcCalNLayer; lay++){
97 m_nTot[lay] = 0;
98 for(bin=0; bin<m_nXtBin; bin++){
99 m_nEntries[lay][bin] = 0;
100 }
101 }
102}
103
105 IMessageSvc* msgSvc;
106 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
107 MsgStream log(msgSvc, "PreXtMdcCalib");
108 log << MSG::DEBUG << "PreXtMdcCalib::fillHist()" << endreq;
109
110 int bin;
111 int lay;
112 int digiId;
113 unsigned fgOverFlow;
114 double tdc;
115 double adc;
116 double traw;
117 double trec;
118 Identifier id;
119
120 IDataProviderSvc* eventSvc = NULL;
121 Gaudi::svcLocator()->service("EventDataSvc", eventSvc);
122
123 double xtpar[8];
124 if(! m_fgIniTm){
125 for(lay=0; lay<MdcCalNLayer; lay++){
126 m_t0[lay] = m_mdcFunSvc -> getT0(lay, 0);
127 m_mdcFunSvc->getXtpar(lay, 0, 0, xtpar);
128 m_tm[lay] = xtpar[6];
129 }
130 m_fgIniTm = true;
131 }
132
133 // get EsTimeCol
134 double tes = -9999.0;
135 int esTimeflag = -1;
136 SmartDataPtr<RecEsTimeCol> aevtimeCol(eventSvc,"/Event/Recon/RecEsTimeCol");
137 if( (!aevtimeCol) || (aevtimeCol->size()==0) ){
138 tes = -9999.0;
139 esTimeflag = -1;
140 }else{
141 RecEsTimeCol::iterator iter_evt = aevtimeCol->begin();
142 for(; iter_evt!=aevtimeCol->end(); iter_evt++){
143 tes = (*iter_evt)->getTest();
144 esTimeflag = (*iter_evt)->getStat();
145 }
146 }
147 bool flagTes = false;
148 for(int iEs=0; iEs<m_param.nEsFlag; iEs++){
149 if(esTimeflag == m_param.esFlag[iEs]){
150 flagTes = true;
151 break;
152 }
153 }
154 if( (!flagTes) || (tes < m_param.tesMin) || (tes > m_param.tesMax) ) return -1;
155
156 // retrieve Mdc digi
157 SmartDataPtr<MdcDigiCol> mdcDigiCol(eventSvc,"/Event/Digi/MdcDigiCol");
158 if (!mdcDigiCol) {
159 log << MSG::FATAL << "Could not find event" << endreq;
160 }
161
162 MdcDigiCol::iterator iter = mdcDigiCol->begin();
163 digiId = 0;
164 for(; iter != mdcDigiCol->end(); iter++, digiId++) {
165 MdcDigi *aDigi = (*iter);
166 id = (aDigi)->identify();
167
168 lay = MdcID::layer(id);
169
170 tdc = (aDigi) -> getTimeChannel();
171 adc = (aDigi) -> getChargeChannel();
172 fgOverFlow = (aDigi) -> getOverflow();
173 traw = tdc * MdcCalTdcCnv;
174 trec = traw - tes - m_t0[lay];
175
176// if ( !((fgOverFlow == 0)||(fgOverFlow ==12)) ||
177// (aDigi->getTimeChannel() == 0x7FFFFFFF) ||
178// (aDigi->getChargeChannel() == 0x7FFFFFFF) ) continue;
179 if ( ((fgOverFlow & 1) !=0 ) || ((fgOverFlow & 12) != 0) ||
180 (aDigi->getTimeChannel() == 0x7FFFFFFF) ||
181 (aDigi->getChargeChannel() == 0x7FFFFFFF) ) continue;
182
183
184 /* integrated time Spectrum for determination X-T */
185 if(trec > 0){
186 if(trec < m_tm[lay]){
187 m_nTot[lay]++;
188 m_htrec[lay]->Fill(trec);
189 m_nhitTot->Fill(lay);
190 }
191 for(bin=0; bin<m_nXtBin; bin++){
192 if(trec < m_tbin[bin]){
193 m_nEntries[lay][bin]++;
194 m_nhitBin[lay]->Fill(m_tbin[bin]);
195 }
196 }
197 }
198 }
199 return 1;
200}
201
203 IMessageSvc* msgSvc;
204 Gaudi::svcLocator() -> service("MessageSvc", msgSvc);
205 MsgStream log(msgSvc, "PreXtMdcCalib");
206 log << MSG::DEBUG << "PreXtMdcCalib::updateConst()" << endreq;
207
208 int lay;
209 int bin;
210 int iLR;
211 int iEntr;
212 int ord;
213
214 double layRad;
215 double ncel;
216 double pi = 3.141592653;
217
218 double dm;
219 double dist[40];
220 double xtpar[6];
221 char hname[200];
222
223 TF1* funXt = new TF1("funXt", xtfun, 0, 300, 6);
224 funXt -> FixParameter(0, 0.0);
225 funXt -> SetParameter(1, 0.03);
226 funXt -> SetParameter(2, 0.0);
227 funXt -> SetParameter(3, 0.0);
228 funXt -> SetParameter(4, 0.0);
229 funXt -> SetParameter(5, 0.0);
230
231 ofstream fxtlog("preXtpar.dat");
232 for(lay=0; lay<MdcCalNLayer; lay++){
233 sprintf(hname, "grPreXt%02d", lay);
234 m_grXt[lay] = new TGraph();
235 m_grXt[lay] -> SetName(hname);
236 m_grXt[lay] -> SetMarkerStyle(20);
237 m_fdPreXt -> Add(m_grXt[lay]);
238
239 layRad = m_mdcGeomSvc -> Layer(lay) -> Radius();
240 ncel = m_mdcGeomSvc -> Layer(lay) -> NCell();
241 dm = pi * layRad / ncel;
242
243 fxtlog << "layer " << lay << endl;
244 for(bin=0; bin<m_nXtBin; bin++){
245 dist[bin] = dm * m_nEntries[lay][bin] / m_nTot[lay];
246 m_grXt[lay] -> SetPoint(bin, m_tbin[bin], dist[bin]);
247 fxtlog << setw(4) << bin << setw(15) << m_tbin[bin]
248 << setw(15) << dist[bin] << setw(15) << dm
249 << setw(10) << m_nEntries[lay][bin]
250 << setw(10) << m_nTot[lay] << endl;
251
252 if(m_tbin[bin] >= m_tm[lay]) break;
253 }
254
255 if(1 == m_param.fgCalib[lay]){
256 m_grXt[lay] -> Fit(funXt, "Q", "", 0.0, m_tm[lay]);
257 for(ord=0; ord<6; ord++){
258 xtpar[ord] = funXt -> GetParameter(ord);
259 }
260
261 for(iEntr=0; iEntr<MdcCalNENTRXT; iEntr++){
262 for(iLR=0; iLR<MdcCalLR; iLR++){
263 for(ord=0; ord<6; ord++){
264 calconst -> resetXtpar(lay, iEntr, iLR, ord, xtpar[ord]);
265 }
266 }
267 }
268 } else{
269 for(ord=0; ord<6; ord++) xtpar[ord] = calconst->getXtpar(lay, 0, 0, ord);
270 }
271
272 for(ord=0; ord<6; ord++) fxtlog << setw(14) << xtpar[ord];
273 fxtlog << setw(10) << m_tm[lay] << " 0" << endl;
274 } // end of layer loop
275 fxtlog.close();
276 cout << "preXt.dat was written." << endl;
277
278 delete funXt;
279 return 1;
280}
281
282Double_t PreXtMdcCalib::xtfun(Double_t *x, Double_t *par){
283 Double_t val = 0.0;
284 for(Int_t ord=0; ord<6; ord++){
285 val += par[ord] * pow(x[0], ord);
286 }
287 return val;
288}
289
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)
void Fit()
Definition: Eangle1D/Fit.cxx:3
*******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
const int MdcCalNLayer
Definition: MdcCalParams.h:6
const double MdcCalTdcCnv
Definition: MdcCalParams.h:24
const int MdcCalNENTRXT
Definition: MdcCalParams.h:12
const int MdcCalLR
Definition: MdcCalParams.h:10
IMessageSvc * msgSvc()
#define NULL
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
virtual void getXtpar(int layid, int entr, int lr, double par[]) const =0
int fgCalib[MdcCalNLayer]
Definition: MdcCalParams.h:75
double tesMax
Definition: MdcCalParams.h:44
int esFlag[50]
Definition: MdcCalParams.h:41
double tesMin
Definition: MdcCalParams.h:43
double getXtpar(int lay, int entr, int lr, int order)
static int layer(const Identifier &id)
Values of different levels (failure returns 0)
Definition: MdcID.cxx:49
int updateConst(MdcCalibConst *calconst)
int fillHist(MdcCalEvent *event)
void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc, IMdcUtilitySvc *mdcUtilitySvc)
unsigned int getChargeChannel() const
Definition: RawData.cxx:45
unsigned int getTimeChannel() const
Definition: RawData.cxx:40
Definition: Event.h:21
const float pi
Definition: vector3.h:133