LCIO  02.17
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
test_simtrackerhit.cc
Go to the documentation of this file.
1 // test lcio::SimTrackerHit
4 
5 #include "tutil.h"
6 #include "lcio.h"
7 
8 #include "EVENT/LCIO.h"
9 #include "IO/LCReader.h"
10 #include "IO/LCWriter.h"
11 #include "IMPL/LCEventImpl.h"
12 #include "IMPL/LCCollectionVec.h"
13 #include "IMPL/SimTrackerHitImpl.h"
14 #include "IMPL/LCFlagImpl.h"
15 //#include "UTIL/Operators.h"
16 
17 #include "UTIL/CellIDEncoder.h"
18 #include "UTIL/CellIDDecoder.h"
19 #include "UTIL/ILDConf.h"
20 #include "UTIL/LCTrackerConf.h"
21 
22 //#include <iostream>
23 
24 using namespace std ;
25 using namespace lcio ;
26 
27 //static const int NRUN = 10 ;
28 static const int NEVENT = 10 ; // events
29 static const int NHITS = 1000 ; // tracker hits per event
30 
31 static string FILEN = "simtrackerhits.slcio" ;
32 
33 // replace mytest with the name of your test
34 const static string testname="test_simtrackerhit";
35 
36 //=============================================================================
37 
38 int main(int /*argc*/, char** /*argv*/ ){
39 
40  // this should be the first line in your test
41  TEST MYTEST=TEST( testname, std::cout );
42 
43  try{
44 
45 
46  MYTEST.LOG( " writing SimTrackerHits " );
47 
48  // create sio writer
49  LCWriter* lcWrt = LCFactory::getInstance()->createLCWriter() ;
50 
51  lcWrt->open( FILEN , LCIO::WRITE_NEW ) ;
52 
53  // EventLoop - create some events and write them to the file
54  for(int i=0;i<NEVENT;i++){
55 
56  // we need to use the implementation classes here
57  LCEventImpl* evt = new LCEventImpl() ;
58 
59 
60  evt->setRunNumber( 4711 ) ;
61  evt->setEventNumber( i ) ;
62 
63  LCCollectionVec* trkHits = new LCCollectionVec( LCIO::SIMTRACKERHIT ) ;
64 
65  ILDCellIDEncoder<SimTrackerHitImpl> idEnc( "readoutUnit:8,daqChannel:16," , trkHits ) ;
66  // this is effectively the same as:
67  // CellIDEncoder<SimTrackerHitImpl> idEnc( LCTrackerCellID::encoding_string() + ",readoutUnit:8,daqChannel:16" , trkHits ) ;
68 
69 
70 
71  for(int j=0;j<NHITS;j++){
72 
73  SimTrackerHitImpl* trkHit = new SimTrackerHitImpl ;
74 
75  idEnc.reset() ;
76 
77  idEnc[ LCTrackerCellID::subdet() ] = ILDDetID::FTD ;
78  idEnc[ LCTrackerCellID::layer() ] = j % 100 ;
79  idEnc[ LCTrackerCellID::side() ] = ILDDetID::bwd ;
80  idEnc[ LCTrackerCellID::module() ] = j / 100 + 1 ;
81  idEnc[ LCTrackerCellID::sensor() ] = j % 4 ;
82  idEnc["daqChannel"] = j*8 ;
83 
84  idEnc.setCellID( trkHit ) ;
85 
86  trkHit->setEDep( i*j*117. ) ;
87 
88  if( (i+j) % 3 == 0 ) {
89  trkHit->setProducedBySecondary( true );
90  }
91 
92  if( (i+j) % 5 == 0 ) {
93  trkHit->setOverlay( true );
94  }
95 
96  double pos[3] = { 1.*i, 1.*j, 1.*i*j } ;
97  trkHit->setPosition( pos ) ;
98 
99  trkHits->addElement( trkHit ) ;
100  }
101  evt->addCollection( trkHits , "SimTrackerHits") ;
102 
103  lcWrt->writeEvent(evt) ;
104 
105  delete evt ;
106  }
107 
108 
109  lcWrt->close() ;
110 
111  MYTEST.LOG(" reading back SimTrackerHits from file " ) ;
112 
113  // create sio reader
114  LCReader* lcRdr = LCFactory::getInstance()->createLCReader() ;
115 
116  lcRdr->open( FILEN ) ;
117 
118  for(int i=0;i<NEVENT;i++){
119 
120  //std::cout << " testing event " << i << std::endl ;
121 
122  LCEvent* evt = lcRdr->readNextEvent() ;
123 
124  MYTEST( evt->getRunNumber() , 4711 , " run number " ) ;
125 
126  MYTEST( evt->getEventNumber() , i , " event number " ) ;
127 
128  LCCollection* trkHits = evt->getCollection( "SimTrackerHits") ;
129 
130  CellIDDecoder<SimTrackerHit> idDec( trkHits ) ;
131 
132  for(int j=0;j<NHITS;j++) {
133 
134  //std::cout << " testing hit " << j << std::endl ;
135 
136  SimTrackerHit* trkHit = dynamic_cast<SimTrackerHit*>(trkHits->getElementAt(j)) ;
137 
138 
139  MYTEST( idDec(trkHit)[ LCTrackerCellID::subdet() ] , ILDDetID::FTD , " cellID(trkHit) == ( ILDDetID::FTD ) " ) ;
140  MYTEST( idDec(trkHit)[ LCTrackerCellID::layer() ] , j % 100 , " cellID(trkHit) == ( j % 100 ) " ) ;
141  MYTEST( idDec(trkHit)[ LCTrackerCellID::side() ] , ILDDetID::bwd , " cellID(trkHit) == ( ILDDetID::bwd ) " ) ;
142  MYTEST( idDec(trkHit)[ LCTrackerCellID::module() ] , j / 100 + 1 , " cellID(trkHit) == ( j / 100 + 1 ) " ) ;
143  MYTEST( idDec(trkHit)[ LCTrackerCellID::sensor() ] , j % 4 , " cellID(trkHit) == ( j % 4 ) " ) ;
144  MYTEST( idDec(trkHit)["daqChannel"] , j*8 , " cellID(trkHit) == ( j*8 ) " ) ;
145 
146  MYTEST( trkHit->getEDep() , i*j*117. , "EDep" ) ;
147 
148  MYTEST( trkHit->isProducedBySecondary() , bool((i+j)%3==0) , "Secondary" ) ;
149  MYTEST( trkHit->isOverlay() , bool((i+j)%5==0) , "Overlay" ) ;
150 
151  const double* pos = trkHit->getPosition() ;
152 
153  MYTEST( pos[0] , i , " pos[0] " ) ;
154  MYTEST( pos[1] , j , " pos[1] " ) ;
155  MYTEST( pos[2] , i*j , " pos[2] " ) ;
156 
157  }
158  }
159  lcRdr->close() ;
160 
161 
162  } catch( Exception &e ){
163  MYTEST.FAILED( e.what() );
164  }
165 
166  return 0;
167 }
168 
169 //=============================================================================
static string FILEN
Definition: tutil.h:7
static const string testname
void FAILED(const std::string &msg)
Definition: tutil.h:42
int main(int argc, char **argv)
Simple program that opens existing LCIO files and appends the records needed for direct access - if t...
static const int NEVENT
static const int NHITS
void LOG(const std::string &msg)
Definition: tutil.h:21