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

#include <PndGemFindClusters.h>

Inheritance diagram for PndGemFindClusters:

Public Member Functions

 PndGemFindClusters ()
 
 PndGemFindClusters (Int_t iVerbose)
 
 PndGemFindClusters (const char *name, Int_t iVerbose=0)
 
virtual ~PndGemFindClusters ()
 
virtual void Exec (Option_t *opt)
 

Private Member Functions

virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
Int_t SortDigis ()
 
Int_t CreateClusters ()
 
Int_t WriteClusters ()
 
Bool_t CompareDigiToClustersDigis (Int_t digiNumber)
 
void SortClusters ()
 
void SortCluster (Int_t clus)
 
void PrintClusters ()
 
void PrintCluster (Int_t clus)
 
void JoinTwoClusters (Int_t clus1, Int_t clus2)
 
void AddDigiToCluster (Int_t digiNr, Int_t clusNr)
 
virtual void Finish ()
 
 ClassDef (PndGemFindClusters, 1)
 

Private Attributes

PndGemDigiParfDigiPar
 
TClonesArray * fDigis
 
TClonesArray * fClusters
 
std::map< Double_t, Int_t > fTimeOrderedDigis
 
std::vector< DigiClusterfDigiClusters
 
Int_t fTNofEvents
 
Int_t fTNofDigis
 
Int_t fTNofClusters
 
TStopwatch fTimer
 

Detailed Description

Definition at line 51 of file PndGemFindClusters.h.

Constructor & Destructor Documentation

PndGemFindClusters::PndGemFindClusters ( )

Default constructor

Definition at line 50 of file PndGemFindClusters.cxx.

References fClusters, fDigiPar, fDigis, fTNofClusters, fTNofDigis, and fTNofEvents.

51  : FairTask("GEMFindClusters", 1) {
52  fDigiPar = NULL;
53  fDigis = NULL;
54  fClusters = NULL;
55 
56  fTNofEvents = 0;
57  fTNofDigis = 0;
58  fTNofClusters = 0;
59 }
PndGemDigiPar * fDigiPar
TClonesArray * fClusters
TClonesArray * fDigis
PndGemFindClusters::PndGemFindClusters ( Int_t  iVerbose)

Standard constructor

Definition at line 65 of file PndGemFindClusters.cxx.

References fClusters, fDigiPar, fDigis, fTNofClusters, fTNofDigis, and fTNofEvents.

66  : FairTask("GEMFindClusters", iVerbose) {
67  fDigiPar = NULL;
68  fDigis = NULL;
69  fClusters = NULL;
70 
71  fTNofEvents = 0;
72  fTNofDigis = 0;
73  fTNofClusters = 0;
74 }
PndGemDigiPar * fDigiPar
TClonesArray * fClusters
TClonesArray * fDigis
Int_t iVerbose
PndGemFindClusters::PndGemFindClusters ( const char *  name,
Int_t  iVerbose = 0 
)

Constructor with task name

Definition at line 80 of file PndGemFindClusters.cxx.

References fClusters, fDigiPar, fDigis, fTNofClusters, fTNofDigis, and fTNofEvents.

81  : FairTask(name, iVerbose) {
82  fDigiPar = NULL;
83  fDigis = NULL;
84  fClusters = NULL;
85 
86  fTNofEvents = 0;
87  fTNofDigis = 0;
88  fTNofClusters = 0;
89 }
PndGemDigiPar * fDigiPar
TClonesArray * fClusters
TClonesArray * fDigis
TString name
Int_t iVerbose
PndGemFindClusters::~PndGemFindClusters ( )
virtual

Destructor

Definition at line 95 of file PndGemFindClusters.cxx.

References fClusters.

95  {
96  if ( fClusters ) {
97  fClusters->Delete();
98  delete fClusters;
99  }
100 }
TClonesArray * fClusters

Member Function Documentation

void PndGemFindClusters::AddDigiToCluster ( Int_t  digiNr,
Int_t  clusNr 
)
private

Definition at line 948 of file PndGemFindClusters.cxx.

References digi, Double_t, fDigiClusters, fDigiPar, fDigis, fVerbose, PndGemDigi::GetChannelNr(), PndGemDigi::GetCharge(), PndGemDigi::GetDetectorId(), PndGemSensor::GetMeanChannel(), PndGemDigiPar::GetSensor(), PndGemDigi::GetSensorNr(), PndGemDigi::GetSide(), PndGemDigi::GetStationNr(), and sensor.

Referenced by CompareDigiToClustersDigis(), and JoinTwoClusters().

948  {
949  if ( fVerbose > 0 ) cout << "PndGemFindClusters::AddDigiToCluster()" << endl;
950 
951  PndGemDigi* digi = (PndGemDigi*) fDigis->At(digiNr);
952  Int_t stationNr = digi->GetStationNr();
953  Int_t sensorNr = digi->GetSensorNr();
954  Int_t iSide = digi->GetSide();
955 
956  PndGemSensor* sensor = fDigiPar->GetSensor(stationNr, sensorNr);
957 
958  if ( fDigiClusters[clusNr].digiNr.size() == 0 ) {
959  // cout << "NEW CLUSTER " << clusNr << " CREATED WITH AddDigiToCluster(" << digiNr << "," << clusNr << ") IN DETECTOR " << digi->GetDetectorId() << endl;
960  fDigiClusters[clusNr].detId = digi->GetDetectorId();
961  fDigiClusters[clusNr].cluPos = digi->GetChannelNr();
962  fDigiClusters[clusNr].cluADC = digi->GetCharge();
963  fDigiClusters[clusNr].cluTDC = digi->GetTimeStamp();
964  fDigiClusters[clusNr].cluPMn = digi->GetChannelNr();
965  fDigiClusters[clusNr].cluPMx = digi->GetChannelNr();
966  fDigiClusters[clusNr].cluMVl = digi->GetCharge();
967  fDigiClusters[clusNr].cluMPs = digi->GetChannelNr();
968  }
969  else {
970  fDigiClusters[clusNr].cluTDC = (fDigiClusters[clusNr].cluTDC*((Double_t)(fDigiClusters[clusNr].sigADC.size()-1.))+digi->GetTimeStamp())/((Double_t)(fDigiClusters[clusNr].sigADC.size()));
971  fDigiClusters[clusNr].cluPos = sensor->GetMeanChannel(iSide,fDigiClusters[clusNr].cluPos,fDigiClusters[clusNr].cluADC,digi->GetChannelNr(),digi->GetCharge());
972  if ( fDigiClusters[clusNr].cluPMn > digi->GetChannelNr() ) fDigiClusters[clusNr].cluPMn = digi->GetChannelNr();
973  if ( fDigiClusters[clusNr].cluPMx < digi->GetChannelNr() ) fDigiClusters[clusNr].cluPMx = digi->GetChannelNr();
974  fDigiClusters[clusNr].cluADC = fDigiClusters[clusNr].cluADC+digi->GetCharge();
975  if ( digi->GetCharge() > fDigiClusters[clusNr].cluMVl ) {
976  fDigiClusters[clusNr].cluMVl = digi->GetCharge();
977  fDigiClusters[clusNr].cluMPs = digi->GetChannelNr();
978  }
979  }
980 
981  fDigiClusters[clusNr].digiNr.push_back(digiNr);
982  fDigiClusters[clusNr].chanNr.push_back((Int_t)digi->GetChannelNr());
983  fDigiClusters[clusNr].sigADC.push_back(digi->GetCharge());
984 }
PndGemDigiPar * fDigiPar
int fVerbose
Definition: poormantracks.C:24
TClonesArray * digi
Int_t GetStationNr() const
Definition: PndGemDigi.h:84
TClonesArray * fDigis
TGeoVolume * sensor
Int_t GetDetectorId() const
Definition: PndGemDigi.h:79
Int_t GetSensorNr() const
Definition: PndGemDigi.h:86
PndGemSensor * GetSensor(Int_t stationNr, Int_t sensorNr)
Double_t
std::vector< DigiCluster > fDigiClusters
Double_t GetChannelNr() const
Definition: PndGemDigi.h:80
Double_t GetCharge() const
Definition: PndGemDigi.h:91
Double_t GetMeanChannel(Int_t iSide, Double_t chan1, Double_t weight1, Double_t chan2, Double_t weight2)
Int_t GetSide() const
Definition: PndGemDigi.h:88
PndGemFindClusters::ClassDef ( PndGemFindClusters  ,
 
)
private
Bool_t PndGemFindClusters::CompareDigiToClustersDigis ( Int_t  digiNumber)
private

Compare digi to clusters Compare digi to digis

Definition at line 877 of file PndGemFindClusters.cxx.

