StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StjTrgSoftGetAdcEt.cxx
1 // $Id: StjTrgSoftGetAdcEt.cxx,v 1.5 2008/10/16 20:39:22 tai Exp $
2 // Copyright (C) 2008 Tai Sakuma <sakuma@bnl.gov>
3 #include "StjTrgSoftGetAdcEt.h"
4 
5 #include "StjTrgBEMCJetPatchTowerIdMap.h"
6 
7 #include "StjTowerEnergyCutBemcStatus.h"
8 #include "StjTowerEnergyCutEnergy.h"
9 #include "StjTowerEnergyCutEt.h"
10 #include "StjTowerEnergyCutAdc.h"
11 #include "StjTowerEnergyCutTowerId.h"
12 
13 #include "StjTrg.h"
14 
15 #include <TVector3.h>
16 #include <TMath.h>
17 
18 #include <iostream>
19 
20 ClassImp(StjTrgSoftGetAdcEt)
21 
22 using namespace std;
23 
25  : _bemc(bemc), _bemcJpTowerMap(bemcJpTowerMap), _trg(0)
26  , _runNumber(-1), _eventId(-1)
27 {
28  _cut.addCut( new StjTowerEnergyCutEnergy(0.0) );
29  _cut.addCut( new StjTowerEnergyCutBemcStatus(1) );
30  _cut.addCut( new StjTowerEnergyCutAdc(0, 2.0) );
31  _cut.addCut( new StjTowerEnergyCutEt(0.2) );
32 }
33 
34 bool StjTrgSoftGetAdcEt::isNewEvent()
35 {
36  if(_runNumber != _trg->runNumber()) return true;
37  if(_eventId != _trg->eventId()) return true;
38  return false;
39 }
40 
41 void StjTrgSoftGetAdcEt::read()
42 {
43  _towerAdc.clear();
44  _towerEnergy.clear();
45  _towerEt.clear();
46 
47  _jetPatchAdc.clear();
48  _jetPatchEnergy.clear();
49  _jetPatchEt.clear();
50 
51  vector<int> towers = _trg->towers();
52  vector<int> jetPatches = _trg->jetPatches();
53 
54  StjTowerEnergyList energyList = _bemc->getEnergyList();
55  energyList = _cut(energyList);
56 
57  for(vector<int>::const_iterator tower = towers.begin(); tower != towers.end(); ++tower) {
58  bool found = false;
59  for(StjTowerEnergyList::const_iterator it = energyList.begin(); it != energyList.end(); ++it) {
60  if((*it).towerId == (*tower)) {
61  _towerAdc.push_back((*it).adc);
62  _towerEnergy.push_back((*it).energy);
63  TVector3 vec3;
64  vec3.SetPtEtaPhi((*it).towerR, (*it).towerEta, (*it).towerPhi);
65  _towerEt.push_back(((*it).energy)*TMath::Sin(vec3.Theta()));
66  found = true;
67  break;
68  }
69  }
70  if(!found) {
71  _towerAdc.push_back(0);
72  _towerEnergy.push_back(0);
73  _towerEt.push_back(0);
74  }
75  }
76 
77  for(vector<int>::const_iterator jp = jetPatches.begin(); jp != jetPatches.end(); ++jp) {
78  unsigned int adc = 0;
79  double energy = 0;
80  double et = 0;
81  for(StjTowerEnergyList::const_iterator it = energyList.begin(); it != energyList.end(); ++it) {
82  if(*jp == _bemcJpTowerMap->getJetPatchIdForTower((*it).towerId)) {
83  adc += (*it).adc;
84  energy += (*it).energy;
85  TVector3 vec3;
86  vec3.SetPtEtaPhi((*it).towerR, (*it).towerEta, (*it).towerPhi);
87  et += ((*it).energy)*TMath::Sin(vec3.Theta());
88  }
89  }
90  _jetPatchAdc.push_back(adc);
91  _jetPatchEnergy.push_back(energy);
92  _jetPatchEt.push_back(et);
93  }
94 
95 
96 }
97 
98 vector<unsigned int> StjTrgSoftGetAdcEt::towerAdc()
99 {
100  if(isNewEvent()) read();
101  return _towerAdc;
102 }
103 
104 vector<double> StjTrgSoftGetAdcEt::towerEnergy()
105 {
106  if(isNewEvent()) read();
107  return _towerEnergy;
108 }
109 
110 vector<double> StjTrgSoftGetAdcEt::towerEt()
111 {
112  if(isNewEvent()) read();
113  return _towerEt;
114 }
115 
116 vector<unsigned int> StjTrgSoftGetAdcEt::jetPatchAdc()
117 {
118  if(isNewEvent()) read();
119  return _jetPatchAdc;
120 }
121 
122 vector<double> StjTrgSoftGetAdcEt::jetPatchEnergy()
123 {
124  if(isNewEvent()) read();
125  return _jetPatchEnergy;
126 }
127 
128 vector<double> StjTrgSoftGetAdcEt::jetPatchEt()
129 {
130  if(isNewEvent()) read();
131  return _jetPatchEt;
132 }