FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndMvdConvertApv Class Reference

Convert meassured data into PndSdsDigiStrip. More...

#include <PndMvdConvertApv.h>

Public Member Functions

 PndMvdConvertApv ()
 
 PndMvdConvertApv (const TString &CalibFileName, const TString &HitFileName)
 
virtual ~PndMvdConvertApv ()
 
 PndMvdConvertApv (const PndMvdConvertApv &)=delete
 
PndMvdConvertApvoperator= (const PndMvdConvertApv &)=delete
 
long int GetNofEvents ()
 
std::vector< PndSdsDigiStripReadAll ()
 
std::vector< PndSdsDigiStripReadNext ()
 
Bool_t Init ()
 
void SetFakePair (Int_t TopModuleID, Int_t BottomModuleID)
 

Private Member Functions

void LoadCalibration (TString CalibFileName)
 
std::vector< PndSdsDigiStripCalc (std::vector< PndMvdApvHit > hitlist)
 
 ClassDef (PndMvdConvertApv, 1)
 Gives Access to the Path info of a hit. More...
 

Private Attributes

std::map< Int_t, std::map
< Int_t, double > > 
fCalibPars
 Calib Parameters: <moduleID<FE<channel> > > More...
 
long int fNofEvents
 stored last eventID of the hitfile More...
 
long int fEvent
 stored current eventID More...
 
long int fLastEvent
 last eventID More...
 
bool fNoCalib
 true if calib was succesfull loaded More...
 
TString fHitFileName
 name of hitfile More...
 
std::ifstream fDataFile
 fstream of hitfile More...
 
std::vector< PndMvdApvHitfhitlist
 stored readed events More...
 
Int_t fTopModuleID
 stored fake module top side More...
 
Int_t fBottomModuleID
 stored fake module bottom side More...
 
Bool_t fFake
 knows if fake is allowed More...
 
TFile * f
 
TTree * t
 
PndMvdTsEventtsEv
 
TClonesArray * arr
 
PndGeoHandlingfGeoH
 

Detailed Description

Convert meassured data into PndSdsDigiStrip.

PndMvdConvertApv.h

Author
L.Ackermann lars..nosp@m.acke.nosp@m.rmann.nosp@m.@phy.nosp@m.sik.t.nosp@m.u-dr.nosp@m.esden.nosp@m..de

Load calibration parameter and keeps strean to meassured data. Event wise read in and converting to PndSdsDigiStrip

Author
Lars Ackermann
Date
11.03.2009

Definition at line 39 of file PndMvdConvertApv.h.

Constructor & Destructor Documentation

PndMvdConvertApv::PndMvdConvertApv ( )
inline

default constructor

Definition at line 43 of file PndMvdConvertApv.h.

43  :
44  fCalibPars(),
45  fNofEvents(0),
46  fEvent(-1),
47  fLastEvent(0),
48  fNoCalib(kFALSE),
49  fHitFileName(""),
50  fDataFile(),
51  fhitlist(),
52  fTopModuleID(0),
53  fBottomModuleID(0),
54  fFake(kFALSE),
55  f(NULL),
56  t(NULL),
57  tsEv(NULL),
58  arr(NULL),
59  fGeoH(NULL)
60  {};
Int_t fTopModuleID
stored fake module top side
std::map< Int_t, std::map< Int_t, double > > fCalibPars
Calib Parameters: &lt;moduleID&lt;FE&lt;channel&gt; &gt; &gt;
Bool_t fFake
knows if fake is allowed
long int fEvent
stored current eventID
TString fHitFileName
name of hitfile
long int fNofEvents
stored last eventID of the hitfile
PndMvdTsEvent * tsEv
long int fLastEvent
last eventID
std::ifstream fDataFile
fstream of hitfile
Int_t fBottomModuleID
stored fake module bottom side
std::vector< PndMvdApvHit > fhitlist
stored readed events
TClonesArray * arr
bool fNoCalib
true if calib was succesfull loaded
PndGeoHandling * fGeoH
PndMvdConvertApv::PndMvdConvertApv ( const TString CalibFileName,
const TString HitFileName 
)