References AddDigiToCluster(), digi, Double_t, fDigiClusters, fDigiPar, fDigis, fVerbose, PndGemDigi::GetChannelNr(), PndGemDigi::GetDetectorId(), PndGemSensor::GetDistance2(), PndGemDigiPar::GetSensor(), PndGemDigi::GetSensorNr(), PndGemDigi::GetSide(), PndGemDigi::GetStationNr(), JoinTwoClusters(), and sensor.

Referenced by CreateClusters().

877  {
878  if ( fVerbose > 0 ) cout << "PndGemFindClusters::CompareDigiToClustersDigis" << endl;
879 
881  PndGemDigi* digi;
882  digi = (PndGemDigi*) fDigis->At(digiNumber);
883 
884  Int_t iDigiUsed = -1;
885  for ( size_t idc = 0 ; idc < fDigiClusters.size() ; idc++ ) {
886  if ( fDigiClusters[idc].cluPos < -0.5 ) continue;
887  if ( digi->GetDetectorId() != fDigiClusters[idc].detId ) continue;
888 
889  Int_t stationNr = digi->GetStationNr();
890  Int_t sensorNr = digi->GetSensorNr();
891  Int_t iSide = digi->GetSide();
892 
893  sensor = fDigiPar->GetSensor(stationNr, sensorNr);
894 
895  Double_t digiDigiDist=999999999, testDist;
896  for ( size_t id = 0 ; id < fDigiClusters[idc].digiNr.size() ; id++ ) {
897  testDist = sensor->GetDistance2(iSide,digi->GetChannelNr(),fDigiClusters[idc].chanNr[id]);
898  if ( testDist < digiDigiDist ) digiDigiDist=testDist;
899  }
900 
901  if ( digiDigiDist < 0. ) continue;
902  if ( digiDigiDist > 1 ) continue;//1 for next strip to accept
903 
904  if ( fVerbose > 1 )
905  cout << "digi " << digiNumber << " ( " << digi->GetDetectorId() << " / " << digi->GetChannelNr() << " ) could attach to cluster "
906  << idc << " ( " << fDigiClusters[idc].cluPos << " ) with distance " << digiDigiDist << " . " << fDigiClusters[idc].digiNr.size() << " digis in cluster" << flush;
907 
908  AddDigiToCluster(digiNumber,idc);
909 
910  if ( fVerbose > 1 )
911  cout << " >>> mean channel = " << fDigiClusters[idc].cluPos << endl;
912 
913  if ( iDigiUsed >= 0 ) {
914  //cout << "DIGI " << digiNumber << " MATCHES TO CLUSTER " << iDigiUsed << " and " << idc << endl;
915  JoinTwoClusters(idc,iDigiUsed);
916  }
917 
918  iDigiUsed = idc;
919  }
920  return (iDigiUsed==-1?kFALSE:kTRUE);
921 
922 }
PndGemDigiPar * fDigiPar
int fVerbose
Definition: poormantracks.C:24
TClonesArray * digi
Int_t GetStationNr() const
Definition: PndGemDigi.h:84
TClonesArray * fDigis
TGeoVolume * sensor
Int_t GetDetectorId() const
Definition: PndGemDigi.h:79
Int_t GetSensorNr() const
Definition: PndGemDigi.h:86
PndGemSensor * GetSensor(Int_t stationNr, Int_t sensorNr)
Double_t
std::vector< DigiCluster > fDigiClusters
void AddDigiToCluster(Int_t digiNr, Int_t clusNr)
Double_t GetChannelNr() const
Definition: PndGemDigi.h:80
void JoinTwoClusters(Int_t clus1, Int_t clus2)
Double_t GetDistance2(Int_t iSide, Double_t chan1, Double_t chan2)
Int_t GetSide() const
Definition: PndGemDigi.h:88
Int_t PndGemFindClusters::CreateClusters ( )
private

Find clusters Look for clusters

Definition at line 322 of file PndGemFindClusters.cxx.

References DigiCluster::chanNr, DigiCluster::cluADC, DigiCluster::cluMPs, DigiCluster::cluMVl, DigiCluster::cluPMn, DigiCluster::cluPMx, DigiCluster::cluPos, DigiCluster::cluTDC, CompareDigiToClustersDigis(), DigiCluster::detId, digi, DigiCluster::digiNr, fDigiClusters, fDigis, fVerbose, PndGemDigi::GetChannelNr(), PndGemDigi::GetCharge(), PndGemDigi::GetDetectorId(), and DigiCluster::sigADC.

Referenced by Exec().

