BOSS 8.0.0
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
destinationoperator (const source &) const
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
virtual destinationconvert (const source &) const =0
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
typedef Ty1 source
typedef Ty2 destination
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 22 of file RootDedxCalibDataCnv.h.

Constructor & Destructor Documentation

◆ RootDedxCalibDataCnv()

RootDedxCalibDataCnv::RootDedxCalibDataCnv ( ISvcLocator * svc)

Definition at line 34 of file RootDedxCalibDataCnv.cxx.

RootCalBaseCnv(ISvcLocator *svc, const CLID &clid)

◆ ~RootDedxCalibDataCnv()

virtual RootDedxCalibDataCnv::~RootDedxCalibDataCnv ( )
inlinevirtual

Definition at line 31 of file RootDedxCalibDataCnv.h.

31{};

Member Function Documentation

◆ classID()

const CLID & RootDedxCalibDataCnv::classID ( )
static

Definition at line 39 of file RootDedxCalibDataCnv.cxx.

39{ return CLID_Calib_DedxCal; }

Referenced by CalibRootCnvSvc::createConverter().

◆ 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 291 of file RootDedxCalibDataCnv.cxx.

292 {
293
294 MsgStream log( msgSvc(), "RootDedxCalibDataCnv" );
295
296 // Open the file, create the branch
297 StatusCode sc = openWrite( fname );
298 if ( !sc ) { log << MSG::ERROR << "unable to open files" << endmsg; }
299 // write the Data in the TCDS to RootFile
300 int i;
301 CalibData::DedxCalibData* tmpObject = dynamic_cast<CalibData::DedxCalibData*>( pTDSObj );
302 // write rungcalib ------------------------------------------------------------
303 double rungain;
304 double runmean;
305 double runno;
306 double runresol;
307 TTree* rungcalibtree = new TTree( "rungcalib", "rungCalib" );
308 rungcalibtree->Branch( "rungain", &rungain, "S/D" );
309 rungcalibtree->Branch( "runmean", &runmean, "S/D" );
310 rungcalibtree->Branch( "runno", &runno, "S/D" );
311 rungcalibtree->Branch( "runresol", &runresol, "S/D" );
312 int N = tmpObject->getrunNO();
313 for ( i = 0; i < N; i++ )
314 {
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 // write ddgcalib----------------------------------------------------------------
323 double ddg0;
324 double ddg1;
325 double ddg2;
326 double ddg3;
327 TTree* ddgcalibtree = new TTree( "ddgcalib", "ddgCalib" );
328 ddgcalibtree->Branch( "ddg0", &ddg0, "S/D" );
329 ddgcalibtree->Branch( "ddg1", &ddg1, "S/D" );
330 ddgcalibtree->Branch( "ddg2", &ddg2, "S/D" );
331 ddgcalibtree->Branch( "ddg3", &ddg3, "S/D" );
332 for ( i = 0; i < 43; i++ )
333 {
334 ddg0 = tmpObject->getddg( 0, i );
335 ddg1 = tmpObject->getddg( 1, i );
336 ddg2 = tmpObject->getddg( 2, i );
337 ddg3 = tmpObject->getddg( 3, i );
338 ddgcalibtree->Fill();
339 }
340
341 // write ggscalib----------------------------------------------------------------
342 double ggs0;
343 double ggs1;
344 double ggs2;
345 double ggs3;
346 TTree* ggscalibtree = new TTree( "ggscalib", "ggsCalib" );
347 ggscalibtree->Branch( "ggs0", &ggs0, "S/D" );
348 ggscalibtree->Branch( "ggs1", &ggs1, "S/D" );
349 ggscalibtree->Branch( "ggs2", &ggs2, "S/D" );
350 ggscalibtree->Branch( "ggs3", &ggs3, "S/D" );
351 for ( i = 0; i < 43; i++ )
352 {
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 {
372 zdep0 = tmpObject->getzdep( 0, i );
373 zdep1 = tmpObject->getzdep( 1, i );
374 zdep2 = tmpObject->getzdep( 2, i );
375 zdep3 = tmpObject->getzdep( 3, i );
376 zdepcalibtree->Fill();
377 }
378
379 // write entacalib----------------------------------------------------------------
380 double enta0;
381 double enta1;
382 double enta2;
383 double enta3;
384 TTree* entacalibtree = new TTree( "entacalib", "entaCalib" );
385 entacalibtree->Branch( "enta0", &enta0, "S/D" );
386 entacalibtree->Branch( "enta1", &enta1, "S/D" );
387 entacalibtree->Branch( "enta2", &enta2, "S/D" );
388 entacalibtree->Branch( "enta3", &enta3, "S/D" );
389 for ( i = 0; i < 43; i++ )
390 {
391 enta0 = tmpObject->getenta( 0, i );
392 enta1 = tmpObject->getenta( 1, i );
393 enta2 = tmpObject->getenta( 2, i );
394 enta3 = tmpObject->getenta( 3, i );
395 entacalibtree->Fill();
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 // write wiregcalib----------------------------------------------------------------
413 double wireg;
414 TTree* wiregcalibtree = new TTree( "wiregcalib", "wiregCalib" );
415 wiregcalibtree->Branch( "wireg", &wireg, "S/D" );
416 for ( i = 0; i < 6796; i++ )
417 {
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 {
428 layerg = tmpObject->getlayerg( i );
429 layergcalibtree->Fill();
430 }
431
432 // write all the trees
433 rungcalibtree->Write();
434 ddgcalibtree->Write();
435 ggscalibtree->Write();
436 zdepcalibtree->Write();
437 entacalibtree->Write();
438 gaintree->Write();
439 resoltree->Write();
440 wiregcalibtree->Write();
441 layergcalibtree->Write();
442
443 delete rungcalibtree;
444 delete ddgcalibtree;
445 delete ggscalibtree;
446 delete zdepcalibtree;
447 delete entacalibtree;
448 delete gaintree;
449 delete resoltree;
450 delete wiregcalibtree;
451 delete layergcalibtree;
452 closeWrite();
453 log << MSG::INFO << "successfully create RootFile" << endmsg;
454 return sc;
455}
IMessageSvc * msgSvc()
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 41 of file RootDedxCalibDataCnv.cxx.

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

◆ objType()

const CLID & RootDedxCalibDataCnv::objType ( ) const

Definition at line 37 of file RootDedxCalibDataCnv.cxx.

37{ return CLID_Calib_DedxCal; }

◆ repSvcType()

virtual long RootDedxCalibDataCnv::repSvcType ( ) const
inlinevirtual

Definition at line 35 of file RootDedxCalibDataCnv.h.

◆ 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: