CGEM BOSS 6.6.5.g
BESIII Offline Software System
Loading...
Searching...
No Matches
RootMdcCalibDataCnv.cxx
Go to the documentation of this file.
1// $Header: /bes/bes/BossCvs/Calibration/CalibSvc/CalibROOTCnv/src/cnv/RootMdcCalibDataCnv.cxx,v 1.9 2011/12/20 22:33:05 zhangy Exp $
2#include "GaudiKernel/MsgStream.h"
7
8#include "TFile.h"
9#include "TTree.h"
10#include "TDirectory.h"
11#include "TObject.h"
12#include "TBufferFile.h"
13#include "TObjArray.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<RootMdcCalibDataCnv> MdcCal_factory;
30//const ICnvFactory& RootMdcCalibDataCnvFactory = MdcCal_factory;
31
32
33
36
37}
38
39
40const CLID& RootMdcCalibDataCnv::objType() const {
41 return CLID_Calib_MdcCal;
42}
43
45 return CLID_Calib_MdcCal;
46}
47
48StatusCode RootMdcCalibDataCnv::i_createObj(const std::string& fname,
49 DataObject*& refpObject) {
50
51 MsgStream log(msgSvc(), "RootMdcCalibDataCnv");
52 log<<MSG::DEBUG<<"SetProperty"<<endreq;
53
54 StatusCode sc = openRead(fname);
55 if(!sc)
56 { log<<MSG::ERROR<<"unable to open files"<<endreq;
57 }
58
60 // Read in our object
61 int i;
62 int nentries;
63 // read xttree ------------------------------------------------------------
64 double xtpar;
65 int xtkey;
66 TTree *xttree = (TTree*)m_inFile -> Get("XtTree");
67 xttree -> SetBranchAddress("xtpar", &xtpar);
68 xttree -> SetBranchAddress("xtkey", &xtkey);
69 nentries = xttree -> GetEntries();
70 for(i=0; i<nentries; i++){
71 xttree -> GetEntry(i);
72 tmpObject -> setXtpar(xtkey,xtpar);
73 }
74
75 // read newxttree ------------------------------------------------------------
76 TObjArray newXtTrees;
77 if(NULL!=m_inFile->Get("trNewXt00_00_0")){
78 for(int layid=0; layid<43; layid++){
79 for(int entr=0; entr<18; entr++){
80 for(int lr=0; lr<2; lr++){
81 char newXtTreeName[20];
82 sprintf(newXtTreeName,"trNewXt%02d_%02d_%d",layid,entr,lr);
83 TTree* newXtTree = ((TTree*)m_inFile->Get(newXtTreeName));
84 if(!newXtTree) break;
85 newXtTrees.Add(newXtTree->CloneTree());
86 }//end lr
87 }//end entr
88 }//end layid
89 if((43*18*2)==newXtTrees.GetEntries())tmpObject->setNewXtpar(&newXtTrees);
90 }
91
92
93 // read r2ttree ------------------------------------------------------------
94 TObjArray r2tTrees;
95 if(NULL!=m_inFile->Get("r2t00")){
96 for(int layid=0; layid<43; layid++){
97 char r2tTreeName[20];
98 sprintf(r2tTreeName,"r2t%02d",layid);
99 TTree* r2tTree = ((TTree*)m_inFile->Get(r2tTreeName));
100 if(!r2tTree) break;
101 r2tTrees.Add(r2tTree->CloneTree());
102 }
103 if(43==r2tTrees.GetEntries()) tmpObject->setR2tpar(&r2tTrees);
104 }
105
106 // read t0tree ------------------------------------------------------------
107 double t0;
108 double delt0;
109 TTree *t0tree = (TTree*)m_inFile -> Get("T0Tree");
110 t0tree -> SetBranchAddress("t0", &t0);
111 t0tree -> SetBranchAddress("delt0", &delt0);
112 nentries = t0tree -> GetEntries();
113 for(i=0; i<nentries; i++){
114 t0tree -> GetEntry(i);
115 tmpObject -> setT0(t0);
116 tmpObject -> setDelT0(delt0);
117 }
118
119 // read qttree ------------------------------------------------------------
120 double qtpar0;
121 double qtpar1;
122 TTree *qttree = (TTree*)m_inFile -> Get("QtTree");
123 qttree -> SetBranchAddress("qtpar0", &qtpar0);
124 qttree -> SetBranchAddress("qtpar1", &qtpar1);
125 nentries = qttree -> GetEntries();
126 for(i=0; i<nentries; i++){
127 qttree -> GetEntry(i);
128 tmpObject -> setQtpar0(qtpar0);
129 tmpObject -> setQtpar1(qtpar1);
130 }
131
132 // read Sdtree ---------------------------------------------------------
133 double sdpar;
134 int sdkey;
135 TTree *sdtree = (TTree*)m_inFile -> Get("SdTree");
136 sdtree -> SetBranchAddress("sdpar", &sdpar);
137 sdtree -> SetBranchAddress("sdkey", &sdkey);
138 nentries = sdtree -> GetEntries();
139
140 for(i=0; i<nentries; i++){
141 sdtree -> GetEntry(i);
142 tmpObject -> setSdpar(sdkey,sdpar);
143 }
144
145 refpObject=tmpObject;
146 return StatusCode::SUCCESS;
147}
148
149StatusCode RootMdcCalibDataCnv::createRoot(const std::string& fname,
150 CalibData::CalibBase1* pTDSObj) {
151
152 MsgStream log(msgSvc(), "RootMdcCalibDataCnv");
153 // Open the file, create the branch
154 StatusCode sc = openWrite(fname);
155 if(!sc)
156 { log<<MSG::ERROR<<"unable to open files"<<endreq;
157 }
158 // write the Data in the TCDS to RootFile
159 CalibData::MdcCalibData* tmpObject = dynamic_cast<CalibData::MdcCalibData*>(pTDSObj);
160 int tmpNo;
161 int key;
162 double xtpar;
163 double t0;
164 double delt0;
165 double qtpar[2];
166 double sdpar;
167 int i;
168
169 //xttree------------------------------------------------------------------
170 TTree *xttree = new TTree("XtTree", "XtTree");
171 xttree -> Branch("xtkey", &key, "key/I");
172 xttree -> Branch("xtpar", &xtpar, "xtpar/D");
173 tmpObject -> setXtBegin();
174 while( tmpObject -> getNextXtpar(key, xtpar) ){
175 xttree -> Fill();
176 }
177 //newxttree------------------------------------------------------------------
178 //FIXME
179 /*
180 TObjArray newxttrees_1;
181 //TBufferFile newxttrees_buf(TBuffer::kWrite);
182 int nTree=0;
183 for(int t_layer=0; t_layer<43; t_layer++){
184 for(int t_bin=0; t_bin<18; t_bin++){
185 for(int t_lr=0; t_lr<2; t_lr++){
186 char newXtTreeName[20];
187 sprintf(newXtTreeName,"trNewXt%02d_%02d_%d",t_layer,t_bin,t_lr);
188 TTree *newXttree = ((TTree*)f->Get(newXtTreeName))->CloneTree();
189 newxttrees_1.Add(newXttree);
190 nTree++;
191 }
192 }
193 }
194 newxttrees_1.Streamer(newxttrees_buf);
195 */
196
197
198 //t0tree-------------------------------------------------------------------
199 TTree *t0tree = new TTree("T0Tree", "T0Tree");
200 t0tree -> Branch("t0", &t0, "t0/D");
201 t0tree -> Branch("delt0", &delt0, "delt0/D");
202 tmpNo = tmpObject -> gett0No();
203 for(i=0; i<tmpNo; i++){
204 t0 = tmpObject -> getT0(i);
205 delt0 = tmpObject -> getDelT0(i);
206 t0tree -> Fill();
207 }
208
209 //qttree--------------------------------------------------------------------
210 TTree *qttree = new TTree("QtTree", "QtTree");
211 qttree -> Branch("qtpar0", &(qtpar[0]), "qtpar0/D");
212 qttree -> Branch("qtpar1", &(qtpar[1]), "qtpar1/D");
213 tmpNo = tmpObject -> getqtparNo();
214 for(i=0; i<tmpNo; i++){
215 qtpar[0] = tmpObject -> getQtpar0(i);
216 qtpar[1] = tmpObject -> getQtpar1(i);
217 qttree -> Fill();
218 }
219
220 //sdtree--------------------------------------------------------------------
221 TTree *sdtree = new TTree("SdTree", "SdTree");
222 sdtree -> Branch("sdkey", &key, "key/I");
223 sdtree -> Branch("sdpar", &sdpar, "sdpar/D");
224 tmpObject -> setSdBegin();
225 while( tmpObject -> getNextSdpar(key, sdpar) ){
226 sdtree -> Fill();
227 }
228
229 xttree -> Write();
230 t0tree -> Write();
231 qttree -> Write();
232 sdtree -> Write();
233
234 delete xttree;
235 delete t0tree;
236 delete qttree;
237 delete sdtree;
238
239 closeWrite();
240 log<<MSG::INFO<<"successfully create RootFile"<<endreq;
241 return sc;
242}
const CLID CLID_Calib_MdcCal
Definition: CalibModel.h:41
data SetBranchAddress("time",&time)
data GetEntry(0)
Int_t nentries
IMessageSvc * msgSvc()
*************DOUBLE PRECISION m_pi *DOUBLE PRECISION m_HvecTau2 DOUBLE PRECISION m_HvClone2 DOUBLE PRECISION m_gamma1 DOUBLE PRECISION m_gamma2 DOUBLE PRECISION m_thet1 DOUBLE PRECISION m_thet2 INTEGER m_IFPHOT *COMMON c_Taupair $ !Spin Polarimeter vector first Tau $ !Spin Polarimeter vector second Tau $ !Clone Spin Polarimeter vector first Tau $ !Clone Spin Polarimeter vector second Tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning nd tau $ !Random Euler angle for cloning nd tau $ !Random Euler angle for cloning nd tau $ !phi of HvecTau1 $ !theta of HvecTau1 $ !phi of HvecTau2 $ !theta of HvecTau2 $ !super key
Definition: Taupair.h:42
void setR2tpar(TObjArray *r2tTrees)
void setNewXtpar(TObjArray *newXtTrees)
StatusCode openRead(const std::string &fname)
StatusCode closeWrite()
virtual StatusCode openWrite(const std::string &fname)
const CLID & objType() const
virtual StatusCode createRoot(const std::string &fname, CalibData::CalibBase1 *pTDSObj)
virtual StatusCode i_createObj(const std::string &fname, DataObject *&refpObject)
static const CLID & classID()
RootMdcCalibDataCnv(ISvcLocator *svc)