322  {
323  if ( fVerbose > 0 ) {
324  cout << "-I- PndGemFindClusters::CreateClusters()" <<endl;
325  }
326 
327  fDigiClusters.clear();
328 
329  PndGemDigi* digi = NULL;
330  for ( Int_t idigi = 0 ; idigi < fDigis->GetEntries() ; idigi++ ) {
331  digi = (PndGemDigi*) fDigis->At(idigi);
332 
333  if ( fVerbose ) {
334  // cout << "digi " << idigi << " in " << digi->GetDetectorId() << ", there are already " << fDigiClusters.size() << " clusters" << endl;
335  // cout << idigi << " ---> " << digi->GetDetectorId() << " @ " << digi->GetChannelNr() << " @ " << digi->GetTimeStamp() << endl;
336  }
337  if ( CompareDigiToClustersDigis(idigi) == kTRUE ) continue;
338 
339  DigiCluster tempDC;
340  tempDC.detId = digi->GetDetectorId();
341  tempDC.digiNr.push_back(idigi);
342  tempDC.chanNr.push_back((Int_t)digi->GetChannelNr());
343  tempDC.sigADC.push_back(digi->GetCharge());
344  tempDC.cluTDC = digi->GetTimeStamp();
345  tempDC.cluADC = digi->GetCharge();
346  tempDC.cluPos = digi->GetChannelNr();
347  tempDC.cluPMn = digi->GetChannelNr();
348  tempDC.cluPMx = digi->GetChannelNr();
349  tempDC.cluMPs = digi->GetChannelNr();
350  tempDC.cluMVl = digi->GetCharge();
351  fDigiClusters.push_back(tempDC);
352 
353  }
354  return fDigiClusters.size();
355 }
std::vector< Int_t > digiNr
int fVerbose
Definition: poormantracks.C:24
Bool_t CompareDigiToClustersDigis(Int_t digiNumber)
TClonesArray * digi
std::vector< Int_t > chanNr
TClonesArray * fDigis
Int_t GetDetectorId() const
Definition: PndGemDigi.h:79
std::vector< DigiCluster > fDigiClusters
Double_t GetChannelNr() const
Definition: PndGemDigi.h:80
std::vector< Double_t > sigADC
Double_t GetCharge() const
Definition: PndGemDigi.h:91
void PndGemFindClusters::Exec ( Option_t *  opt)
virtual

Execution

Definition at line 106 of file PndGemFindClusters.cxx.

References CreateClusters(), fClusters, fDigiClusters, fTimer, fTNofClusters, fTNofEvents, fVerbose, PrintClusters(), SortClusters(), SortDigis(), and WriteClusters().

106  {
107  if ( fVerbose > 0 ) cout <<"-I- PndGemFindClusters::Exec()"<<endl;
108 
109  fTNofEvents++;
110 
111  fTimer.Start();
112  //Bool_t warn = kFALSE; //[R.K. 01/2017] unused variable?
113 
114  // Clear output array
115  fClusters->Clear();
116 
117  //Int_t nDigis = 0; //[R.K.03/2017] unused variable
118  Int_t nClusters = 0;
119 
120  // Sort digis with respect to time
121  SortDigis();//nDigis = //[R.K.03/2017] unused variable
122 
123  //FindClusters();
124 
125  CreateClusters();
126 
127  if ( fVerbose > 1 ) {
128  Int_t nofC2W = 0;
129  for ( size_t idc = 0 ; idc < fDigiClusters.size() ; idc++ )
130  if ( fDigiClusters[idc].cluPos > -0.5 )
131  nofC2W++;
132  cout << ">>>> created " << fDigiClusters.size() << " clusters, " << nofC2W << " to write" << endl;
133  }
134 
135  //Skip Analyze... the join is already done in the FindClusters()..
136  //AnalyzeClusters();
137 
138 // if ( fVerbose ) {
139 // Int_t nofC2W = 0;
140 // for ( Int_t idc = 0 ; idc < fDigiClusters.size() ; idc++ )
141 // if ( fDigiClusters[idc].cluPos > -0.5 )
142 // nofC2W++;
143 // cout << ">>>> analyzed " << fDigiClusters.size() << " clusters, " << nofC2W << " to write" << endl;
144 // }
145 
146  //ClearClusters();
147  //ClearClusters2();//sort and print only..
148  SortClusters();
149 
150  if ( fVerbose > 1 ) {
151  PrintClusters();
152  }
153 
154 // if ( fVerbose ) {
155 // Int_t nofC2W = 0;
156 // for ( Int_t idc = 0 ; idc < fDigiClusters.size() ; idc++ )
157 // if ( fDigiClusters[idc].cluPos > -0.5 )
158 // nofC2W++;
159 // cout << ">>>> cleared " << fDigiClusters.size() << " clusters, " << nofC2W << " to write" << endl;
160 // }
161 
162  nClusters = WriteClusters();
163 
164  fTimer.Stop();
165 
166  // cout << "CREATED " << nClusters << " CLUSTERS OUT OF " << nDigis << " DIGIS IN " << fTimer.RealTime() << "s" << endl;
167 
168  fTNofClusters += nClusters;
169 }
int fVerbose
Definition: poormantracks.C:24
TClonesArray * fClusters
std::vector< DigiCluster > fDigiClusters
void PndGemFindClusters::Finish ( )
privatevirtual