main constructor, call all function to be ready for converting hits from hitfile

Parameters
CalibFileNamename of file where the calibration is stored
HitFileNamename of hitfile

Definition at line 20 of file PndMvdConvertApv.cxx.

References arr, f, fGeoH, fNofEvents, PndGeoHandling::Instance(), LoadCalibration(), t, and tsEv.

20  :
21  fCalibPars(),
22  fNofEvents(0),
23  fEvent(-1),
24  fLastEvent(0),
25  fNoCalib(kFALSE),
26  fHitFileName(""),
27  fDataFile(),
28  fhitlist(),
29  fTopModuleID(0),
30  fBottomModuleID(0),
31  fFake(kFALSE),
32  f(NULL),
33  t(NULL),
34  tsEv(NULL),
35  arr(NULL),
36  fGeoH(NULL)
37 {
38 
39  f = new TFile(HitFileName);
41  t = (TTree*) f->Get("T");
42  tsEv = new PndMvdTsEvent();
43  arr = new TClonesArray("PndMvdSiHit");
44 
45  cout << "---------------------------------------" << endl;
46  cout << "Number of events: " << t->GetEntries() << endl;
47  cout << "---------------------------------------" << endl;
48 
49  t->SetBranchAddress("events",&tsEv);
50  cout << t->GetEntries() << " events in File" << endl;
51 
52  //LoadCalibration(CalibFileName, fes);
54 
55  fNofEvents=t->GetEntries();
56 
57  cout<<"** end of PndMvdConvertApv::PndMvdConvertApv(const TString& , const TString&) **"<<endl;
58 }
Int_t fTopModuleID
stored fake module top side
std::map< Int_t, std::map< Int_t, double > > fCalibPars
Calib Parameters: &lt;moduleID&lt;FE&lt;channel&gt; &gt; &gt;
Bool_t fFake
knows if fake is allowed
long int fEvent
stored current eventID
TString fHitFileName
name of hitfile
TString HitFileName
long int fNofEvents
stored last eventID of the hitfile
PndMvdTsEvent * tsEv
static PndGeoHandling * Instance()
long int fLastEvent
last eventID
std::ifstream fDataFile
fstream of hitfile
Int_t fBottomModuleID
stored fake module bottom side
std::vector< PndMvdApvHit > fhitlist
stored readed events
TClonesArray * arr
void LoadCalibration(TString CalibFileName)
bool fNoCalib
true if calib was succesfull loaded
PndGeoHandling * fGeoH
TString CalibFileName
virtual PndMvdConvertApv::~PndMvdConvertApv ( )
inlinevirtual

Destructor

Definition at line 70 of file PndMvdConvertApv.h.

References fDataFile.

71  { fDataFile.close(); };
std::ifstream fDataFile
fstream of hitfile
PndMvdConvertApv::PndMvdConvertApv ( const PndMvdConvertApv )
delete

Member Function Documentation

std::vector< PndSdsDigiStrip > PndMvdConvertApv::Calc ( std::vector< PndMvdApvHit hitlist)
private

Convert the readed hit and store them in the PndSdsDigiStrip

Parameters
hitlistvector of PndMvdApvHit hits of one event from hitfile
Returns
list of converted PndSdsDigiStrip

Definition at line 114 of file PndMvdConvertApv.cxx.

References Double_t, fCalibPars, fFake, fNoCalib, and kMVDHitsStrip.

Referenced by ReadNext().

