SND@LHC Software
Loading...
Searching...
No Matches
ecalModule Class Reference

#include <ecalModule.h>

Inheritance diagram for ecalModule:
Collaboration diagram for ecalModule:

Public Member Functions

 ecalModule (char type=1, Int_t cellnumber=-1, Float_t x1=0, Float_t y1=0, Float_t x2=0, Float_t y2=0, Int_t mc=0, Float_t energy=0)
 
ecalCellLocate (Int_t x, Int_t y) const
 
ecalCellAt (Int_t x, Int_t y) const
 
ecalCellFindCell (Float_t x, Float_t y) const
 
void AddEnergy (Float_t x, Float_t y, Float_t energy)
 
Float_t GetEnergy (Float_t x, Float_t y) const
 
void ResetModule ()
 
Float_t GetDX () const
 
Float_t GetDY () const
 
std::vector< ecalCell * > GetCells () const
 
std::list< ecalCell * > GetCellsX (Float_t x) const
 
std::list< ecalCell * > GetCellsY (Float_t y) const
 
- Public Member Functions inherited from ecalCell
 ecalCell (Int_t cellnumber, Float_t x1=0, Float_t y1=0, Float_t x2=0, Float_t y2=0, Char_t type=0, Float_t energy=0)
 
Bool_t IsInside (Float_t x, Float_t y)
 
Char_t GetType () const
 
Float_t X1 () const
 
Float_t Y1 () const
 
Float_t X2 () const
 
Float_t Y2 () const
 
Float_t GetX1 () const
 
Float_t GetY1 () const
 
Float_t GetX2 () const
 
Float_t GetY2 () const
 
Float_t GetCenterX () const
 
Float_t GetCenterY () const
 
Short_t ADC () const
 
Short_t GetADC () const
 
Int_t GetCellNumber () const
 
Float_t GetEnergy () const
 
Float_t GetTime () const
 
void SetTime (Float_t time)
 
void GetNeighborsList (std::list< ecalCell * > &neib) const
 
void SetNeighborsList (std::list< ecalCell * > &neib)
 
void Get5x5Cluster (std::list< ecalCell * > &cls)
 
void SetEnergy (Float_t energy)
 
void SetADC (Short_t adc)
 
void ResetEnergyFast ()
 
void AddEnergy (Float_t energy)
 
void GetClusterEnergy (Float_t &EcalEnergy)
 
void SetCoord (Float_t x1, Float_t y1, Float_t x2, Float_t y2)
 
void SetType (Char_t type)
 
Int_t CountNeighbors (const std::list< ecalCell * > &lst) const
 

Private Member Functions

 ClassDef (ecalModule, 1)
 

Private Attributes

Float_t fDx
 
Float_t fDy
 
std::vector< ecalCell * > fCells
 

Detailed Description

ecalModule.h

Author
Mikhail Prokudin

ECAL module structure, consisting of cells Useless if we have only modules with one lightisolated cell

Definition at line 16 of file ecalModule.h.

Constructor & Destructor Documentation

◆ ecalModule()

ecalModule::ecalModule ( char  type = 1,
Int_t  cellnumber = -1,
Float_t  x1 = 0,
Float_t  y1 = 0,
Float_t  x2 = 0,
Float_t  y2 = 0,
Int_t  mc = 0,
Float_t  energy = 0 
)

Definition at line 18 of file ecalModule.cxx.

19 : ecalCell(cellnumber, x1,y1,x2,y2, type, energy),
20 fDx(x2-x1),
21 fDy(y2-y1),
22 fCells()
23{
24 if (GetType()<1) return;
25
26 Int_t i;
27 Int_t j;
28 Int_t mt;
29 Int_t num;
30
31 mt=type;
32 fCells.resize(mt*mt,NULL);
33
34 num=cellnumber/100;
35
36 if (mc==0)
37 for(i=0;i<mt;i++)
38 for(j=0;j<mt;j++)
39 fCells[i*mt+j]=new ecalCell((num*10+i+1)*10+j+1, x1+j*fDx/mt, y1+i*fDy/mt, x1+(j+1)*fDx/mt, y1+(i+1)*fDy/mt, type);
40 else
41 for(i=0;i<mt;i++)
42 for(j=0;j<mt;j++)
43 fCells[i*mt+j]=new ecalCellMC((num*10+i+1)*10+j+1, x1+j*fDx/mt, y1+i*fDy/mt, x1+(j+1)*fDx/mt, y1+(i+1)*fDy/mt, type);
44}
Char_t GetType() const
Definition ecalCell.h:26
std::vector< ecalCell * > fCells
Definition ecalModule.h:52
Float_t fDx
Definition ecalModule.h:48
Float_t fDy
Definition ecalModule.h:50
int i
Definition ShipAna.py:86