Finish at the end of each event

Definition at line 1063 of file PndGemFindClusters.cxx.

References Double_t, fClusters, fTNofClusters, fTNofDigis, and fTNofEvents.

1063  {
1064  if ( fClusters ) fClusters->Clear();
1065 
1066  cout << "-------------------- " << fName.Data() << " : Summary -----------------------" << endl;
1067  cout << " Events: " << setw(10) << fTNofEvents << endl;
1068  cout << " Digis: " << setw(10) << fTNofDigis << " ( " << (Double_t)fTNofDigis /((Double_t)fTNofEvents) << " per event )" << endl;
1069  cout << " Clusters: " << setw(10) << fTNofClusters << " ( " << (Double_t)fTNofClusters/((Double_t)fTNofEvents) << " per event )" << endl;
1070  cout << " --> ( " << (Double_t)fTNofDigis /((Double_t)fTNofClusters) << " digis per cluster )" << endl;
1071  cout << "---------------------------------------------------------------------" << endl;
1072 
1073 }
TClonesArray * fClusters
Double_t
InitStatus PndGemFindClusters::Init ( )
privatevirtual

Intialisation

Definition at line 193 of file PndGemFindClusters.cxx.

References fClusters, fDigiPar, fDigis, PndGemDigiPar::GetNStations(), PndGemSensor::GetOuterRadius(), PndGemStation::GetSensor(), PndGemDigiPar::GetStation(), and sensor.

193  {
194 
195  // Get input array
196  FairRootManager* ioman = FairRootManager::Instance();
197  if ( ! ioman ) Fatal("Init", "No FairRootManager");
198  fDigis = (TClonesArray*) ioman->GetObject("GEMDigi");
199 
200  // Register output array
201  //fClusters = new TClonesArray("PndGemDigi", 1000);
202  fClusters = new TClonesArray("PndGemCluster", 1000);
203  ioman->Register("GEMCluster", "Cluster in GEM", fClusters, kTRUE);
204 
205  cout << "-I- " << fName.Data() << "::Init(). There are " << fDigiPar->GetNStations() << " GEM stations." << endl;
206  cout << "-I- " << fName.Data() << "::Init(). Initialization succesfull." << endl;
207 
209  PndGemSensor* sensor = (PndGemSensor*)station->GetSensor(0);
210  cout << "sensor out rad is " << sensor->GetOuterRadius() << endl;
211 
212  return kSUCCESS;
213 }
PndGemDigiPar * fDigiPar
TClonesArray * fClusters
TClonesArray * fDigis
TGeoVolume * sensor
PndGemSensor * GetSensor(Int_t iSensor)
Definition: PndGemStation.h:63
Int_t GetNStations()
Definition: PndGemDigiPar.h:45
Double_t GetOuterRadius() const
Definition: PndGemSensor.h:103
PndGemStation * GetStation(Int_t iStation)
void PndGemFindClusters::JoinTwoClusters ( Int_t  clus1,
Int_t  clus2 
)
private

Definition at line 926 of file PndGemFindClusters.cxx.

References AddDigiToCluster(), fDigiClusters, and fVerbose.

Referenced by CompareDigiToClustersDigis().

926  {
927  if ( fVerbose > 0 ) cout << "PndGemFindClusters::JointTwoClusters()" << endl;
928 
929  //PndGemDigi* digi; //[R.K. 01/2017] unused variable?
930 
931  for ( size_t id2 = 0 ; id2 < fDigiClusters[clus2].digiNr.size() ; id2++ ) {
932  Int_t digi2 = fDigiClusters[clus2].digiNr[id2];
933  Int_t chan2 = fDigiClusters[clus2].chanNr[id2];
934 
935  Int_t addDigi = kTRUE;
936  for ( size_t id1 = 0 ; id1 < fDigiClusters[clus1].digiNr.size() ; id1++ ) {
937  if ( chan2 == fDigiClusters[clus1].chanNr[id1] ) addDigi = kFALSE;
938  }
939 
940  if ( addDigi ) AddDigiToCluster(digi2,clus1);
941  }
942  fDigiClusters[clus2].cluPos = -1.;
943 }
int fVerbose
Definition: poormantracks.C:24
std::vector< DigiCluster > fDigiClusters
void AddDigiToCluster(Int_t digiNr, Int_t clusNr)
void PndGemFindClusters::PrintCluster ( Int_t  clus)
private

Definition at line 1033 of file PndGemFindClusters.cxx.

References Double_t, fDigiClusters, and fVerbose.

Referenced by PrintClusters().

1033  {
1034  if ( fVerbose > 0 ) cout << "PndGemFindClusters::PrintCluster()" << endl;
1035 
1036  cout << "CLUSTER " << clus << " IN " << fDigiClusters[clus].detId << " AT " << fDigiClusters[clus].cluPos
1037  << " CREATED FROM " << fDigiClusters[clus].digiNr.size() << " DIGIS, FROM "
1038  << fDigiClusters[clus].cluPMn << " TO " << fDigiClusters[clus].cluPMx << " WITH MAXIMUM OF "
1039  << fDigiClusters[clus].cluMVl << " AT " << fDigiClusters[clus].cluMPs << endl;
1040 
1041  for ( Int_t icol = 0 ; icol < fDigiClusters[clus].cluPMx-fDigiClusters[clus].cluPMn+3 ; icol++ ) cout << "-" << flush; cout<<endl;
1042  for ( Double_t fh = 0.99 ; fh >= 0. ; fh-=0.247475 ) {
1043  cout << "|" << flush;
1044  Int_t thisDigi = 0;
1045  Double_t thisADC = 0.;
1046  for ( Int_t ichan = fDigiClusters[clus].cluPMn ; ichan <= fDigiClusters[clus].cluPMx ; ichan++ ) {
1047  thisADC = 0.;
1048  if ( fDigiClusters[clus].chanNr[thisDigi] == ichan ) {
1049  thisADC = fDigiClusters[clus].sigADC[thisDigi];
1050  thisDigi++;
1051  }
1052  if ( thisADC > fh*fDigiClusters[clus].cluMVl ) cout << "*" << flush;
1053  else cout << " " << flush;
1054  }
1055  cout << "|" << endl;
1056  }
1057  for ( Int_t icol = 0 ; icol < fDigiClusters[clus].cluPMx-fDigiClusters[clus].cluPMn+3 ; icol++ ) cout << "-" << flush; cout<<endl;
1058 }
int fVerbose
Definition: poormantracks.C:24
Double_t
std::vector< DigiCluster > fDigiClusters
void PndGemFindClusters::PrintClusters ( )
private

Definition at line 1021 of file PndGemFindClusters.cxx.

References fDigiClusters, fVerbose, and PrintCluster().

Referenced by Exec().

1021  {
1022  if ( fVerbose > 0 ) cout << "PndGemFindClusters::PrintClusters()" << endl;
1023 
1024  for ( size_t idc = 0 ; idc < fDigiClusters.size() ; idc++ ) {
1025  // if ( fDigiClusters[idc].detId != 568329024 ) continue;
1026  if ( fDigiClusters[idc].cluPos < -0.5 ) continue;
1027  PrintCluster(idc);
1028  }
1029 }
int fVerbose
Definition: poormantracks.C:24
void PrintCluster(Int_t clus)
std::vector< DigiCluster > fDigiClusters
InitStatus PndGemFindClusters::ReInit ( )
privatevirtual

Reinitialisation

Definition at line 220 of file PndGemFindClusters.cxx.

220  {
221 
222  // Create sectorwise digi sets
223  // MakeSets();
224 
225  return kSUCCESS;
226 }
void PndGemFindClusters::SetParContainers ( )
privatevirtual

Get parameter containers

Definition at line 175 of file PndGemFindClusters.cxx.

References fDigiPar, and run.

