BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
DstReformAlg.cxx
Go to the documentation of this file.
1#include "GaudiKernel/Bootstrap.h"
2#include "GaudiKernel/ISvcLocator.h"
3#include "GaudiKernel/MsgStream.h"
4#include "GaudiKernel/SmartIF.h"
5#include <GaudiKernel/IProperty.h>
6#include <TFile.h>
7#include <TTree.h>
8#include <iostream>
9#include <map>
10#include <vector>
11
12#include "DstReformAlg.h"
13
14#include "CLHEP/Geometry/Point3D.h"
15#ifndef ENABLE_BACKWARDS_COMPATIBILITY
16typedef HepGeom::Point3D<double> HepPoint3D;
17#endif
18
19using namespace CLHEP;
20
21/////////////////////////////////////////////////////////////////////////////
23DstReformAlg::DstReformAlg( const std::string& name, ISvcLocator* pSvcLocator )
24 : Algorithm( name, pSvcLocator ) {
25
26 declareProperty( "OutputDstFile", m_outputFile );
27}
28
29// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
31 MsgStream log( msgSvc(), name() );
32
33 log << MSG::INFO << "in initialize()" << endmsg;
34 m_evtNum = 0;
35
36 Gaudi::svcLocator()->service( "TagFilterSvc", m_tagFilterSvc ).ignore();
37
38 m_dstDataType = m_tagFilterSvc->getDstDataType();
39 std::cout << "init DstReformAlg, m_dstDataType: " << m_dstDataType << std::endl;
40
41 return StatusCode::SUCCESS;
42}
43
45 UInt_t tagData0 = m_tagFilterSvc->getTagData0();
46 UInt_t tagData1 = m_tagFilterSvc->getTagData1();
47 UInt_t tagData2 = m_tagFilterSvc->getTagData2();
48 UInt_t tagData3 = m_tagFilterSvc->getTagData3();
49 UInt_t tagData4 = m_tagFilterSvc->getTagData4();
50
51 // tagData0 and tagData1 is for DST reform, tagData0 is first priority
52
53 // for DTagSetAlg and LTagSetAlg, if tagData0=3, then reform these events with
54 // tagData1(nGoodCharged).
55 int dstDataType = m_tagFilterSvc->getDstDataType();
56 if ( ( dstDataType == 2 || dstDataType == 4 ) && ( tagData0 == 3 ) )
57 { cmap[tagData0 + tagData1].push_back( m_evtNum ); }
58 else cmap[tagData0].push_back( m_evtNum );
59
60 // add tag information to each event
61 if ( m_dstDataType == 1 )
62 {
63 evtmap[m_evtNum].push_back( tagData0 );
64 evtmap[m_evtNum].push_back( tagData1 );
65 evtmap[m_evtNum].push_back( tagData2 );
66 evtmap[m_evtNum].push_back( tagData3 );
67 evtmap[m_evtNum].push_back( tagData4 );
68 }
69 else if ( ( m_dstDataType == 2 ) || ( m_dstDataType == 3 ) || ( m_dstDataType == 4 ) )
70 {
71 UInt_t tagData5 = m_tagFilterSvc->getTagData5();
72 UInt_t tagData6 = m_tagFilterSvc->getTagData6();
73 UInt_t tagData7 = m_tagFilterSvc->getTagData7();
74 UInt_t tagData8 = m_tagFilterSvc->getTagData8();
75 UInt_t tagData9 = m_tagFilterSvc->getTagData9();
76 evtmap[m_evtNum].push_back( tagData1 );
77 evtmap[m_evtNum].push_back( tagData2 );
78 evtmap[m_evtNum].push_back( tagData3 );
79 evtmap[m_evtNum].push_back( tagData4 );
80 evtmap[m_evtNum].push_back( tagData5 );
81 evtmap[m_evtNum].push_back( tagData6 );
82 evtmap[m_evtNum].push_back( tagData7 );
83 evtmap[m_evtNum].push_back( tagData8 );
84 evtmap[m_evtNum].push_back( tagData9 );
85 }
86
87 m_evtNum++;
88
89 return StatusCode::SUCCESS;
90}
91
92// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
94 MsgStream log( msgSvc(), name() );
95
96 std::string fn;
97
98 auto iProp = service<IProperty>( "EventCnvSvc" );
99 if ( iProp )
100 {
101 iProp->getProperty( "digiRootInputFile", fn ).ignore();
102 std::string::size_type p1 = fn.find( '\'' ) + 1;
103 std::string::size_type p2 = fn.find( '\'', p1 );
104 fn = fn.substr( p1, p2 - p1 );
105 std::cout << "input dst file name: " << fn << std::endl;
106 }
107
108 TFile* f1 = TFile::Open( fn.c_str(), "READ" );
109 TFile* f2 = TFile::Open( m_outputFile.c_str(), "RECREATE" );
110 TTree* t1 = (TTree*)f1->Get( "Event" );
111
112 log << MSG::INFO << "in finalize()" << endmsg;
113 TTree* jt1 = (TTree*)( f1->Get( "JobInfoTree" ) );
114 TTree* jt2 = jt1->CloneTree();
115 TTree* t0 = new TTree( "Metadata", "" );
116 Int_t mode = -1, begin = 0, end = 0;
117 t0->Branch( "mode", &mode, "mode/I" );
118 t0->Branch( "begin", &begin, "begin/I" );
119 t0->Branch( "end", &end, "end/I" );
120
121 for ( std::map<int, std::vector<long int>>::iterator it = cmap.begin(); it != cmap.end();
122 ++it )
123 {
124 mode = it->first;
125 begin = end;
126 end += ( it->second ).size();
127 std::cout << "mode: " << mode << " begin: " << begin << " end: " << end << std::endl;
128 t0->Fill();
129 }
130
131 TTree* t2 = t1->CloneTree( 0 );
132 Int_t nEvents = t1->GetEntries();
133 std::cout << "Total nEvents: " << nEvents << std::endl;
134 Int_t origEntry = 0, currEntry = 0;
135 UInt_t tagData0, tagData1, tagData2, tagData3, tagData4;
136 UInt_t tagData5, tagData6, tagData7, tagData8;
137 TTree* t3 = new TTree( "Entries", "" );
138 t3->Branch( "origEntry", &origEntry, "origEntry/I" );
139 t3->Branch( "currEntry", &currEntry, "currEntry/I" );
140 t3->Branch( "tagData0", &tagData0, "tagData0/i" );
141 t3->Branch( "tagData1", &tagData1, "tagData1/i" );
142 t3->Branch( "tagData2", &tagData2, "tagData2/i" );
143 t3->Branch( "tagData3", &tagData3, "tagData3/i" );
144 t3->Branch( "tagData4", &tagData4, "tagData4/i" );
145 if ( ( m_dstDataType == 2 ) || ( m_dstDataType == 3 ) || ( m_dstDataType == 4 ) )
146 {
147 t3->Branch( "tagData5", &tagData5, "tagData5/i" );
148 t3->Branch( "tagData6", &tagData6, "tagData6/i" );
149 t3->Branch( "tagData7", &tagData7, "tagData7/i" );
150 t3->Branch( "tagData8", &tagData8, "tagData8/i" );
151 }
152
153 for ( std::map<int, std::vector<long int>>::iterator it = cmap.begin(); it != cmap.end();
154 ++it )
155 {
156 int nNt = it->second.size();
157 std::cout << "Writing " << nNt << "\t events with " << it->first << " TagData0..."
158 << std::endl;
159 for ( int i = 0; i < nNt; ++i )
160 {
161 origEntry = ( it->second )[i];
162 t1->GetEntry( origEntry );
163 tagData0 = ( evtmap[origEntry] )[0];
164 tagData1 = ( evtmap[origEntry] )[1];
165 tagData2 = ( evtmap[origEntry] )[2];
166 tagData3 = ( evtmap[origEntry] )[3];
167 tagData4 = ( evtmap[origEntry] )[4];
168 if ( ( m_dstDataType == 2 ) || ( m_dstDataType == 3 ) || ( m_dstDataType == 4 ) )
169 {
170 tagData5 = ( evtmap[origEntry] )[5];
171 tagData6 = ( evtmap[origEntry] )[6];
172 tagData7 = ( evtmap[origEntry] )[7];
173 tagData8 = ( evtmap[origEntry] )[8];
174 }
175 t2->Fill();
176 t3->Fill();
177 currEntry++;
178 }
179 ( it->second ).clear();
180 }
181
182 f2->Write();
183 f2->Close();
184 return StatusCode::SUCCESS;
185}
double p2[4]
double p1[4]
DECLARE_COMPONENT(BesBdkRc)
HepGeom::Point3D< double > HepPoint3D
TFile * f1
IMessageSvc * msgSvc()
DstReformAlg(const std::string &name, ISvcLocator *pSvcLocator)
StatusCode finalize()
StatusCode execute()
StatusCode initialize()