BOSS 7.0.6
BESIII Offline Software System
Loading...
Searching...
No Matches
RootEstTofCalibDataCnv.cxx
Go to the documentation of this file.
1// $Header: /bes/bes/BossCvs/Calibration/CalibSvc/CalibROOTCnv/src/cnv/RootEstTofCalibDataCnv.cxx,v 1.9 2019/09/20 07:01:13 sunss Exp $
2#include "GaudiKernel/MsgStream.h"
9
10#include "TFile.h"
11#include "TTree.h"
12#include "TDirectory.h"
13#include "TObject.h"
14
15#include "GaudiKernel/CnvFactory.h"
16#include "GaudiKernel/IOpaqueAddress.h"
17#include "GaudiKernel/DataObject.h"
18#include "GaudiKernel/IAddressCreator.h"
19#include "GaudiKernel/IDataProviderSvc.h"
20#include "GaudiKernel/IConversionSvc.h"
21#include "GaudiKernel/GenericAddress.h"
22
23#include "CalibDataSvc/ICalibRootSvc.h" //maybe
25
26// Temporary. Hope to find a better way to do this
28using namespace CalibData;
29//static CnvFactory<RootEstTofCalibDataCnv> TofCalib_factory;
30//const ICnvFactory& RootEstTofCalibDataCnvFactory = TofCalib_factory;
31
32const unsigned int nBarPar = 10;
33const unsigned int nBarParOff = 20;
34const unsigned int nEndPar = 8;
35const unsigned int nEtfPar = 20;
36const unsigned int nBarOffset = 2;
37
40 }
41
42
45}
46
49}
50
51StatusCode RootEstTofCalibDataCnv::i_createObj(const std::string& fname,
52 DataObject*& refpObject) {
53
54 MsgStream log(msgSvc(), "RootEstTofCalibDataCnv");
55 log<<MSG::DEBUG<<"SetProperty"<<endreq;
56
57 // open the file
58 StatusCode sc = openRead(fname);
59 if(!sc)
60 { log<<MSG::ERROR<<"unable to open files"<<endreq;
61 }
62
68
69 std::vector<CalibData::bTofCalibBase> tmpbTof;
70 std::vector<CalibData::eTofCalibBase> tmpeTof;
71 std::vector<CalibData::etfCalibBase> tmpetf;
72 std::vector<CalibData::bTofCommonCalibBase> tmpbTofCommon;
73 std::vector<CalibData::tofCalibInfoBase> tofinfoCol;
74
75 // Read in the object
76 int cnt;
77 // read btoftree ------------------------------------------------------------
78 double cnvBarPar1[nBarPar];
79 double cnvBarPar2[nBarPar];
80 double cnvBarParOff1_bunch0[nBarParOff];
81 double cnvBarParOff2_bunch0[nBarParOff];
82 double cnvBarParOff1_bunch1[nBarParOff];
83 double cnvBarParOff2_bunch1[nBarParOff];
84 double cnvBarParOff1_bunch2[nBarParOff];
85 double cnvBarParOff2_bunch2[nBarParOff];
86 double cnvBarParOff1_bunch3[nBarParOff];
87 double cnvBarParOff2_bunch3[nBarParOff];
88
89 TTree *btoftree = (TTree*)m_inFile -> Get("BarTofPar");
90
91 char brname[10];
92 for( unsigned int i=0; i<nBarPar; i++ ) {
93 sprintf( brname, "P%i", i );
94 btoftree -> SetBranchAddress( brname, &cnvBarPar1[i] );
95 }
96 for( unsigned int i=0; i<nBarPar; i++ ) {
97 sprintf( brname, "P%i", i+nBarPar );
98 btoftree -> SetBranchAddress( brname, &cnvBarPar2[i] );
99 }
100 for( unsigned int i=0; i<nBarParOff; i++ ) {
101 sprintf( brname, "Bunch0_Poff%i", i );
102 btoftree -> SetBranchAddress( brname, &cnvBarParOff1_bunch0[i] );
103 }
104 for( unsigned int i=0; i<nBarParOff; i++ ) {
105 sprintf( brname, "Bunch0_Poff%i", i+nBarParOff );
106 btoftree -> SetBranchAddress( brname, &cnvBarParOff2_bunch0[i] );
107 }
108 for( unsigned int i=0; i<nBarParOff; i++ ) {
109 sprintf( brname, "Bunch1_Poff%i", i );
110 btoftree -> SetBranchAddress( brname, &cnvBarParOff1_bunch1[i] );
111 }
112 for( unsigned int i=0; i<nBarParOff; i++ ) {
113 sprintf( brname, "Bunch1_Poff%i", i+nBarParOff );
114 btoftree -> SetBranchAddress( brname, &cnvBarParOff2_bunch1[i] );
115 }
116 for( unsigned int i=0; i<nBarParOff; i++ ) {
117 sprintf( brname, "Bunch2_Poff%i", i );
118 btoftree -> SetBranchAddress( brname, &cnvBarParOff1_bunch2[i] );
119 }
120 for( unsigned int i=0; i<nBarParOff; i++ ) {
121 sprintf( brname, "Bunch2_Poff%i", i+nBarParOff );
122 btoftree -> SetBranchAddress( brname, &cnvBarParOff2_bunch2[i] );
123 }
124 for( unsigned int i=0; i<nBarParOff; i++ ) {
125 sprintf( brname, "Bunch3_Poff%i", i );
126 btoftree -> SetBranchAddress( brname, &cnvBarParOff1_bunch3[i] );
127 }
128 for( unsigned int i=0; i<nBarParOff; i++ ) {
129 sprintf( brname, "Bunch3_Poff%i", i+nBarParOff );
130 btoftree -> SetBranchAddress( brname, &cnvBarParOff2_bunch3[i] );
131 }
132
133 for(cnt=0; cnt< btoftree->GetEntries(); cnt++){
134 btoftree -> GetEntry(cnt);
135 bTof.setP1( cnvBarPar1 );
136 bTof.setP2( cnvBarPar2 );
137 bTof.setPoff1_bunch0( cnvBarParOff1_bunch0 );
138 bTof.setPoff2_bunch0( cnvBarParOff2_bunch0 );
139 bTof.setPoff1_bunch1( cnvBarParOff1_bunch1 );
140 bTof.setPoff2_bunch1( cnvBarParOff2_bunch1 );
141 bTof.setPoff1_bunch2( cnvBarParOff1_bunch2 );
142 bTof.setPoff2_bunch2( cnvBarParOff2_bunch2 );
143 bTof.setPoff1_bunch3( cnvBarParOff1_bunch3 );
144 bTof.setPoff2_bunch3( cnvBarParOff2_bunch3 );
145 tmpbTof.push_back( bTof );
146 }
147
148 //read etoftree
149 double cnvEndPar[nEndPar];
150
151 TTree *etoftree = (TTree*)m_inFile -> Get("EndTofPar");
152
153 char ecname[10];
154 for( unsigned int i=0; i<nEndPar; i++ ) {
155 sprintf( ecname, "P%i", i );
156 etoftree -> SetBranchAddress( ecname, &cnvEndPar[i] );
157 }
158
159 for(cnt=0; cnt<etoftree->GetEntries(); cnt++){
160 etoftree->GetEntry(cnt);
161 eTof.setP( cnvEndPar );
162 tmpeTof.push_back(eTof);
163 }
164
165 //read etftree
166 double cnvEtfPar[nEtfPar];
167 double cnvEtfPar1[nEtfPar];
168 double cnvEtfPar2[nEtfPar];
169
170 if( NULL!=m_inFile->Get("EtfTofPar") ) {
171 TTree *etftree = (TTree*)m_inFile -> Get("EtfTofPar");
172
173 char etfname[10];
174 for( unsigned int i=0; i<nEtfPar; i++ ) {
175 sprintf( etfname, "P%i", i );
176 etftree -> SetBranchAddress( etfname, &cnvEtfPar[i] );
177 }
178 for( unsigned int i=0; i<nEtfPar; i++ ) {
179 sprintf( etfname, "P%i", i+nEtfPar );
180 etftree -> SetBranchAddress( etfname, &cnvEtfPar1[i] );
181 }
182 for( unsigned int i=0; i<nEtfPar; i++ ) {
183 sprintf( etfname, "P%i", i+2*nEtfPar );
184 etftree -> SetBranchAddress( etfname, &cnvEtfPar2[i] );
185 }
186
187 for(cnt=0; cnt<etftree->GetEntries(); cnt++){
188 etftree -> GetEntry(cnt);
189 etf.setP( cnvEtfPar );
190 etf.setP1( cnvEtfPar1 );
191 etf.setP2( cnvEtfPar2 );
192 tmpetf.push_back(etf);
193 }
194 }
195
196 //read bTofCommonCalibBase
197 double cnvBarOffset[nBarOffset];
198 TTree *btofcommontree = (TTree*)m_inFile -> Get("BarTofParCommon");
199 for( unsigned int i=0; i<nBarOffset; i++ ) {
200 sprintf( brname, "t0offset%i", i );
201 btofcommontree-> SetBranchAddress( brname, &cnvBarOffset[i]);
202 }
203
204 int entries = btofcommontree->GetEntries();
205 for(cnt=0;cnt<entries;cnt++){
206 btofcommontree->GetEntry(cnt);
207 bTofCommon.setOffset( cnvBarOffset );
208 tmpbTofCommon.push_back(bTofCommon);
209 }
210
211 int m_runFrom, m_runTo, m_eventFrom, m_eventTo;
212
213 TTree *CalibInfo = (TTree*)m_inFile -> Get("CalibInfo");
214
215 if( CalibInfo->GetBranchStatus("runFrom") ) {
216 CalibInfo->SetBranchAddress("runFrom", &m_runFrom );
217 CalibInfo->SetBranchAddress("runTo", &m_runTo );
218 CalibInfo->SetBranchAddress("eventFrom", &m_eventFrom );
219 CalibInfo->SetBranchAddress("eventTo", &m_eventTo );
220 }
221 else {
222 m_runFrom = -1;
223 m_runTo = -1;
224 m_eventFrom = -1;
225 m_eventTo = -1;
226 }
227
228 entries= CalibInfo->GetEntries();
229 for(cnt=0;cnt<entries;cnt++){
230 CalibInfo->GetEntry(cnt);
231 tofinfo.setRunFrom(m_runFrom);
232 tofinfo.setRunTo(m_runTo);
233 tofinfo.setEventFrom(m_eventFrom);
234 tofinfo.setEventTo(m_eventTo);
235 tofinfoCol.push_back(tofinfo);
236 }
237
238 CalibData::TofCalibData *tmpObject = new CalibData::TofCalibData(&tmpbTof,&tmpbTofCommon,&tmpeTof,&tmpetf,&tofinfoCol);
239
240 refpObject=tmpObject;
241
242 return StatusCode::SUCCESS;
243}
244
245StatusCode RootEstTofCalibDataCnv::createRoot(const std::string& fname,
246 CalibData::CalibBase1* pTDSObj) {
247 MsgStream log(msgSvc(), "RootEstTofCalibDataCnv");
248
249 // Open the file, create the branch
250 StatusCode sc = openWrite(fname);
251 if(!sc)
252 { log<<MSG::ERROR<<"unable to open files"<<endreq;
253 }
254 // write the Data in the TCDS to RootFile
255 CalibData::TofCalibData* btof = dynamic_cast<CalibData::TofCalibData*>(pTDSObj);
256
257 // write btoftree----------------------------------------------------------------
258 double cnvBarPar1[nBarPar];
259 double cnvBarPar2[nBarPar];
260
261 char brname[8], ibrname[8];
262 TTree *btoftree = new TTree("BarTofPar", "BarTofPar");
263 for( unsigned int i=0; i<nBarPar; i++ ) {
264 sprintf( brname, "P%i", i );
265 sprintf( ibrname, "P%i/D", i );
266 btoftree -> Branch( brname, &cnvBarPar1[i], ibrname );
267 }
268
269 for( int i=0; i<176; i++ ) {
270 for( int j=0;j<static_cast<int>(nBarPar);j++){
271 cnvBarPar1[j] = btof->getBTofPleft(i,j);
272 cnvBarPar2[j] = btof->getBTofPright(i,j);
273 }
274 btoftree -> Fill();
275 }
276
277 //write etoftree----------------------------------------------------------------
278 double cnvEndPar[nEndPar];
279
280 char ecname[8], iecname[8];
281 TTree *etoftree = new TTree("EndTofPar", "EndTofPar");
282 for( unsigned int i=0; i<nEndPar; i++ ) {
283 sprintf( ecname, "P%i", i );
284 sprintf( iecname, "P%i/D", i );
285 etoftree -> Branch( ecname, &cnvEndPar[i], iecname );
286 }
287
288 for(int i=0; i<96; i++){
289 for(int j=0;j<static_cast<int>(nEndPar);j++){
290 cnvEndPar[j] = btof->getETofP(i,j);
291 }
292 etoftree->Fill();
293 }
294
295 //write etftree----------------------------------------------------------------
296 double cnvEtfPar[nEtfPar];
297 double cnvEtfPar1[nEtfPar];
298 double cnvEtfPar2[nEtfPar];
299
300 char etfname[8], ietfname[8];
301 TTree *etftree = new TTree("EtfTofPar", "EtfTofPar");
302 for( unsigned int i=0; i<nEtfPar; i++ ) {
303 sprintf( etfname, "P%i", i );
304 sprintf( ietfname, "P%i/D", i );
305 etftree -> Branch( etfname, &cnvEtfPar[i], ietfname );
306 }
307 for( unsigned int i=0; i<nEtfPar; i++ ) {
308 sprintf( etfname, "P%i", i+nEtfPar );
309 sprintf( ietfname, "P%i/D", i+nEtfPar );
310 etftree -> Branch( etfname, &cnvEtfPar1[i], ietfname );
311 }
312 for( unsigned int i=0; i<nEtfPar; i++ ) {
313 sprintf( etfname, "P%i", i+2*nEtfPar );
314 sprintf( ietfname, "P%i/D", i+2*nEtfPar );
315 etftree -> Branch( etfname, &cnvEtfPar2[i], ietfname );
316 }
317
318 for( int i=0; i<72; i++ ) {
319 for( int k=0; k<12; k++ ) {
320 for(int j=0;j<static_cast<int>(nEtfPar);j++){
321 cnvEtfPar[j] = btof->getEtfPcombine(i,k,j);
322 cnvEtfPar1[j] = btof->getEtfPleft(i,k,j);
323 cnvEtfPar2[j] = btof->getEtfPright(i,k,j);
324 }
325 etftree -> Fill();
326 }
327 }
328
329 // write all the trees
330 btoftree -> Write();
331 etoftree -> Write();
332 etftree -> Write();
333 delete btoftree;
334 delete etoftree;
335 delete etftree;
336 closeWrite();
337 log<<MSG::INFO<<"successfully create RootFile"<<endreq;
338
339 return sc;
340}
const CLID CLID_Calib_EstTofCal
Definition: CalibModel.h:53
curve Branch("CurveSize",&CurveSize,"CurveSize/I")
curve Fill()
curve Write()
const unsigned int nBarOffset
const unsigned int nBarParOff
const unsigned int nEndPar
const unsigned int nEtfPar
const unsigned int nBarPar
IMessageSvc * msgSvc()
double getETofP(int index, int pardex)
double getEtfPleft(int index, int strip, int pardex)
double getEtfPcombine(int index, int strip, int pardex)
double getEtfPright(int index, int strip, int pardex)
double getBTofPleft(int index, int pardex)
double getBTofPright(int index, int pardex)
void setP2(const double *TofP2)
void setPoff1_bunch3(const double *TofPoff1_bunch3)
void setPoff2_bunch1(const double *TofPoff2_bunch1)
void setPoff1_bunch0(const double *TofPoff1_bunch0)
void setPoff2_bunch2(const double *TofPoff2_bunch2)
void setPoff1_bunch1(const double *TofPoff1_bunch1)
void setPoff1_bunch2(const double *TofPoff1_bunch2)
void setPoff2_bunch3(const double *TofPoff2_bunch3)
void setPoff2_bunch0(const double *TofPoff2_bunch0)
void setP1(const double *TofP1)
void setOffset(const double *offset)
void setP(const double *TofP)
void setP(const double *etfP)
void setP1(const double *etfP1)
void setP2(const double *etfP2)
void setRunFrom(const int runFrom)
void setEventTo(const int eventTo)
void setRunTo(const int runTo)
void setEventFrom(const int eventFrom)
StatusCode openRead(const std::string &fname)
StatusCode closeWrite()
virtual StatusCode openWrite(const std::string &fname)
const CLID & objType() const
RootEstTofCalibDataCnv(ISvcLocator *svc)
virtual StatusCode createRoot(const std::string &fname, CalibData::CalibBase1 *pTDSObj)
static const CLID & classID()
virtual StatusCode i_createObj(const std::string &fname, DataObject *&refpObject)
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)
curve GetEntry(0)
curve SetBranchAddress("CurveSize",&CurveSize)