StFms  0.0.0
FMS software in the STAR framework
StFmsEventClusterer.h
Go to the documentation of this file.
1 // $Id$
2 //
3 // $Log$
13 #ifndef STROOT_STFMSPOINTMAKER_STFMSEVENTCLUSTERER_H_
14 #define STROOT_STFMSPOINTMAKER_STFMSEVENTCLUSTERER_H_
15 
16 #include <vector>
17 
18 #include "TObject.h"
19 
21 
22 namespace FMSCluster { // $NMSPC
23 class StFmsClusterFitter;
24 class StFmsFittedPhoton;
25 class StFmsGeometry;
26 class StFmsEventClusterer: public TObject {
27  public:
36  StFmsEventClusterer(const StFmsGeometry* geometry, Int_t detectorId);
44  int detector() const { return mDetectorId; }
52  Bool_t cluster(std::vector<FMSCluster::StFmsTower>* towers);
53 #ifndef __CINT__ // Hide ClusterList from CINT as it uses C++11
54 
57  const ClusterList& clusters() const { return mClusters; }
58 #endif // __CINT__
59 
60  private:
61 #ifndef __CINT__ // Hide ClusterList from CINT as it uses C++11
62 
63  typedef ClusterList::iterator ClusterIter;
65  typedef ClusterList::const_iterator ClusterConstIter;
66 #endif // __CINT__
67 
85  Int_t fitEvent();
91  Double_t photonEnergyInCluster(Double_t towerWidth,
93  const StFmsFittedPhoton* photon) const;
99  Double_t photonEnergyInTower(Double_t towerWidth, const StFmsTower* tower,
100  const StFmsFittedPhoton* photon) const;
107 #ifndef __CINT__ // Hide Cluster(Const)Iter from CINT as it uses C++11
108  /*
109  Perform a global fit of all photons in an event.
110 
111  Update the (x, y) positions end energies of the photons in each cluster based
112  on a global fit including all photons.
113  Only makes sense when there is more than one photon in the event.
114  Arguments:
115  - nPh, number of photons in the event
116  - nCl, number of clusters containing those photons
117  - first, iterator to the first cluster
118 
119  Returns the &chi;<sup>2</sup> of the fit.
120  */
121  Float_t globalFit(const Int_t, const Int_t, ClusterIter first);
122  /*
123  Special 2-photon fit for a single cluster.
124 
125  Cluster moments must have been calculated first
126 
127  Returns the &chi;<sup>2</sup> of the fit.
128  */
130  /*
131  Run tests on the lower-energy photon in a 2-photon cluster.
132 
133  Return true if the photon passes tests, in which case it is a real photon.
134  Return false if it fails, in which case it is a bogus photon due to some
135  problem in reconstruction - the cluster is actually a 1-photon cluster.
136 
137  Arguments:
138  - clusterIndex: index of cluster to test
139  - nRealClusters: total number of clusters in the event
140  */
143 #endif // __CINT__
146  Int_t mDetectorId;
147  std::vector<FMSCluster::StFmsTower>* mTowers;
149  std::vector<Float_t> mTowerWidthXY;
150  ClassDef(StFmsEventClusterer, 0)
151 };
152 } // namespace FMSCluster
153 #endif // STROOT_STFMSPOINTMAKER_STFMSEVENTCLUSTERER_H_
bool validate2ndPhoton(ClusterConstIter cluster) const
Double_t photonEnergyInCluster(Double_t towerWidth, const StFmsTowerCluster *cluster, const StFmsFittedPhoton *photon) const
Float_t fit2PhotonClust(ClusterIter cluster)
std::vector< Float_t > mTowerWidthXY
Geometry for this sub-detector (cm)
const StFmsGeometry * mGeometry
FMS geometry for current run.
Declaration of StFmsClusterFinder, an FMS tower clustering algorithm.
ClusterList mClusters
List of clusters in this sub-detector/event.
StFmsClusterFinder mClusterFinder
Cluster-finding routine.
Double_t photonEnergyInTower(Double_t towerWidth, const StFmsTower *tower, const StFmsFittedPhoton *photon) const
const ClusterList & clusters() const
Float_t globalFit(const Int_t, const Int_t, ClusterIter first)
std::list< std::unique_ptr< StFmsTowerCluster > > ClusterList
ClusterList::const_iterator ClusterConstIter
std::vector< FMSCluster::StFmsTower > * mTowers
Towers to cluster.
StFmsEventClusterer & operator=(const StFmsEventClusterer &)
StFmsEventClusterer(const StFmsGeometry *geometry, Int_t detectorId)
Float_t fitOnePhoton(StFmsTowerCluster *cluster)
Int_t mDetectorId
ID of this FMS sub-detector.
StFmsClusterFitter * mFitter
Routine for fitting photons to clusters.
Bool_t cluster(std::vector< FMSCluster::StFmsTower > *towers)