BOSS 7.0.6
BESIII Offline Software System
Loading...
Searching...
No Matches
RootMucCalibDataCnv.cxx
Go to the documentation of this file.
1// $Header: /bes/bes/BossCvs/Calibration/CalibSvc/CalibROOTCnv/src/cnv/RootMucCalibDataCnv.cxx,v 1.4 2008/10/22 06:33:46 huangb Exp $
2#include <iostream>
3
4#include "GaudiKernel/MsgStream.h"
8
9#include "TFile.h"
10#include "TTree.h"
11#include "TH1F.h"
12#include "TDirectory.h"
13#include "TObject.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
28//using namespace CalibData;
31
32using namespace std;
33
34//static CnvFactory<RootMucCalibDataCnv> MucCal_factory;
35//const ICnvFactory& RootMucCalibDataCnvFactory = MucCal_factory;
36
37
38
41 m_ptrIdTr = new MucIdTransform();
42}
43
44
45const CLID& RootMucCalibDataCnv::objType() const {
46 return CLID_Calib_MucCal;
47}
48
50 return CLID_Calib_MucCal;
51}
52
53StatusCode RootMucCalibDataCnv::i_createObj(const std::string& fname,
54 DataObject*& refpObject) {
55
56 MsgStream log(msgSvc(), "RootMucCalibDataCnv");
57 log<<MSG::DEBUG<<"SetProperty"<<endreq;
58
59 StatusCode sc = openRead(fname);
60 if(!sc)
61 { log<<MSG::ERROR<<"unable to open files"<<endreq;
62 return StatusCode::FAILURE;
63 }
64
65 MucCalibData *tmpObject = new MucCalibData() ;
66 // Read in our object
67 // int i;
68 // int nentries;
69
70 // read DigiCalibConst ------------------------------------------------------------
71
72 Double_t lay_eff, box_eff, str_eff;
73 Double_t lay_cnt, box_cnt, str_cnt;
74 Double_t lay_nos, box_nos, str_nos;
75 Double_t lay_nos_ratio, box_nos_ratio, str_nos_ratio;
76 lay_eff = box_eff = str_eff = 0.0;
77 lay_cnt = box_cnt = str_cnt = 0.0;
78 lay_nos = box_nos = str_nos = 0.0;
79 lay_nos_ratio = box_nos_ratio = str_nos_ratio = 0.0;
80
81 char name[60];
82
83 TTree* tr_Lvl[3];
84 //TTree* tr_ClstPro[2];
85
86 //TTree *ddgtree = (TTree*)m_inFile -> Get("ddgcalib");
87 tr_Lvl[0] = (TTree*)m_inFile->Get("LayConst");
88 tr_Lvl[0]->SetBranchAddress("layer_eff", &lay_eff);
89 tr_Lvl[0]->SetBranchAddress("layer_cnt", &lay_cnt);
90 tr_Lvl[0]->SetBranchAddress("layer_noise", &lay_nos);
91 tr_Lvl[0]->SetBranchAddress("layer_nosratio", &lay_nos_ratio);
92
93 tr_Lvl[1] = (TTree*)m_inFile->Get("BoxConst");
94 tr_Lvl[1]->SetBranchAddress("box_eff", &box_eff);
95 tr_Lvl[1]->SetBranchAddress("box_cnt", &box_cnt);
96 tr_Lvl[1]->SetBranchAddress("box_noise", &box_nos);
97 tr_Lvl[1]->SetBranchAddress("box_nosratio", &box_nos_ratio);
98
99 tr_Lvl[2] = (TTree*)m_inFile->Get("StrConst");
100 tr_Lvl[2]->SetBranchAddress("strip_eff", &str_eff);
101 tr_Lvl[2]->SetBranchAddress("strip_cnt", &str_cnt);
102 tr_Lvl[2]->SetBranchAddress("strip_noise", &str_nos);
103 tr_Lvl[2]->SetBranchAddress("strip_nosratio", &str_nos_ratio);
104 //tr_ClstPro[0] = (TTree*)m_inFile->Get("LayClstPro");
105 //tr_ClstPro[1] = (TTree*)m_inFile->Get("BoxClstPro");
106
107 int part, segment, layer, strip;
108 part = segment = layer = strip = 0;
109 for(int i=0; i<LAYER_MAX; i++)
110 {
111
112 tr_Lvl[0]->GetEntry(i);
113 tmpObject->setLayerEff(lay_eff, i);
114 tmpObject->setLayerCnt(lay_cnt, i);
115 tmpObject->setLayerNos(lay_nos, i);
116 tmpObject->setLayerNosRatio(lay_nos_ratio, i);
117
118 sprintf(name,"LayClstPro");
119 //tr_ClstPro[0] = (TTree*)m_inFile->Get(name);
120 for(int j=0; j<CLST_MAX; j++) {
121 //if( tr_ClstPro[0] != NULL ) tmpObject->setLayerClstPro(tr_ClstPro[0]->GetBinContent(j),i,j);
122 //else tmpObject->setLayerClstPro(DEFAULT_CLST_PRO[j],i,j);
123 tmpObject->setLayerClstPro(DEFAULT_CLST_PRO[j],i,j);
124 }
125 //log<<MSG::DEBUG<<"layer: " << i << "\t" << lay_eff <<endreq;
126 }
127
128 for(int i=0; i<BOX_MAX; i++)
129 {
130 m_ptrIdTr->SetBoxPos( i, &part, &segment, &layer );
131 tr_Lvl[1]->GetEntry(i);
132 tmpObject->setBoxEff(box_eff, part, segment, layer);
133 tmpObject->setBoxCnt(box_cnt, part, segment, layer);
134 tmpObject->setBoxNos(box_nos, part, segment, layer);
135 tmpObject->setBoxNosRatio(box_nos_ratio, part, segment, layer);
136
137 sprintf(name,"BoxClstPro_B%d",i);
138 //tr_ClstPro[1] = (TTree*)m_inFile->Get(name);
139 for(int j=0; j<CLST_MAX; j++) {
140 //if( tr_ClstPro[1] != NULL ) tmpObject->setBoxClstPro(tr_ClstPro[1]->GetBinContent(j),part,segment,layer,j);
141 //else tmpObject->setBoxClstPro(DEFAULT_CLST_PRO[j],part,segment,layer,j);
142 tmpObject->setBoxClstPro(DEFAULT_CLST_PRO[j],part,segment,layer,j);
143 }
144 //log<<MSG::DEBUG<<"box: " << i << "\t" << box_eff <<endreq;
145 }
146
147 for(int i=0; i<STRIP_MAX; i++)
148 {
149 m_ptrIdTr->SetStripPos( i, &part, &segment, &layer, &strip );
150 tr_Lvl[2]->GetEntry(i);
151 tmpObject->setStripEff(str_eff, part, segment, layer, strip);
152 tmpObject->setStripCnt(str_cnt, part, segment, layer, strip);
153 tmpObject->setStripNos(str_nos, part, segment, layer, strip);
154 tmpObject->setStripNosRatio(str_nos_ratio, part, segment, layer, strip);
155 //log<<MSG::DEBUG<<"strip: " << i << "\t" << str_eff <<endreq;
156 }
157 refpObject=tmpObject;
158
159 return StatusCode::SUCCESS;
160}
161
162StatusCode RootMucCalibDataCnv::createRoot(const std::string& fname,
163 CalibData::CalibBase1* pTDSObj) {
164
165 MsgStream log(msgSvc(), "RootMucCalibDataCnv");
166
167 // Open the file, create the branch
168/* StatusCode sc = openWrite(fname);
169 if(!sc)
170 { log<<MSG::ERROR<<"unable to open files"<<endreq;
171 }
172 // write the Data in the TCDS to RootFile
173 CalibData::MucCalibData* tmpObject = dynamic_cast<CalibData::MucCalibData*>(pTDSObj);
174 int tmpNo;
175 double MucCalibConst;
176 double EnCoeff;
177 double PosCoeff;
178 int i;
179
180 //DigiCalibConst------------------------------------------------------------------
181 TTree *Digitree = new TTree("DigiCalibConst", "DigiCalibConst");
182 Digitree -> Branch("DigiCalibConst", &MucCalibConst, "MucCalibConst/D");
183 tmpNo = tmpObject -> getDigiCalibConstNo();
184 for(i=0; i<tmpNo; i++){
185 MucCalibConst = tmpObject -> getDigiCalibConst(i);
186 Digitree -> Fill();
187 }
188
189
190 //EnCoeff-------------------------------------------------------------------------
191 TTree *Entree = new TTree("EnCoeff", "EnCoeff");
192 Entree -> Branch("EnCoeff", &EnCoeff, "EnCoeff/D");
193 tmpNo = tmpObject -> getEnCoeffNo();
194 for(i=0; i<tmpNo; i++){
195 EnCoeff = tmpObject -> getEnCoeff(i);
196 Entree -> Fill();
197 }
198
199 //PosCoeff-------------------------------------------------------------------------
200 TTree *Postree = new TTree("PosCoeff", "PosCoeff");
201 Postree -> Branch("PosCoeff", &PosCoeff, "PosCoeff/D");
202 tmpNo = tmpObject -> getPosCoeffNo();
203 for(i=0; i<tmpNo; i++){
204 PosCoeff = tmpObject -> getPosCoeff(i);
205 Postree -> Fill();
206 }
207
208
209 Digitree -> Write();
210 Entree -> Write();
211 Postree -> Write();
212
213
214 delete Digitree;
215 delete Entree;
216 delete Postree;
217
218 closeWrite();
219 log<<MSG::INFO<<"successfully create RootFile"<<endreq;
220 */
221 return StatusCode::SUCCESS;
222
223}
const CLID CLID_Calib_MucCal
Definition: CalibModel.h:49
const int STRIP_MAX
Definition: MucMappingAlg.h:24
IMessageSvc * msgSvc()
void setLayerClstPro(const double layerpro, int i, int j)
Definition: MucCalibData.h:95
void setLayerNos(const double layernos, int i)
Definition: MucCalibData.h:92
void setLayerCnt(const double layercnt, int i)
Definition: MucCalibData.h:93
void setBoxCnt(const double boxcnt, int i, int j, int k)
Definition: MucCalibData.h:99
void setStripNos(const double stripnos, int i, int j, int k, int l)
Definition: MucCalibData.h:104
void setBoxNosRatio(const double boxratio, int i, int j, int k)
Definition: MucCalibData.h:100
void setStripEff(const double stripeff, int i, int j, int k, int l)
Definition: MucCalibData.h:103
void setStripNosRatio(const double stripnosratio, int i, int j, int k, int l)
Definition: MucCalibData.h:106
void setStripCnt(const double stripcnt, int i, int j, int k, int l)
Definition: MucCalibData.h:105
void setBoxEff(const double boxeff, int i, int j, int k)
Definition: MucCalibData.h:97
void setBoxNos(const double boxnos, int i, int j, int k)
Definition: MucCalibData.h:98
void setLayerEff(const double layereff, int i)
Definition: MucCalibData.h:91
void setLayerNosRatio(const double layernosratio, int i)
Definition: MucCalibData.h:94
void setBoxClstPro(const double boxpro, int i, int j, int k, int l)
Definition: MucCalibData.h:101
bool SetBoxPos(int boxid, int *part, int *segment, int *layer)
bool SetStripPos(int stripid, int *part, int *segment, int *layer, int *subid)
StatusCode openRead(const std::string &fname)
const CLID & objType() const
virtual StatusCode i_createObj(const std::string &fname, DataObject *&refpObject)
virtual StatusCode createRoot(const std::string &fname, CalibData::CalibBase1 *pTDSObj)
RootMucCalibDataCnv(ISvcLocator *svc)
static const CLID & classID()
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)