CGEM BOSS 6.6.5.h
BESIII Offline Software System
Loading...
Searching...
No Matches
RootDedxCalibDataCnv Class Reference

#include <RootDedxCalibDataCnv.h>

+ Inheritance diagram for RootDedxCalibDataCnv:

Public Member Functions

const CLID & objType () const
 
 RootDedxCalibDataCnv (ISvcLocator *svc)
 
virtual ~RootDedxCalibDataCnv ()
 
virtual StatusCode createRoot (const std::string &fname, CalibData::CalibBase1 *pTDSObj)
 
virtual long repSvcType () const
 
- Public Member Functions inherited from RootCalBaseCnv
virtual ~RootCalBaseCnv ()
 
virtual StatusCode initialize ()
 
virtual StatusCode finalize ()
 
virtual StatusCode createObj (IOpaqueAddress *addr, DataObject *&refpObject)
 
ICalibRootSvcgetCalibRootSvc ()
 
 RootCalBaseCnv (ISvcLocator *svc, const CLID &clid)
 
virtual StatusCode readRootObj (const std::string &treename, const std::string &branch, TObject *&pCalib, unsigned index=0)
 
virtual StatusCode readRootObj (TTree *tree, const std::string &branch, TObject *&pCalib, unsigned index=0)
 
- Public Member Functions inherited from Converter< Ty1, Ty2 >
destinationoperator (const source &) const
 

Static Public Member Functions

static const CLID & classID ()
 
- Static Public Member Functions inherited from RootCalBaseCnv
static const unsigned char storageType ()
 

Protected Member Functions

virtual StatusCode i_createObj (const std::string &fname, DataObject *&refpObject)
 
- Protected Member Functions inherited from RootCalBaseCnv
virtual StatusCode internalCreateObj (const std::string &fname, DataObject *&refpObject, IOpaqueAddress *address)
 
virtual StatusCode i_processObj (DataObject *pObject, IOpaqueAddress *address)
 In case there is additional work to do on the created object.
 
virtual StatusCode fillRoot (CalibData::CalibBase *pTDSObj, TObject *pRootObj)
 
virtual StatusCode openWrite (const std::string &fname)
 
StatusCode closeWrite ()
 
StatusCode openRead (const std::string &fname)
 
StatusCode closeRead ()
 
void setBaseInfo (CalibData::CalibBase1 *pObj)
 Another utility for derived classes to use.
 
- Protected Member Functions inherited from Converter< Ty1, Ty2 >
virtual destinationconvert (const source &) const =0
 

Friends

class CnvFactory< RootDedxCalibDataCnv >
 

Additional Inherited Members

- Public Types inherited from Converter< Ty1, Ty2 >
typedef Ty1 source
 
typedef Ty2 destination
 
- Protected Attributes inherited from RootCalBaseCnv
ICalibRootSvcm_rootSvc
 
ICalibMetaCnvSvcm_metaSvc
 
IInstrumentNamem_instrSvc
 
int m_serNo
 
ITime * m_vstart
 
ITime * m_vend
 
int m_runfrm
 
int m_runto
 
TFile * m_outFile
 
TTree * m_ttree
 
TFile * m_inFile
 
TDirectory * m_saveDir
 

Detailed Description

Definition at line 20 of file RootDedxCalibDataCnv.h.

Constructor & Destructor Documentation

◆ RootDedxCalibDataCnv()

RootDedxCalibDataCnv::RootDedxCalibDataCnv ( ISvcLocator * svc)

Definition at line 35 of file RootDedxCalibDataCnv.cxx.

35 :
37}
const CLID CLID_Calib_DedxCal
Definition CalibModel.h:45
RootCalBaseCnv(ISvcLocator *svc, const CLID &clid)

◆ ~RootDedxCalibDataCnv()

virtual RootDedxCalibDataCnv::~RootDedxCalibDataCnv ( )
inlinevirtual

Definition at line 29 of file RootDedxCalibDataCnv.h.

29{};

Member Function Documentation

◆ classID()

const CLID & RootDedxCalibDataCnv::classID ( )
static

