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