175  {
176 
177  // Get run and runtime database
178  FairRunAna* run = FairRunAna::Instance();
179  if ( ! run ) Fatal("SetParContainers", "No analysis run");
180 
181  FairRuntimeDb* db = run->GetRuntimeDb();
182  if ( ! db ) Fatal("SetParContainers", "No runtime database");
183 
184  // Get GEM digitisation parameter container
185  fDigiPar = (PndGemDigiPar*) db->getContainer("PndGemDetectors");
186 }
PndGemDigiPar * fDigiPar
Int_t run
Definition: autocutx.C:47
Digitization Parameter Class for GEM part.
Definition: PndGemDigiPar.h:31
void PndGemFindClusters::SortCluster ( Int_t  clus)
private

Definition at line 997 of file PndGemFindClusters.cxx.

References Double_t, fDigiClusters, and fVerbose.

Referenced by SortClusters().

997  {
998  if ( fVerbose > 0 ) cout << "PndGemFindClusters::SortCluster()" << endl;
999  Int_t channelPosition = 0;
1000 
1001  for ( Int_t ich = fDigiClusters[clus].cluPMn ; ich <= fDigiClusters[clus].cluPMx ; ich++ ) {
1002  for ( size_t idigi = channelPosition ; idigi < fDigiClusters[clus].digiNr.size() ; idigi++ ) {
1003  if ( fDigiClusters[clus].chanNr[idigi] == ich ) {
1004  Int_t dN = fDigiClusters[clus].digiNr[idigi];
1005  Int_t cN = fDigiClusters[clus].chanNr[idigi];
1006  Double_t sA = fDigiClusters[clus].sigADC[idigi];
1007  fDigiClusters[clus].digiNr[idigi] = fDigiClusters[clus].digiNr[channelPosition];
1008  fDigiClusters[clus].chanNr[idigi] = fDigiClusters[clus].chanNr[channelPosition];
1009  fDigiClusters[clus].sigADC[idigi] = fDigiClusters[clus].sigADC[channelPosition];
1010  fDigiClusters[clus].digiNr[channelPosition] = dN;
1011  fDigiClusters[clus].chanNr[channelPosition] = cN;
1012  fDigiClusters[clus].sigADC[channelPosition] = sA;
1013  channelPosition++;
1014  }
1015  }
1016  }
1017 }
int fVerbose
Definition: poormantracks.C:24
Double_t
std::vector< DigiCluster > fDigiClusters
void PndGemFindClusters::SortClusters ( )
private

Analyze clusters Clear clusters

Definition at line 988 of file PndGemFindClusters.cxx.

References fDigiClusters, fVerbose, and SortCluster().

Referenced by Exec().

988  {
989  if ( fVerbose > 0 ) cout << "PndGemFindClusters::SortClusters()" << endl;
990  for ( size_t idc = 0 ; idc < fDigiClusters.size() ; idc++ ) {
991  SortCluster(idc);
992  }
993 }
int fVerbose
Definition: poormantracks.C:24
void SortCluster(Int_t clus)
std::vector< DigiCluster > fDigiClusters
Int_t PndGemFindClusters::SortDigis ( )
private

Sort digis timewise

Definition at line 231 of file PndGemFindClusters.cxx.

References Double_t, fDigis, fTimeOrderedDigis, fTNofDigis, and fVerbose.

Referenced by Exec().

231  {
232  if ( fVerbose ) cout << "-I- PndGemFindClusters::SortDigis()" <<endl;
233 
234  // Check input array
235  if ( ! fDigis ) {
236  cout << "-E- " << fName << "::SortDigis: No input array!" << endl;
237  return -1;
238  }
239  TRandom2* randNumber = new TRandom2();
240 
241  //PndGemDigi* digi; //[R.K. 01/2017] unused variable?
242  Int_t nDigis = fDigis->GetEntriesFast();
243 
244  fTimeOrderedDigis.clear();
245 
246  fTNofDigis += nDigis;
247 
248  for (Int_t iDigi=0; iDigi<nDigis; iDigi++) {
249  //digi = (PndGemDigi*) fDigis->At(iDigi); //[R.K. 01/2017] unused variable?
250  // Double_t time = digi->GetTimeStamp();
251  Double_t time = randNumber->Gaus(0.,5.);
252 
253  fTimeOrderedDigis[time] = iDigi;
254  }
255 
256  // cout << "have " << nDigis << " digis ( " << fTimeOrderedDigis.size() << " ) " << endl;
257 
258  map<Double_t, Int_t>::iterator mapIt;
259  Int_t nd = 0;
260  for (mapIt=fTimeOrderedDigis.begin(); mapIt!=fTimeOrderedDigis.end(); mapIt++) {
261  // cout << mapIt->first << " --> " << mapIt->second << endl;
262  nd++;
263  }
264  // cout << " ---> " << nd << endl;
265 
266  return nDigis;
267 }
int fVerbose
Definition: poormantracks.C:24
std::map< Double_t, Int_t > fTimeOrderedDigis
TClonesArray * fDigis
Double_t
Int_t PndGemFindClusters::WriteClusters ( )
private

