BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
TreeMdcCalibDataCnv.cxx
Go to the documentation of this file.
2#include "CalibData/Mdc/MdcCalibData.h"
3#include "CalibDataSvc/IInstrumentName.h"
4#include "CalibMySQLCnvSvc/TreeAddress.h"
5#include "GaudiKernel/MsgStream.h"
6#include "TBuffer.h"
7#include "TDirectory.h"
8#include "TFile.h"
9#include "TObjArray.h"
10#include "TObject.h"
11#include "TTree.h"
12
13#include "GaudiKernel/DataObject.h"
14#include "GaudiKernel/GenericAddress.h"
15#include "GaudiKernel/IAddressCreator.h"
16#include "GaudiKernel/IConversionSvc.h"
17#include "GaudiKernel/IDataProviderSvc.h"
18#include "GaudiKernel/IOpaqueAddress.h"
19
20#include "CalibDataSvc/ICalibMetaCnvSvc.h"
21#include "CalibDataSvc/ICalibTreeSvc.h" //maybe
22
23// Temporary. Hope to find a better way to do this
24#include "CalibData/CalibModel.h"
25using namespace CalibData;
26// static CnvFactory<TreeMdcCalibDataCnv> MdcCal_factory;
27// const ICnvFactory& TreeMdcCalibDataCnvFactory = MdcCal_factory;
28
31
32const CLID& TreeMdcCalibDataCnv::objType() const { return CLID_Calib_MdcCal; }
33
35
36StatusCode TreeMdcCalibDataCnv::i_createObj( IOpaqueAddress* addr, DataObject*& refpObject ) {
37
38 MsgStream log( msgSvc(), "TreeMdcCalibDataCnv" );
39 log << MSG::DEBUG << "SetProperty" << endmsg;
41 TreeAddress* add = dynamic_cast<TreeAddress*>( addr );
42 DatabaseRecord* records = add->pp();
43
44 TBufferFile* buf1 = new TBufferFile( TBuffer::kRead );
45 TBufferFile* buf2 = new TBufferFile( TBuffer::kRead );
46 TBufferFile* buf3 = new TBufferFile( TBuffer::kRead );
47 TBufferFile* buf4 = new TBufferFile( TBuffer::kRead );
48 TBufferFile* buf5 = new TBufferFile( TBuffer::kRead );
49 TBufferFile* buf6 = new TBufferFile( TBuffer::kRead );
50
51 buf1->SetBuffer( ( *records )["XtTree"], 512000, kFALSE );
52 buf2->SetBuffer( ( *records )["QtTree"], 512000, kFALSE );
53 buf3->SetBuffer( ( *records )["T0Tree"], 512000, kFALSE );
54 buf4->SetBuffer( ( *records )["SdTree"], 512000, kFALSE );
55 buf5->SetBuffer( ( *records )["NewXtTrees"], 51200000, kFALSE );
56 buf6->SetBuffer( ( *records )["R2tTrees"], 25600000, kFALSE );
57
58 std::cout << " SftVer is " << ( *records )["SftVer"] << std::endl;
59 std::cout << "TreeMdcCalibDataCnv: CalVerSft is " << ( *records )["CalParVer"] << std::endl;
60 std::cout << "TreeMdcCalibDataCnv: Calib file name is " << ( *records )["FileName"]
61 << std::endl;
62
63 bool debug = std::string( ( *records )["FileName"] ) ==
64 std::string( "MdcCalibConst_9810_10878_662b.root" );
65
66 TTree* xttree = new TTree();
67 xttree->Streamer( *buf1 );
68
69 TTree* qttree = new TTree();
70 qttree->Streamer( *buf2 );
71
72 TTree* t0tree = new TTree();
73 t0tree->Streamer( *buf3 );
74
75 TTree* sdtree = new TTree();
76 sdtree->Streamer( *buf4 );
77
78 TObjArray newxttrees;
79 DatabaseRecord::iterator it = ( *records ).find( "NewXtTrees" );
80 if ( it != ( *records ).end() )
81 {
82 if ( ( *it ).second != NULL ) { newxttrees.Streamer( *buf5 ); }
83 }
84
85 TObjArray r2ttrees;
86 it = ( *records ).find( "R2tTrees" );
87 if ( it != ( *records ).end() )
88 {
89 if ( ( *it ).second != NULL ) { r2ttrees.Streamer( *buf6 ); }
90 }
91
92 // Read in the object
93 int i;
94 int nentries;
95
96 // read xttree-----------------------------
97 double xtpar;
98 int xtkey;
99 xttree->SetBranchAddress( "xtpar", &xtpar );
100 xttree->SetBranchAddress( "xtkey", &xtkey );
101 nentries = xttree->GetEntries();
102 for ( i = 0; i < nentries; i++ )
103 {
104 xttree->GetEntry( i );
105 tmpObject->setXtpar( xtkey, xtpar );
106 }
107
108 // read newxttrees-----------------------------
109 if ( ( 43 * 18 * 2 ) == newxttrees.GetEntries() )
110 {
111 tmpObject->setNewXtpar( &newxttrees );
112 for ( int i = 0; i < 43 * 18 * 2; i++ )
113 {
114 TTree* tempTree = (TTree*)newxttrees.At( i );
115 delete tempTree;
116 }
117 }
118 // read r2ttrees-----------------------------
119 if ( 43 == r2ttrees.GetEntries() )
120 {
121 tmpObject->setR2tpar( &r2ttrees );
122 for ( int i = 0; i < 43; i++ )
123 {
124 TTree* tempTree = (TTree*)r2ttrees.At( i );
125 delete tempTree;
126 }
127 }
128
129 // read t0tree ------------------------------------------------------------
130 double t0;
131 double delt0;
132 t0tree->SetBranchAddress( "t0", &t0 );
133 t0tree->SetBranchAddress( "delt0", &delt0 );
134 nentries = t0tree->GetEntries();
135 for ( i = 0; i < nentries; i++ )
136 {
137 t0tree->GetEntry( i );
138 tmpObject->setT0( t0 );
139 tmpObject->setDelT0( delt0 );
140 }
141 // read qttree ------------------------------------------------------------
142 double qtpar0;
143 double qtpar1;
144 qttree->SetBranchAddress( "qtpar0", &qtpar0 );
145 qttree->SetBranchAddress( "qtpar1", &qtpar1 );
146 nentries = qttree->GetEntries();
147 for ( i = 0; i < nentries; i++ )
148 {
149 qttree->GetEntry( i );
150 tmpObject->setQtpar0( qtpar0 );
151 tmpObject->setQtpar1( qtpar1 );
152 }
153
154 // read Sdtree ---------------------------------------------------------
155 double sdpar;
156 int sdkey;
157 sdtree->SetBranchAddress( "sdpar", &sdpar );
158 sdtree->SetBranchAddress( "sdkey", &sdkey );
159 nentries = sdtree->GetEntries();
160
161 for ( i = 0; i < nentries; i++ )
162 {
163 sdtree->GetEntry( i );
164 tmpObject->setSdpar( sdkey, sdpar );
165 }
166
167 refpObject = tmpObject;
168 delete xttree;
169 delete qttree;
170 delete t0tree;
171 delete sdtree;
172
173 delete buf1;
174 delete buf2;
175 delete buf3;
176 delete buf4;
177 delete buf5;
178 delete buf6;
179
180 return StatusCode::SUCCESS;
181}
Int_t nentries
IMessageSvc * msgSvc()
void setR2tpar(TObjArray *r2tTrees)
void setQtpar1(double val)
void setXtpar(int xtkey, double val)
void setNewXtpar(TObjArray *newXtTrees)
void setSdpar(int sdkey, double val)
void setDelT0(double val)
void setQtpar0(double val)
TreeCalBaseCnv(ISvcLocator *svc, const CLID &clid)
static const CLID & classID()
TreeMdcCalibDataCnv(ISvcLocator *svc)
const CLID & objType() const
virtual StatusCode i_createObj(IOpaqueAddress *address, DataObject *&refpObject)