115 {
116  std::vector<PndSdsDigiStrip> result;
117  for(UInt_t hitnumber=0;hitnumber<hitlist.size();hitnumber++)
118  {
119  Double_t q=0.;
120  if(fNoCalib)
121  {
122  q=1.*hitlist[hitnumber].GetADC(); // no calib adc -> e !!!
123  }else{
124  if (fCalibPars[hitlist[hitnumber].GetModuleID()].size())
125  {
126  if (fCalibPars[hitlist[hitnumber].GetModuleID()][hitlist[hitnumber].GetChannel()])
127  {
128  q=fCalibPars[hitlist[hitnumber].GetModuleID()][hitlist[hitnumber].GetChannel()]*(hitlist[hitnumber].GetADC())*1000.; // in electrons
129  }
130  }
131  }
132 
133  //TString detPath="Module";
134  //Int_t modId=-1; //[R.K.02/2017] Unused variable?
135  //if(fFake) //[R.K.02/2017] Unused variable?
136  //{ //[R.K.02/2017] Unused variable?
137  //if(fTopModuleID==hitlist[hitnumber].GetModuleID() || fBottomModuleID==hitlist[hitnumber].GetModuleID())
138  //modId = 1; //[R.K.02/2017] Unused variable?
139  //if(fBottomModuleID==hitlist[hitnumber].GetModuleID()) hitlist[hitnumber].SetFeID(hitlist[hitnumber].GetFeID()+3); //[R.K.02/2017] Unused variable?
140  //}else{ //[R.K.02/2017] Unused variable?
141  //modId = hitlist[hitnumber].GetModuleID(); //[R.K.02/2017] Unused variable?
142  //} //[R.K.02/2017] Unused variable?
143  /* detPath+=modId;
144  detPath+="Rect";
145  // std::cout<<detPath.Data()<<" | "<<modId<<std::endl;
146  //TGeoVolume* Vol=gGeoManager->FindVolumeFast(detPath);
147  TGeoVolume* Vol=0;
148  if(Vol!=0) {
149  // std::cout<<Vol->GetName()<<std::endl;
150  // Vol->GetNode(-1)->cd();
151  detPath="/SiliconTestStation_1/DummysensorAss_0/";
152  detPath+=Vol->GetName();
153  detPath+="_0";
154  gGeoManager->cd(detPath.Data());
155  }
156  else {
157  // std::cout<<" -E- PndMvdConvertApv::Calc(): "<<detPath.Data()<<" does not exist"<<std::endl;
158  detPath+="_nonexistent";
159  }*/
160  // std::cout<<detPath.Data()<<std::endl;
161  // std::cout<<gGeoManager->GetPath()<<std::endl;
162  // if (0==Vol) std::cout<<"0";
163  // else std::cout<<"1";
164  // std::cout<<Vol->GetName();
165 
166  // detPath = Vol->GetName();
167  // std::cout << "write Digi with "<< detPath.Data()<<" ( "<<fGeoH->GetID(detPath)<<" )"<<std::endl;
168  if(fFake)
169  {
170  PndSdsDigiStrip DigiHit(hitlist[hitnumber].GetEventID(), // index
171  kMVDHitsStrip, // panda detID
172  -1, // No SensorID (from geopath)
173  hitlist[hitnumber].GetFeID(), // fe
174  hitlist[hitnumber].GetChannel(), // chan
175  q/*/1000/1000*/, // charge
176  hitlist[hitnumber].GetTimestamp()// timestamp
177  );
178  DigiHit.SetTimeStampError(hitlist[hitnumber].GetTriggerID());
179 
180  // cout << "Ev. " << hitlist[hitnumber].GetEventID() << ", FE: " << hitlist[hitnumber].GetFeID() << ", ch: " << hitlist[hitnumber].GetChannel() << ", FAKE" << endl;
181 
182  result.push_back(DigiHit);
183  }else{
184  PndSdsDigiStrip DigiHit(hitlist[hitnumber].GetEventID(),
186  hitlist[hitnumber].GetModuleID(),
187  hitlist[hitnumber].GetFeID(),
188  hitlist[hitnumber].GetChannel(),
189  q/*/1000/1000*/,
190  hitlist[hitnumber].GetTimestamp()
191  );
192  DigiHit.SetTimeStampError(hitlist[hitnumber].GetTriggerID());
193 
194 
195  // cout << "Ev. " << hitlist[hitnumber].GetEventID() << ", FE: " << hitlist[hitnumber].GetChannel()/128 << ", ch: " << hitlist[hitnumber].GetChannel()%128 << ", sens" << hitlist[hitnumber].GetModuleID() << endl;
196  result.push_back(DigiHit);
197  }
198  }
199  return result;
200 }
std::map< Int_t, std::map< Int_t, double > > fCalibPars
Calib Parameters: &lt;moduleID&lt;FE&lt;channel&gt; &gt; &gt;
Class for digitised strip hits.
Bool_t fFake
knows if fake is allowed
Double_t
bool fNoCalib
true if calib was succesfull loaded
PndMvdConvertApv::ClassDef ( PndMvdConvertApv  ,
 
)
private