Write clusters

Definition at line 790 of file PndGemFindClusters.cxx.

References fDigiClusters, and fVerbose.

Referenced by Exec().

790  {
791  if ( fVerbose > 0 ) cout << "-I- PndGemFindClusters::WriteClusters()" <<endl;
792 
793  Int_t nClusters = 0;
794  //PndGemDigi* digi; //[R.K. 01/2017] unused variable?
795 
796  for ( size_t idc = 0 ; idc < fDigiClusters.size() ; idc++ ) {
797 
798  // if ( fDigiClusters[idc].cluADC < 1. ) continue;
799  if ( fDigiClusters[idc].cluPos < -0.5 ) continue;
800 
801  if ( fVerbose > 1 ) {
802  cout << "cluster at " << fDigiClusters[idc].cluPos << ", h = " << fDigiClusters[idc].cluADC << ", for detId = " << fDigiClusters[idc].detId << " from " << fDigiClusters[idc].digiNr.size() << " digis (from " << fDigiClusters[idc].cluPMn << " to " << fDigiClusters[idc].cluPMx << "): " << endl;
803 
804  for ( size_t idigi = 0 ; idigi < fDigiClusters[idc].digiNr.size() ; idigi++ ) {
805  cout <<" digiNr=" << fDigiClusters[idc].digiNr[idigi]
806  <<" chanNr=" << fDigiClusters[idc].chanNr[idigi]
807  <<" sigADC=" << fDigiClusters[idc].sigADC[idigi] << endl;
808  }
809  }
810 
811  new ((*fClusters)[nClusters]) PndGemCluster(fDigiClusters[idc].detId,
812  fDigiClusters[idc].cluPos,
813  fDigiClusters[idc].cluPMn,
814  fDigiClusters[idc].cluPMx,
815  fDigiClusters[idc].cluADC,
816  fDigiClusters[idc].cluTDC,
817  fDigiClusters[idc].digiNr);
818 
819  if ( fVerbose > 1 )
820  cout << "creating cluster at " << fDigiClusters[idc].cluPos << " with height of " << fDigiClusters[idc].cluADC << " /// compare maximum at " << fDigiClusters[idc].cluMPs << endl;
821  nClusters++;
822  }
823  return nClusters;
824 }
int fVerbose
Definition: poormantracks.C:24
std::vector< DigiCluster > fDigiClusters

Member Data Documentation

TClonesArray* PndGemFindClusters::fClusters
private

Input array of PndGemDigi

Definition at line 82 of file PndGemFindClusters.h.

Referenced by Exec(), Finish(), Init(), PndGemFindClusters(), and ~PndGemFindClusters().

std::vector<DigiCluster> PndGemFindClusters::fDigiClusters
private
PndGemDigiPar* PndGemFindClusters::fDigiPar
private
TClonesArray* PndGemFindClusters::fDigis
private

Digitisation parameters

Definition at line 81 of file PndGemFindClusters.h.

Referenced by AddDigiToCluster(), CompareDigiToClustersDigis(), CreateClusters(), Init(), PndGemFindClusters(), and SortDigis().

std::map<Double_t, Int_t> PndGemFindClusters::fTimeOrderedDigis
private

Output array of PndGemDigi

Definition at line 84 of file PndGemFindClusters.h.

Referenced by SortDigis().

TStopwatch PndGemFindClusters::fTimer
private

Definition at line 91 of file PndGemFindClusters.h.

Referenced by Exec().

Int_t PndGemFindClusters::fTNofClusters
private

Definition at line 89 of file PndGemFindClusters.h.

Referenced by Exec(), Finish(), and PndGemFindClusters().

Int_t PndGemFindClusters::fTNofDigis
private

Definition at line 88 of file PndGemFindClusters.h.

Referenced by Finish(), PndGemFindClusters(), and SortDigis().

Int_t PndGemFindClusters::fTNofEvents
private

Definition at line 87 of file PndGemFindClusters.h.

Referenced by Exec(), Finish(), and PndGemFindClusters().


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