BOSS 7.1.2
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.10 2023/10/06 22:57:09 maqm 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
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 //reduce Memerey use
102 r2tTree->SetCacheSize(2500);
103 r2tTrees.Add(r2tTree->CloneTree());
104 }
105 if(43==r2tTrees.GetEntries()) tmpObject->setR2tpar(&r2tTrees);
106 }
107
108 // read t0tree ------------------------------------------------------------
109 double t0;
110 double delt0;
111 TTree *t0tree = (TTree*)m_inFile -> Get("T0Tree");
112 t0tree -> SetBranchAddress("t0", &t0);
113 t0tree -> SetBranchAddress("delt0", &delt0);
114 nentries = t0tree -> GetEntries();
115 for(i=0; i<nentries; i++){
116 t0tree -> GetEntry(i);
117 tmpObject -> setT0(t0);
118 tmpObject -> setDelT0(delt0);
119 }
120
121 // read qttree ------------------------------------------------------------
122 double qtpar0;
123 double qtpar1;
124 TTree *qttree = (TTree*)m_inFile -> Get("QtTree");
125 qttree -> SetBranchAddress("qtpar0", &qtpar0);
126 qttree -> SetBranchAddress("qtpar1", &qtpar1);
127 nentries = qttree -> GetEntries();
128 for(i=0; i<nentries; i++){
129 qttree -> GetEntry(i);
130 tmpObject -> setQtpar0(qtpar0);
131 tmpObject -> setQtpar1(qtpar1);
132 }
133
134 // read Sdtree ---------------------------------------------------------
135 double sdpar;
136 int sdkey;
137 TTree *sdtree = (TTree*)m_inFile -> Get("SdTree");
138 sdtree -> SetBranchAddress("sdpar", &sdpar);
139 sdtree -> SetBranchAddress("sdkey", &sdkey);
140 nentries = sdtree -> GetEntries();
141
142 for(i=0; i<nentries; i++){
143 sdtree -> GetEntry(i);
144 tmpObject -> setSdpar(sdkey,sdpar);
145 }
146
147 refpObject=tmpObject;
148 return StatusCode::SUCCESS;
149}
150
151StatusCode RootMdcCalibDataCnv::createRoot(const std::string& fname,
152 CalibData::CalibBase1* pTDSObj) {
153
154 MsgStream log(msgSvc(), "RootMdcCalibDataCnv");
155 // Open the file, create the branch
156 StatusCode sc = openWrite(fname);
157 if(!sc)
158 { log<<MSG::ERROR<<"unable to open files"<<endreq;
159 }
160 // write the Data in the TCDS to RootFile
161 CalibData::MdcCalibData* tmpObject = dynamic_cast<CalibData::MdcCalibData*>(pTDSObj);
162 int tmpNo;
163 int key;
164 double xtpar;
165 double t0;
166 double delt0;
167 double qtpar[2];
168 double sdpar;
169 int i;
170
171 //xttree------------------------------------------------------------------
172 TTree *xttree = new TTree("XtTree", "XtTree");
173 xttree -> Branch("xtkey", &key, "key/I");
174 xttree -> Branch("xtpar", &xtpar, "xtpar/D");
175 tmpObject -> setXtBegin();
176 while( tmpObject -> getNextXtpar(key, xtpar) ){
177 xttree -> Fill();
178 }
179 //newxttree------------------------------------------------------------------
180 //FIXME
181 /*
182 TObjArray newxttrees_1;
183 //TBufferFile newxttrees_buf(TBuffer::kWrite);
184 int nTree=0;
185 for(int t_layer=0; t_layer<43; t_layer++){
186 for(int t_bin=0; t_bin<18; t_bin++){
187 for(int t_lr=0; t_lr<2; t_lr++){
188 char newXtTreeName[20];
189 sprintf(newXtTreeName,"trNewXt%02d_%02d_%d",t_layer,t_bin,t_lr);
190 TTree *newXttree = ((TTree*)f->Get(newXtTreeName))->CloneTree();
191 newxttrees_1.Add(newXttree);
192 nTree++;
193 }
194 }
195 }
196 newxttrees_1.Streamer(newxttrees_buf);
197 */
198
199
200 //t0tree-------------------------------------------------------------------
201 TTree *t0tree = new TTree("T0Tree", "T0Tree");
202 t0tree -> Branch("t0", &t0, "t0/D");
203 t0tree -> Branch("delt0", &delt0, "delt0/D");
204 tmpNo = tmpObject -> gett0No();
205 for(i=0; i<tmpNo; i++){
206 t0 = tmpObject -> getT0(i);
207 delt0 = tmpObject -> getDelT0(i);
208 t0tree -> Fill();
209 }
210
211 //qttree--------------------------------------------------------------------
212 TTree *qttree = new TTree("QtTree", "QtTree");
213 qttree -> Branch("qtpar0", &(qtpar[0]), "qtpar0/D");
214 qttree -> Branch("qtpar1", &(qtpar[1]), "qtpar1/D");
215 tmpNo = tmpObject -> getqtparNo();
216 for(i=0; i<tmpNo; i++){
217 qtpar[0] = tmpObject -> getQtpar0(i);
218 qtpar[1] = tmpObject -> getQtpar1(i);
219 qttree -> Fill();
220 }
221
222 //sdtree--------------------------------------------------------------------
223 TTree *sdtree = new TTree("SdTree", "SdTree");
224 sdtree -> Branch("sdkey", &key, "key/I");
225 sdtree -> Branch("sdpar", &sdpar, "sdpar/D");
226 tmpObject -> setSdBegin();
227 while( tmpObject -> getNextSdpar(key, sdpar) ){
228 sdtree -> Fill();
229 }
230
231 xttree -> Write();
232 t0tree -> Write();
233 qttree -> Write();
234 sdtree -> Write();
235
236 delete xttree;
237 delete t0tree;
238 delete qttree;
239 delete sdtree;
240
241 closeWrite();
242 log<<MSG::INFO<<"successfully create RootFile"<<endreq;
243 return sc;
244}
const CLID CLID_Calib_MdcCal
Definition CalibModel.h:41
curve Branch("CurveSize",&CurveSize,"CurveSize/I")
curve Fill()
curve Write()
sprintf(cut,"kal_costheta0_em>-0.93&&kal_costheta0_em<0.93&&kal_pxy0_em>=0.05+%d*0.1&&kal_pxy0_em<0.15+%d*0.1&&NGch>=2", j, j)
data GetEntry(0)
data SetBranchAddress("time",&time)
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
#define NULL
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)