Gives Access to the Path info of a hit.

long int PndMvdConvertApv::GetNofEvents ( )
Returns
long integer of last eventID from hitfile

Definition at line 204 of file PndMvdConvertApv.cxx.

References fNofEvents.

205 {
206  return fNofEvents;
207 }
long int fNofEvents
stored last eventID of the hitfile
Bool_t PndMvdConvertApv::Init ( )

Initialize global geometry

Returns
success

Definition at line 60 of file PndMvdConvertApv.cxx.

Referenced by PndMvdConvertApvTask::Init().

61 {
62  return kTRUE;
63 }
void PndMvdConvertApv::LoadCalibration ( TString  CalibFileName)
private

Read the calibration of the modules

Parameters
CalibFileNameFilename where the calibration of the modules are listed
modulesvector of moduleIDs of modules in the hitfile
Returns
void load calibration from database

Definition at line 69 of file PndMvdConvertApv.cxx.

References c, fCalibPars, and fNoCalib.

Referenced by PndMvdConvertApv().

70 {
71 
72  std::ifstream calibfile(CalibFileName);
73  if(!calibfile)
74  {
75  cout<<"Calibration file not found"<<endl;
76  fNoCalib=true;
77  return;
78  }
79  while (!calibfile.eof()) // read data
80  {
81  char c;
82  calibfile>>c;
83  if (calibfile.eof()) break;
84  if (!isdigit(c))
85  {
86  char str[256];
87  calibfile.getline(str,256);
88  continue;
89  }
90  calibfile.putback(c);
91 
92  int boxID, channel;
93  double value;
94 
95  calibfile >> boxID >> channel >> value;
96 
97  //for(unsigned int vec=0;vec<fes.size();vec++)
98  //{
99  //if(feID==fes[vec])
100  //{
101  fCalibPars[boxID][channel]=value;
102  //}
103  //}
104  }
105  calibfile.close();
106  fNoCalib=false;
107  cout<<"Calibration succesfully read"<<endl;
108  return;
109 
110 }
std::map< Int_t, std::map< Int_t, double > > fCalibPars
Calib Parameters: &lt;moduleID&lt;FE&lt;channel&gt; &gt; &gt;
bool fNoCalib
true if calib was succesfull loaded
TString CalibFileName
PndMvdConvertApv& PndMvdConvertApv::operator= ( const PndMvdConvertApv )
delete
std::vector< PndSdsDigiStrip > PndMvdConvertApv::ReadAll ( )

read all events from hitfile

Returns
vector of PndSdsDigiStrip of the event

Definition at line 279 of file PndMvdConvertApv.cxx.

References fEvent, fNofEvents, i, and ReadNext().

280 {
281  std::vector<PndSdsDigiStrip> result;
282  while(fEvent!=fNofEvents)
283  {
284  std::vector<PndSdsDigiStrip> dummy=ReadNext();
285  for(unsigned int i=0;i<dummy.size();++i) result.push_back(dummy[i]);
286  }
287  return result;
288 }
Int_t i
Definition: run_full.C:25
long int fEvent
stored current eventID
std::vector< PndSdsDigiStrip > ReadNext()
long int fNofEvents
stored last eventID of the hitfile
std::vector< PndSdsDigiStrip > PndMvdConvertApv::ReadNext ( )

read the next event from hitfile

Returns
vector of PndSdsDigiStrip of the event

Definition at line 211 of file PndMvdConvertApv.cxx.

