SND@LHC Software
Loading...
Searching...
No Matches
ecalAnalysisSimple.cxx
Go to the documentation of this file.
2
3#include "FairRootManager.h"
4
5#include "ecalStructure.h"
6#include "ecalCell.h"
7#include "ecalPoint.h"
8
9#include "TTree.h"
10#include "TClonesArray.h"
11
12#include <iostream>
13#include <fstream>
14#include <list>
15
16using namespace std;
17
19void ecalAnalysisSimple::Exec(Option_t* option)
20{
21 Int_t n=fTracks->GetEntries();
22 Int_t i;
23 ecalPoint* t;
24 ecalCell* cell;
25 ecalCell* mcell;
26 TVector3 m;
27 list<ecalCell*> cells;
28 list<ecalCell*>::const_iterator p;
29// TVector3 m1;
30
31 fEv++;
32 InitTree();
33
34 for(i=0;i<n;i++)
35 {
36 t=(ecalPoint*)fTracks->At(i);
37 fX=t->GetX();
38 fY=t->GetY();
39 t->Momentum(m);
40 fP=m.Mag();
41 fPX=m.Px();
42 fPY=m.Py();
43 fPZ=m.Pz();
44
45// m1=m.Unit();
46 cell=fStr->GetCell(fX, fY);
47 if (!cell) continue;
48 mcell=cell;
49 cell->GetNeighborsList(cells);
50 for(p=cells.begin();p!=cells.end();++p)
51 if ((*p)->GetEnergy()>mcell->GetEnergy())
52 mcell=(*p);
53
54 mcell->GetNeighborsList(cells);
55 for(p=cells.begin();p!=cells.end();++p)
56 if ((*p)->GetEnergy()>mcell->GetEnergy())
57 break;
58
59 if (p!=cells.end()) continue;
60
61 fCX=mcell->GetCenterX();
62 fCY=mcell->GetCenterY();
63 fCE=mcell->GetEnergy();
64 fCellNum=mcell->GetCellNumber();
65 fADC=mcell->ADC();
66 fOE=fCE;
67 for(p=cells.begin();p!=cells.end();++p)
68 fOE+=(*p)->GetEnergy();
69
70 fTree->Fill();
71 }
72
73}
74
76{
77 if (fTree) return;
78 fTree=new TTree("calib", "calib");
79 fTree->Branch("px", &fPX, "px/D");
80 fTree->Branch("py", &fPY, "py/D");
81 fTree->Branch("pz", &fPZ, "pz/D");
82 fTree->Branch("p" , &fP , "p/D");
83 fTree->Branch("x" , &fX , "x/D");
84 fTree->Branch("y" , &fY , "y/D");
85 fTree->Branch("cx", &fCX, "cx/D");
86 fTree->Branch("cy", &fCY, "cy/D");
87 fTree->Branch("ce", &fCE, "ce/D");
88 fTree->Branch("oe", &fOE, "oe/D");
89 fTree->Branch("ev", &fEv, "ev/I");
90 fTree->Branch("cn", &fCellNum, "cn/I");
91 fTree->Branch("adc", &fADC, "adc/I");
92}
93
94ecalAnalysisSimple::ecalAnalysisSimple(const char* name, const Int_t iVerbose)
95 : FairTask(name, iVerbose),
96 fTree(NULL),
97 fX(0.),
98 fY(0.),
99 fCX(0.),
100 fCY(0.),
101 fP(0.),
102 fCE(0.),
103 fOE(0.),
104 fPX(0.),
105 fPY(0.),
106 fPZ(0.),
107 fEv(0),
108 fCellNum(0),
109 fADC(0),
110 fStr(NULL),
111 fTracks(NULL)
112{
113}
114
116 : FairTask(),
117 fTree(NULL),
118 fX(0.),
119 fY(0.),
120 fCX(0.),
121 fCY(0.),
122 fP(0.),
123 fCE(0.),
124 fOE(0.),
125 fPX(0.),
126 fPY(0.),
127 fPZ(0.),
128 fEv(0),
129 fCellNum(0),
130 fADC(0),
131 fStr(NULL),
132 fTracks(NULL)
133{
134}
135
138{
139 FairRootManager* fManager=FairRootManager::Instance();
140 fStr=(ecalStructure*)fManager->GetObject("EcalStructure");
141 if (!fStr)
142 {
143 Fatal("Init()", "Can't find calorimeter structure. ");
144 return kFATAL;
145 }
146 fTracks=(TClonesArray*)fManager->GetObject("EcalPoint");
147 if (!fTracks)
148 {
149 Fatal("Init()", "Can't find array of reconstructed tracks. ");
150 return kFATAL;
151 }
152 return kSUCCESS;
153}
154
157{
158 if (fTree)
159 fTree->Write();
160}
161
Double_t m
virtual InitStatus Init()
virtual void Exec(Option_t *option)
Int_t GetCellNumber() const
Definition ecalCell.h:40
void GetNeighborsList(std::list< ecalCell * > &neib) const
Definition ecalCell.h:48
Float_t GetCenterX() const
Definition ecalCell.h:35
Short_t ADC() const
Definition ecalCell.h:37
Float_t GetCenterY() const
Definition ecalCell.h:36
Float_t GetEnergy() const
Definition ecalCell.h:42
ecalCell * GetCell(Float_t x, Float_t y) const
ClassImp(ecalContFact) ecalContFact