BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
DedxSaveWireGain.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2
3#include "TCanvas.h"
4#include "TF1.h"
5#include "TFile.h"
6#include "TH1F.h"
7#include "TStyle.h"
8#include "TTree.h"
9#include <sstream>
10#include <string>
11
14using namespace std;
15
16DedxSaveWireGain::DedxSaveWireGain( const std::string& name, ISvcLocator* pSvcLocator )
17 : DedxCalib( name, pSvcLocator ) {
18 declareProperty( "PMin", pMin = 0.2 );
19 declareProperty( "PMax", pMax = 2.3 );
20 f_save = new TFile( "wiregain.root", "recreate" );
21 m_wiregain = new TTree( "wiregain", "wiregain" );
22}
23
25 MsgStream log( msgSvc(), name() );
26 log << MSG::INFO << "DedxSaveWireGain::BookHists()" << endmsg;
27
29}
30
32 MsgStream log( msgSvc(), name() );
33 log << MSG::INFO << "DedxSaveWireGain::FillHists()" << endmsg;
34 int run, evt;
35 double gain;
36 m_wiregain->Branch( "run", &run, "run/I" );
37 m_wiregain->Branch( "evt", &evt, "evt/I" );
38 m_wiregain->Branch( "gain", &gain, "gain/D" );
39
40 TFile* f;
41 TTree* n102;
42 string runlist;
43
44 double dedx = 0;
45 float runNO = 0, evtNO = 0, runFlag = 0, pathlength = 0, wid = 0, layid = 0, dd_in = 0,
46 driftdist = 0, eangle = 0, zhit = 0, costheta = 0, tes = 0, ptrk = 0;
47 cout << "m_recFileVector.size()= " << m_recFileVector.size() << endl;
48 for ( int i = 0; i < m_recFileVector.size(); i++ )
49 {
50 runlist = m_recFileVector[i];
51 cout << "runlist: " << runlist.c_str() << endl;
52 f = new TFile( runlist.c_str() );
53 n102 = (TTree*)f->Get( "n102" );
54 n102->SetBranchAddress( "adc_raw", &dedx );
55 n102->SetBranchAddress( "path_rphi", &pathlength );
56 n102->SetBranchAddress( "wire", &wid );
57 n102->SetBranchAddress( "layer", &layid );
58 n102->SetBranchAddress( "doca_in", &dd_in );
59 n102->SetBranchAddress( "driftdist", &driftdist );
60 n102->SetBranchAddress( "eangle", &eangle );
61 n102->SetBranchAddress( "zhit", &zhit );
62 n102->SetBranchAddress( "runNO", &runNO );
63 n102->SetBranchAddress( "evtNO", &evtNO );
64 n102->SetBranchAddress( "runFlag", &runFlag );
65 n102->SetBranchAddress( "costheta1", &costheta );
66 n102->SetBranchAddress( "t01", &tes );
67 n102->SetBranchAddress( "ptrk1", &ptrk );
68
69 cout << "entries in this file" << n102->GetEntries() << endl;
70 for ( int j = 0; j < n102->GetEntries(); j++ )
71 {
72 n102->GetEntry( j );
73 if ( tes > 1400 ) continue;
74 if ( ptrk > pMax || ptrk < pMin ) continue;
75 if ( layid < 8 )
76 {
77 if ( dedx < MinADCValuecut || dedx > MaxADCValuecut || pathlength > RPhi_PathMaxCut ||
78 pathlength < Iner_RPhi_PathMinCut || driftdist > Iner_DriftDistCut )
79 continue;
80 }
81 else
82 {
83 if ( dedx < MinADCValuecut || dedx > MaxADCValuecut || pathlength > RPhi_PathMaxCut ||
84 pathlength < Out_RPhi_PathMinCut || driftdist > Out_DriftDistCut )
85 continue;
86 }
87 dedx = exsvc->StandardHitCorrec( 0, runFlag, 2, runNO, evtNO, pathlength, wid, layid,
88 dedx, dd_in, eangle, costheta );
89 dedx = exsvc->StandardTrackCorrec( 0, runFlag, 2, runNO, evtNO, dedx, costheta, tes );
90 run = (int)runNO;
91 evt = (int)evtNO;
92 gain = dedx;
93 // cout << "run " << run << " evt " << evt << " gain " << gain << endl;
94 m_wiregain->Fill();
95 } // end of loop in n102 entries
96 } // end of loop in rec files
97} // end of this function
98
100 f_save->cd();
101 m_wiregain->Write();
102 f_save->Save();
103 f_save->Close();
104}
DECLARE_COMPONENT(BesBdkRc)
TFile f("ana_bhabha660a_dqa_mcPat_zy_old.root")
IMessageSvc * msgSvc()
DedxCalib(const std::string &name, ISvcLocator *pSvcLocator)
Definition DedxCalib.cxx:12
void ReadRecFileList()
Definition DedxCalib.cxx:84
IDedxCorrecSvc * exsvc
Definition DedxCalib.h:29
vector< string > m_recFileVector
Definition DedxCalib.h:50
DedxSaveWireGain(const std::string &name, ISvcLocator *pSvcLocator)