References arr, Calc(), ev, PndMvdSiHit::fAdc, PndMvdSiHit::fBox, PndMvdSiHit::fChannel, fe, fEvent, fhitlist, fLastEvent, fNofEvents, PndMvdSiHit::fNumFrames, PndMvdTsEvent::GetEventId(), PndMvdTsEvent::GetExtClockTimeStamp(), PndMvdTsEvent::GetSiHitList(), hit, t, and tsEv.

Referenced by PndMvdConvertApvTask::Exec(), for(), and ReadAll().

212 {
213 
214 
215  std::vector<PndSdsDigiStrip> digiList;
216  //bool work=true;
217 
218  //int triggID=0;
219  int fe=0;
220  //int ts=0;
221  int frame=0;
222  int ch=0;
223  int l=0;
224  int moduleID=0;
225  double q=0.;
226  long int ev=0;
227  UInt_t ClockReset = 0.;
228  ULong64_t ClockCounts = 0.;
229 
230  if (fEvent >= -1 && fEvent <= fNofEvents)
231  {
232 
233  t->GetEvent(fEvent);
234 
235  ClockCounts = tsEv->GetExtClockTimeStamp(ClockReset);
236  // Storing ClockReset in fTriggerID
237  // Storing ClockCounts in fTimestamp
238 
239  arr = tsEv->GetSiHitList();
240  fhitlist.clear();
241 
242  for (Int_t kk = 0 ; kk < arr->GetEntries() ; kk++)
243  {
244 
245  PndMvdSiHit *hit = (PndMvdSiHit*) arr->At(kk);
246  ev = tsEv->GetEventId();
247  //fe = (Int_t) (hit->fChannel)/128;
248  //ch = (Int_t) (hit->fChannel)%128;
249  ch = (Int_t) (hit->fChannel);
250  q = hit->fAdc;
251  l = hit->fNumFrames;
252  moduleID = hit->fBox;
253 
254  //if(fhitlist.size()>20) fhitlist.clear();
255  // digiList = Calc(fhitlist);
256 
257  fLastEvent=ev;
258  //PndMvdApvHit Hit(ev, moduleID, fe, triggID, ts, frame, ch, q, l);
259  PndMvdApvHit Hit(ev, moduleID, fe, ClockReset, ClockCounts, frame, ch, q, l);
260  //cout << "Ev. " << fEvent << ", trigg: " << triggID << ", sens: " << moduleID << ", ch: " << ch << ", fe: " << (Int_t) (hit->fChannel)/128 << ", channel: " << (Int_t) (hit->fChannel)%128 << endl;
261  fhitlist.push_back(Hit);
262  }
263 
264  //if(fhitlist.size()>20) fhitlist.clear();
265  digiList = Calc(fhitlist);
266  fEvent++;
267 
268  }
269 
270  if(!(fEvent%10000)) cout<<"[ "<<(fEvent*100)/fNofEvents<<" %] "<<fEvent<<" events converted ..."<<endl;
271  if(fEvent==fNofEvents) cout<<"[100 %] "<<fEvent<<" events converted"<<endl;
272 
273  return digiList;
274 
275 }
ULong64_t GetExtClockTimeStamp(UInt_t &resetCount)
Definition: PndMvdTsEvent.h:38
UInt_t GetEventId()
Definition: PndMvdTsEvent.h:29
std::vector< PndSdsDigiStrip > Calc(std::vector< PndMvdApvHit > hitlist)
int ev
UShort_t fBox
Definition: PndMvdSiHit.h:12
long int fEvent
stored current eventID
UShort_t fChannel
Definition: PndMvdSiHit.h:13
Class to store data of Apv-Sensors.
Definition: PndMvdApvHit.h:18
long int fNofEvents
stored last eventID of the hitfile
PndMvdTsEvent * tsEv
long int fLastEvent
last eventID
std::vector< PndMvdApvHit > fhitlist
stored readed events
Int_t fAdc
Definition: PndMvdSiHit.h:14
TClonesArray * arr
int fe
Definition: anaLmdDigi.C:67
TClonesArray * GetSiHitList()
Definition: PndMvdTsEvent.h:30
PndSdsMCPoint * hit
Definition: anasim.C:70
UShort_t fNumFrames
Definition: PndMvdSiHit.h:15
void PndMvdConvertApv::SetFakePair ( Int_t  TopModuleID,
Int_t  BottomModuleID 
)