Definition at line 44 of file RootDedxCalibDataCnv.cxx.

44 {
45 return CLID_Calib_DedxCal;
46}

◆ createRoot()

StatusCode RootDedxCalibDataCnv::createRoot ( const std::string & fname,
CalibData::CalibBase1 * pTDSObj )
virtual

Create ROOT file corresponding to TDS object input.
Default implementation is to return an error. Must be separately implemented for each calibration type.

Parameters
fnameFilename for output file
pTDSObjPointer to tds object to be converted

Reimplemented from RootCalBaseCnv.

Definition at line 290 of file RootDedxCalibDataCnv.cxx.

291 {
292
293 MsgStream log(msgSvc(), "RootDedxCalibDataCnv");
294
295 // Open the file, create the branch
296 StatusCode sc = openWrite(fname);
297 if(!sc)
298 { log<<MSG::ERROR<<"unable to open files"<<endreq;
299 }
300 // write the Data in the TCDS to RootFile
301 int i;
302 CalibData::DedxCalibData* tmpObject = dynamic_cast<CalibData::DedxCalibData*>(pTDSObj);
303 // write rungcalib ------------------------------------------------------------
304 double rungain;
305 double runmean;
306 double runno;
307 double runresol;
308 TTree *rungcalibtree = new TTree("rungcalib", "rungCalib");
309 rungcalibtree -> Branch("rungain", &rungain, "S/D");
310 rungcalibtree -> Branch("runmean", &runmean, "S/D");
311 rungcalibtree -> Branch("runno", &runno, "S/D");
312 rungcalibtree -> Branch("runresol", &runresol, "S/D");
313 int N = tmpObject->getrunNO();
314 for(i=0; i<N; i++){
315 rungain = tmpObject->getrung(0,i);
316 runmean = tmpObject->getrung(1,i);
317 runno = tmpObject->getrung(2,i);
318 runresol = tmpObject->getrung(3,i);
319 rungcalibtree -> Fill();
320 }
321
322
323
324 // write ddgcalib----------------------------------------------------------------
325 double ddg0;
326 double ddg1;
327 double ddg2;
328 double ddg3;
329 TTree *ddgcalibtree = new TTree("ddgcalib", "ddgCalib");
330 ddgcalibtree -> Branch("ddg0", &ddg0, "S/D");
331 ddgcalibtree -> Branch("ddg1", &ddg1, "S/D");
332 ddgcalibtree -> Branch("ddg2", &ddg2, "S/D");
333 ddgcalibtree -> Branch("ddg3", &ddg3, "S/D");
334 for(i=0; i<43; i++){
335 ddg0 = tmpObject->getddg(0,i);
336 ddg1 = tmpObject->getddg(1,i);
337 ddg2 = tmpObject->getddg(2,i);
338 ddg3 = tmpObject->getddg(3,i);
339 ddgcalibtree -> Fill();
340 }
341
342 // write ggscalib----------------------------------------------------------------
343 double ggs0;
344 double ggs1;
345 double ggs2;
346 double ggs3;
347 TTree *ggscalibtree = new TTree("ggscalib", "ggsCalib");
348 ggscalibtree -> Branch("ggs0", &ggs0, "S/D");
349 ggscalibtree -> Branch("ggs1", &ggs1, "S/D");
350 ggscalibtree -> Branch("ggs2", &ggs2, "S/D");
351 ggscalibtree -> Branch("ggs3", &ggs3, "S/D");
352 for(i=0; i<43; i++){
353 ggs0 = tmpObject->getggs(0,i);
354 ggs1 = tmpObject->getggs(1,i);
355 ggs2 = tmpObject->getggs(2,i);
356 ggs3 = tmpObject->getggs(3,i);
357 ggscalibtree -> Fill();
358 }
359
360// write zdepcalib----------------------------------------------------------------
361 double zdep0;
362 double zdep1;
363 double zdep2;
364 double zdep3;
365 TTree *zdepcalibtree = new TTree("zdepcalib", "zdepCalib");
366 zdepcalibtree -> Branch("zdep0", &zdep0, "S/D");
367 zdepcalibtree -> Branch("zdep1", &zdep1, "S/D");
368 zdepcalibtree -> Branch("zdep2", &zdep2, "S/D");
369 zdepcalibtree -> Branch("zdep3", &zdep3, "S/D");
370 for(i=0; i<43; i++){
371 zdep0 = tmpObject->getzdep(0,i);
372 zdep1 = tmpObject->getzdep(1,i);
373 zdep2 = tmpObject->getzdep(2,i);
374 zdep3 = tmpObject->getzdep(3,i);
375 zdepcalibtree -> Fill();
376 }
377
378// write entacalib----------------------------------------------------------------
379 double enta0;
380 double enta1;
381 double enta2;
382 double enta3;
383 TTree *entacalibtree = new TTree("entacalib", "entaCalib");
384 entacalibtree -> Branch("enta0", &enta0, "S/D");
385 entacalibtree -> Branch("enta1", &enta1, "S/D");
386 entacalibtree -> Branch("enta2", &enta2, "S/D");
387 entacalibtree -> Branch("enta3", &enta3, "S/D");
388 for(i=0; i<43; i++){
389 enta0 = tmpObject->getenta(0,i);
390 enta1 = tmpObject->getenta(1,i);
391 enta2 = tmpObject->getenta(2,i);
392 enta3 = tmpObject->getenta(3,i);
393 entacalibtree -> Fill();
394 }
395
396
397
398// write gain---------------------------------------------------------------------
399 double gain;
400 TTree *gaintree = new TTree("gaincalib", "gainCalib");
401 gaintree -> Branch("gain", &gain, "S/D");
402 gain = tmpObject->getgain();
403 gaintree -> Fill();
404
405// write resol---------------------------------------------------------------------
406 double resol;
407 TTree *resoltree = new TTree("resolcalib", "resolCalib");
408 resoltree -> Branch("resol", &resol, "S/D");
409 resol = tmpObject->getresol();
410 resoltree -> Fill();
411
412
413// write wiregcalib----------------------------------------------------------------
414 double wireg;
415 TTree *wiregcalibtree = new TTree("wiregcalib", "wiregCalib");
416 wiregcalibtree -> Branch("wireg", &wireg, "S/D");
417 for(i=0; i<6796; i++){
418 wireg = tmpObject->getwireg(i);
419 wiregcalibtree -> Fill();
420 }
421
422// write layergcalib----------------------------------------------------------------
423 double layerg;
424 TTree *layergcalibtree = new TTree("layergcalib", "layergCalib");
425 layergcalibtree -> Branch("layerg", &layerg, "S/D");
426 for(i=0; i<43; i++){
427 layerg = tmpObject->getlayerg(i);
428 layergcalibtree -> Fill();
429 }
430
431// write all the trees
432 rungcalibtree -> Write();
433 ddgcalibtree -> Write();
434 ggscalibtree -> Write();
435 zdepcalibtree -> Write();
436 entacalibtree -> Write();
437 gaintree -> Write();
438 resoltree -> Write();
439 wiregcalibtree -> Write();
440 layergcalibtree -> Write();
441
442 delete rungcalibtree;
443 delete ddgcalibtree;
444 delete ggscalibtree;
445 delete zdepcalibtree;
446 delete entacalibtree;
447 delete gaintree;
448 delete resoltree;
449 delete wiregcalibtree;
450 delete layergcalibtree;
451 closeWrite();
452 log<<MSG::INFO<<"successfully create RootFile"<<endreq;
453 return sc;
454
455}
IMessageSvc * msgSvc()
double getenta(int i, int j) const
double getzdep(int i, int j) const
double getrung(int i, int j) const
double getlayerg(int i) const
double getwireg(int i) const
double getddg(int i, int j) const
double getggs(int i, int j) const
StatusCode closeWrite()
virtual StatusCode openWrite(const std::string &fname)

◆ i_createObj()

StatusCode RootDedxCalibDataCnv::i_createObj ( const std::string & fname,
DataObject *& refpObject )
protectedvirtual

Create the transient representation of an object, given an opaque address. This and the following update method comprise the core functionality of calibration converters. Convenience routine used by most CAL calibration types, which have a <dimension> element describing how the remainder of the data is laid out. Read from TDS; store information internally in protected members.
Given a pointer to a TDS object which can be cast to "our" type, fill in corresponding information in the corresponding root class

Parameters
pTDSObjPointer to tds object to be converted
pRootObjPointer to destination root object Read in object from specified branch. Don't need tree name; it's always Calib

Reimplemented from RootCalBaseCnv.

Definition at line 48 of file RootDedxCalibDataCnv.cxx.

49 {
50
51 MsgStream log(msgSvc(), "RootDedxCalibDataCnv");
52 log<<MSG::DEBUG<<"SetProperty"<<endreq;
53 StatusCode sc = openRead(fname);
54 if(!sc)
55 { log<<MSG::ERROR<<"unable to open files"<<endreq;
56 }
57
59 // Read in our object
60 int i;
61 // read rungcalib ------------------------------------------------------------
62 double rungain;
63 double runmean;
64 int runno;
65 double runresol;
66 TTree *rungtree = (TTree*)m_inFile -> Get("runcalib");
67 rungtree -> SetBranchAddress("rungain", &rungain);
68 rungtree -> SetBranchAddress("runmean", &runmean);
69 rungtree -> SetBranchAddress("runno", &runno);
70 rungtree -> SetBranchAddress("runresol", &runresol);
71 int N = rungtree -> GetEntries();
72 tmpObject -> setrunNO(N);
73 //std::cout<<"N = "<<N<<std::endl;
74 for(i=0; i<N; i++){
75 rungtree -> GetEntry(i);
76 //std::cout<<rungain<<" "<<runmean<<" "<<runno<<" "<<runresol<<std::endl;
77 tmpObject -> setrung(rungain,0,i);
78 tmpObject -> setrung(runmean,1,i);
79 tmpObject -> setrung(runno,2,i);
80 tmpObject -> setrung(runresol,3,i);
81 }
82
83
84
85 // read ddgcalib ------------------------------------------------------------
86 double ddg0[43];
87 double ddg1[43];
88 double ddg2[43];
89 double ddg3[43];
90 double id_doca[1600];
91 double iner_chi[1600];
92 double iner_gain[1600];
93 double iner_hits[1600];
94 double ip_eangle[1600];
95 double out_chi[1600];
96 double out_gain[1600];
97 double out_hits[1600];
98
99 TTree *ddgtree = (TTree*)m_inFile -> Get("ddgcalib");
100 ddgtree -> SetBranchAddress("ddg0", ddg0);
101 ddgtree -> SetBranchAddress("ddg1", ddg1);
102 ddgtree -> SetBranchAddress("ddg2", ddg2);
103 ddgtree -> SetBranchAddress("ddg3", ddg3);
104 TBranch *bbb = ddgtree->FindBranch("Id_doca");
105 if(bbb){
106 ddgtree -> SetBranchAddress("Id_doca", id_doca);
107 ddgtree -> SetBranchAddress("Iner_chi", iner_chi);
108 ddgtree -> SetBranchAddress("Iner_gain", iner_gain);
109 ddgtree -> SetBranchAddress("Iner_hits", iner_hits);
110 ddgtree -> SetBranchAddress("Ip_eangle", ip_eangle);
111 ddgtree -> SetBranchAddress("Out_chi", out_chi);
112 ddgtree -> SetBranchAddress("Out_gain", out_gain);
113 ddgtree -> SetBranchAddress("Out_hits", out_hits);
114 }
115 ddgtree -> GetEntry(0);
116 for(i=0; i<43; i++){
117 tmpObject -> setddg(ddg0[i],0,i);
118 tmpObject -> setddg(ddg1[i],1,i);
119 tmpObject -> setddg(ddg2[i],2,i);
120 tmpObject -> setddg(ddg3[i],3,i);
121 }
122
123 for(i=0; i<1600; i++){
124 if(!bbb){
125 id_doca[i]=0;
126 iner_chi[i]=0;
127 iner_gain[i]=0;
128 iner_hits[i]=0;
129 ip_eangle[i]=0;
130 out_chi[i]=0;
131 out_gain[i]=0;
132 out_hits[i]=0;
133 }
134 tmpObject -> set_id_doca(id_doca[i],i);
135 tmpObject -> set_iner_chi(iner_chi[i],i);
136 tmpObject -> set_iner_gain(iner_gain[i],i);
137 tmpObject -> set_iner_hits(iner_hits[i],i);
138 tmpObject -> set_ip_eangle(ip_eangle[i],i);
139 tmpObject -> set_out_chi(out_chi[i],i);
140 tmpObject -> set_out_gain(out_gain[i],i);
141 tmpObject -> set_out_hits(out_hits[i],i);
142 }
143
144// read entratree---------------------------------------------
145 double entra0[43];
146 double entra1[43];
147 double entra2[43];
148 double entra3[43];
149 double engle[100];
150 int engle_no;
151 TTree *entratree = (TTree*)m_inFile -> Get("entracalib");
152 entratree -> SetBranchAddress("entra0", entra0);
153 entratree -> SetBranchAddress("entra1", entra1);
154 entratree -> SetBranchAddress("entra2", entra2);
155 entratree -> SetBranchAddress("entra3", entra3);
156 entratree -> SetBranchAddress("1denangle", engle);
157 entratree -> SetBranchAddress("1denangle_entry", &engle_no);
158 entratree -> GetEntry(0);
159 for(i=0; i<43; i++){
160 tmpObject -> setenta(entra0[i],0,i);
161 tmpObject -> setenta(entra1[i],1,i);
162 tmpObject -> setenta(entra2[i],2,i);
163 tmpObject -> setenta(entra3[i],3,i);
164 }
165 tmpObject -> set_enanglesize(engle_no);
166 for(i=0; i<engle_no; i++){
167 tmpObject -> set_enangle(engle[i],i);
168 }
169
170
171
172// read ggscalib ------------------------------------------------------------
173 double ggs0[43];
174 double ggs1[43];
175 double ggs2[43];
176 double ggs3[43];
177 double gcostheta[80];
178 int hadron_entry;
179 double hadron[20];
180 TTree *ggstree = (TTree*)m_inFile -> Get("ggscalib");
181 ggstree -> SetBranchAddress("ggs0", ggs0);
182 ggstree -> SetBranchAddress("ggs1", ggs1);
183 ggstree -> SetBranchAddress("ggs2", ggs2);
184 ggstree -> SetBranchAddress("ggs3", ggs3);
185 ggstree -> SetBranchAddress("hadron", hadron);
186 ggstree -> SetBranchAddress("hadronNo", &hadron_entry);
187 if(bbb){
188 ggstree -> SetBranchAddress("costheta", gcostheta);
189 }
190 ggstree -> GetEntry(0);
191 for(i=0; i<43;i++){
192 tmpObject -> setggs(ggs0[i],0,i);
193 tmpObject -> setggs(ggs1[i],1,i);
194 tmpObject -> setggs(ggs2[i],2,i);
195 tmpObject -> setggs(ggs3[i],3,i);
196 }
197
198 for(i=0; i<80;i++){
199 if(!bbb) gcostheta[i]=0;
200 tmpObject ->set_costheta(gcostheta[i],i);
201 }
202
203 if(hadron_entry>20){
204 log<<MSG::FATAL<<"hadron entry is larger than 20, larger than designed"<<endreq;
205 return StatusCode::FAILURE;
206 }
207
208 tmpObject->set_hadronNo(hadron_entry);
209
210 for(i=0;i<hadron_entry;i++){
211 tmpObject->set_hadron(hadron[i],i);
212 }
213
214
215// read zdepcalib ------------------------------------------------------------
216 double zdep0[43];
217 double zdep1[43];
218 double zdep2[43];
219 double zdep3[43];
220 TTree *zdeptree = (TTree*)m_inFile -> Get("zdepcalib");
221 zdeptree -> SetBranchAddress("zdep0", zdep0);
222 zdeptree -> SetBranchAddress("zdep1", zdep1);
223 zdeptree -> SetBranchAddress("zdep2", zdep2);
224 zdeptree -> SetBranchAddress("zdep3", zdep3);
225 zdeptree -> GetEntry(0);
226
227 for(i=0; i<43;i++){
228 tmpObject -> setzdep(zdep0[i],0,i);
229 tmpObject -> setzdep(zdep1[i],1,i);
230 tmpObject -> setzdep(zdep2[i],2,i);
231 tmpObject -> setzdep(zdep3[i],3,i);
232 }
233
234// read gain ------------------------------------------------------------
235 double gain;
236 double gt0[35],gdedx[35];
237 TTree *gaintree = (TTree*)m_inFile -> Get("gaincalib");
238 gaintree -> SetBranchAddress("gain", &gain);
239 if(bbb){
240 gaintree -> SetBranchAddress("t0", gt0);
241 gaintree -> SetBranchAddress("dedx", gdedx);
242 }
243 gaintree -> GetEntry(0);
244 tmpObject -> setgain(gain);
245 for(i=0; i<35;i++){
246 if(!bbb){
247 gt0[i]=0;
248 gdedx[i]=0;
249 }
250 tmpObject->set_t0(gt0[i],i);
251 tmpObject->set_dedx(gdedx[i],i);
252 }
253
254// read resol----------------------------------------------------------------
255 double resol;
256 TTree *resoltree = (TTree*)m_inFile -> Get("resolcalib");
257 resoltree -> SetBranchAddress("resol", &resol);
258 resoltree -> GetEntry(0);
259 tmpObject -> setresol(resol);
260
261// read wireg----------------------------------------------------------------
262 double wireg[6796];
263 TTree *wiregtree = (TTree*)m_inFile -> Get("wiregcalib");
264 wiregtree -> SetBranchAddress("wireg",wireg);
265 wiregtree -> GetEntry(0);
266 for(i=0;i<6796;i++){
267 tmpObject -> setwireg(wireg[i],i);
268 }
269
270// read layerg----------------------------------------------------------------
271 double layerg[43];
272 TTree *layergtree = (TTree*)m_inFile -> Get("layergcalib");
273 layergtree -> SetBranchAddress("layerg",layerg);
274 layergtree -> GetEntry(0);
275
276 for(i=0;i<43;i++){
277 tmpObject -> setlayerg(layerg[i],i);
278 }
279
280 refpObject=tmpObject;
281 // delete pObj;
282 // std::cout<<"pObj.m_t0[Cell]="<<*(pObj->gett0()+Cell-1);
283 // dynamic_cast<CalibData::Mdct0& >(refpObject);
284 //refpObject = pObj;
285 // std::cout<<"refObject.m_t0[Cell]="<<*(refpObject->gett0()+Cell-1);
286 // refpObject = new Mdct0(&t0_tmp[0]);
287 return StatusCode::SUCCESS;
288}
data SetBranchAddress("time",&time)
data GetEntry(0)
void set_costheta(const double aa, int i)
void set_hadron(const double aa, int i)
void set_t0(const double aa, int i)
void set_dedx(const double aa, int i)
StatusCode openRead(const std::string &fname)

◆ objType()

const CLID & RootDedxCalibDataCnv::objType ( ) const

Definition at line 40 of file RootDedxCalibDataCnv.cxx.

40 {
41 return CLID_Calib_DedxCal;
42}

◆ repSvcType()

virtual long RootDedxCalibDataCnv::repSvcType ( ) const
inlinevirtual

Definition at line 34 of file RootDedxCalibDataCnv.h.

34 {
36 }
unsigned const char CALIBROOT_StorageType

Friends And Related Symbol Documentation

◆ CnvFactory< RootDedxCalibDataCnv >

friend class CnvFactory< RootDedxCalibDataCnv >
friend

Definition at line 1 of file RootDedxCalibDataCnv.h.


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