Member Function Documentation

◆ AddEnergy()

void ecalModule::AddEnergy ( Float_t  x,
Float_t  y,
Float_t  energy 
)

Definition at line 64 of file ecalModule.cxx.

65{
66 ecalCell* tmp=FindCell(x,y);
67 if (!tmp) return;
68 tmp->AddEnergy(energy);
69 ecalCell::AddEnergy(energy);
70}
void AddEnergy(Float_t energy)
Definition ecalCell.h:68
ecalCell * FindCell(Float_t x, Float_t y) const

◆ At()

ecalCell * ecalModule::At ( Int_t  x,
Int_t  y 
) const
inline

Definition at line 26 of file ecalModule.h.

◆ ClassDef()

ecalModule::ClassDef ( ecalModule  ,
 
)
private

◆ FindCell()

ecalCell * ecalModule::FindCell ( Float_t  x,
Float_t  y 
) const

Definition at line 54 of file ecalModule.cxx.

55{
56 Int_t ix=(Int_t)TMath::Floor( (x-GetX1())/GetDX()*GetType() );
57 Int_t iy=(Int_t)TMath::Floor( (y-GetY1())/GetDY()*GetType() );
58 if (ix<0) ix=0; if (ix>=GetType()) ix=GetType()-1;
59 if (iy<0) iy=0; if (iy>=GetType()) iy=GetType()-1;
60 return At(ix,iy);
61}
Float_t GetY1() const
Definition ecalCell.h:32
Float_t GetX1() const
Definition ecalCell.h:31
Float_t GetDX() const
Definition ecalModule.h:37
ecalCell * At(Int_t x, Int_t y) const
Definition ecalModule.h:26
Float_t GetDY() const
Definition ecalModule.h:38

◆ GetCells()

std::vector< ecalCell * > ecalModule::GetCells ( ) const
inline

Definition at line 39 of file ecalModule.h.

39{return fCells;}

◆ GetCellsX()

list< ecalCell * > ecalModule::GetCellsX ( Float_t  x) const

Definition at line 73 of file ecalModule.cxx.

74{
75 list<ecalCell*> tmp;
76 vector<ecalCell*>::const_iterator p;
77
78 for(p=fCells.begin();p!=fCells.end();++p)
79 if (x>(*p)->GetX1()&&x<(*p)->GetX2()) tmp.push_back(*p);
80 return tmp;
81}

◆ GetCellsY()

list< ecalCell * > ecalModule::GetCellsY ( Float_t  y) const

Definition at line 84 of file ecalModule.cxx.

85{
86 list<ecalCell*> tmp;
87 vector<ecalCell*>::const_iterator p;
88
89 for(p=fCells.begin();p!=fCells.end();++p)
90 if (y>(*p)->GetY1()&&y<(*p)->GetY2()) tmp.push_back(*p);
91 return tmp;
92}

◆ GetDX()

Float_t ecalModule::GetDX ( ) const
inline

Definition at line 37 of file ecalModule.h.

37{return fDx;}

◆ GetDY()

Float_t ecalModule::GetDY ( ) const
inline

Definition at line 38 of file ecalModule.h.

38{return fDy;}

◆ GetEnergy()

Float_t ecalModule::GetEnergy ( Float_t  x,
Float_t  y 
) const
inline

Definition at line 29 of file ecalModule.h.

30 {
31 ecalCell* tmp=FindCell(x,y);
32 if (tmp) return tmp->GetEnergy();
33 return -1;
34 }

◆ Locate()

ecalCell * ecalModule::Locate ( Int_t  x,
Int_t  y 
) const

Definition at line 47 of file ecalModule.cxx.

48{
49 if (x<0||y<0||x>=GetType()||y>=GetType()) return NULL;
50 return fCells[y*GetType()+x];
51}

◆ ResetModule()

void ecalModule::ResetModule ( )
inline

Definition at line 57 of file ecalModule.h.

58{
60 for(UInt_t i=0;i<fCells.size();i++) fCells[i]->ResetEnergyFast();
61}
void ResetEnergyFast()
Definition ecalCell.h:109

Member Data Documentation

◆ fCells

std::vector<ecalCell*> ecalModule::fCells
private

list of cells contained in a module

Definition at line 52 of file ecalModule.h.

◆ fDx

Float_t ecalModule::fDx
private

module x-size

Definition at line 48 of file ecalModule.h.

◆ fDy

Float_t ecalModule::fDy
private

module y-size

Definition at line 50 of file ecalModule.h.


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