Set two moduleIDs to merge them as one double sided sensor

Parameters
TopModuleIDmoduleID for top side of fake sensor
BottomModuleIDmoduleID for bottom side of fake sensor
Returns
void

Definition at line 290 of file PndMvdConvertApv.cxx.

References fBottomModuleID, fFake, and fTopModuleID.

291 {
292  fFake=true;
293  fTopModuleID=TopModuleID;
294  fBottomModuleID=BottomModuleID;
295  return;
296 }
Int_t fTopModuleID
stored fake module top side
Bool_t fFake
knows if fake is allowed
Int_t fBottomModuleID
stored fake module bottom side

Member Data Documentation

TClonesArray* PndMvdConvertApv::arr
private

Definition at line 184 of file PndMvdConvertApv.h.

Referenced by PndMvdConvertApv(), and ReadNext().

TFile* PndMvdConvertApv::f
private

Definition at line 180 of file PndMvdConvertApv.h.

Referenced by PndMvdConvertApv().

Int_t PndMvdConvertApv::fBottomModuleID
private

stored fake module bottom side

Definition at line 175 of file PndMvdConvertApv.h.

Referenced by SetFakePair().

std::map<Int_t, std::map<Int_t, double> > PndMvdConvertApv::fCalibPars
private

Calib Parameters: <moduleID<FE<channel> > >

Definition at line 148 of file PndMvdConvertApv.h.

Referenced by Calc(), and LoadCalibration().

std::ifstream PndMvdConvertApv::fDataFile
private

fstream of hitfile

Definition at line 166 of file PndMvdConvertApv.h.

Referenced by ~PndMvdConvertApv().

long int PndMvdConvertApv::fEvent
private

stored current eventID

Definition at line 154 of file PndMvdConvertApv.h.

Referenced by ReadAll(), and ReadNext().

Bool_t PndMvdConvertApv::fFake
private

knows if fake is allowed

Definition at line 178 of file PndMvdConvertApv.h.

Referenced by Calc(), and SetFakePair().

PndGeoHandling* PndMvdConvertApv::fGeoH
private

Definition at line 186 of file PndMvdConvertApv.h.

Referenced by PndMvdConvertApv().

TString PndMvdConvertApv::fHitFileName
private

name of hitfile

Definition at line 163 of file PndMvdConvertApv.h.

std::vector<PndMvdApvHit> PndMvdConvertApv::fhitlist
private

stored readed events

Definition at line 169 of file PndMvdConvertApv.h.

Referenced by ReadNext().

long int PndMvdConvertApv::fLastEvent
private

last eventID

Definition at line 157 of file PndMvdConvertApv.h.

Referenced by ReadNext().

bool PndMvdConvertApv::fNoCalib
private

true if calib was succesfull loaded

Definition at line 160 of file PndMvdConvertApv.h.

Referenced by Calc(), and LoadCalibration().

long int PndMvdConvertApv::fNofEvents
private

stored last eventID of the hitfile

Definition at line 151 of file PndMvdConvertApv.h.

Referenced by GetNofEvents(), PndMvdConvertApv(), ReadAll(), and ReadNext().

Int_t PndMvdConvertApv::fTopModuleID
private

stored fake module top side

Definition at line 172 of file PndMvdConvertApv.h.

Referenced by SetFakePair().

TTree* PndMvdConvertApv::t
private

Definition at line 181 of file PndMvdConvertApv.h.

Referenced by PndMvdConvertApv(), and ReadNext().

PndMvdTsEvent* PndMvdConvertApv::tsEv
private

Definition at line 182 of file PndMvdConvertApv.h.

Referenced by PndMvdConvertApv(), and ReadNext().


The documentation for this class was generated from the following files: