StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StEEmcPointMaker.cxx
1 #include "StEEmcPointMaker.h"
2 #include "StEEmcA2EMaker.h"
3 #include "StEEmcClusterMaker.h"
4 
5 #include "StEEmcUtil/EEmcGeom/EEmcGeomDefs.h"
6 #include "StEEmcUtil/EEmcGeom/EEmcGeomSimple.h"
7 #include "StEEmcUtil/StEEmcSmd/EEmcSmdGeom.h"
8 
9 #include <iostream>
10 #include <algorithm>
11 #include <map>
12 
13 #include "eeTowerFunction.h"
14 
15 /* StEvent stuff */
16 #include "StEvent/StEvent.h"
17 #include "StEvent/StEmcCollection.h"
18 #include "StEvent/StEmcPoint.h"
19 
20 /* Root's linear algebra package */
21 #include "TMatrixF.h"
22 
23 //#define DEBUG_buildPoints
24 
25 ClassImp(StEEmcPointMaker);
26 
27 // ----------------------------------------------------------------------------
28 StEEmcPointMaker::StEEmcPointMaker(const Char_t *name):StMaker(name)
29 {
30  std::cout << "StEEmcPointMaker("<<name<<")" << std::endl;
31 
35  mEEtow=new EEmcGeomSimple();
36  mEEsmd=EEmcSmdGeom::instance();
37  mEEmap=EEmcSmdMap::instance();
38 
39  mTowerThreshold=0.;
40  mFillStEvent=false;
41  mEnergyMode=1;
42  mLimit=10;
43 
44 }
45 
46 // ----------------------------------------------------------------------------
48 {
51  assert(mEEanalysis);
52  assert(mEEclusters);
53  return StMaker::Init();
54 }
55 
56 // ----------------------------------------------------------------------------
58 {
59 
66 
67  // Loop over all 12 EEMC sectors
68  for ( Int_t sector=0; sector<12; sector++ )
69  {
70 
72  StEEmcSmdClusterVec_t uclusters=mEEclusters->smdclusters(sector,0);
73  StEEmcSmdClusterVec_t vclusters=mEEclusters->smdclusters(sector,1);
74 
76  std::sort( uclusters.begin(), uclusters.end(), inner );
77  std::sort( vclusters.begin(), vclusters.end(), inner );
78 
79  findPoints( sector, uclusters, vclusters, mPoints );
80 
81  }
82 
83 
90  for ( UInt_t i=0;i<mPoints.size();i++ )
91  {
92  mPoints[i].energy( mPoints[i].energy()/0.007/2. );
93  }
94 
95 
96 
97 
98  StEEmcPointVec_t orgpoints=mPoints;
99 
101  // shareEnergy(); // <<<<<<< leads to negative point energies... rethink
102 
103  if ( mEnergyMode == 1 )
105  else
106  shareEnergySmd();
107 
108 
110  countRelatives();
111 
112 
113  if ( mFillStEvent )
114  {
115  fillStEvent();
116  verifyStEvent();
117  }
118 
119 
120  return kStOK;
121 }
122 
123 // ----------------------------------------------------------------------------
124 
125 StEEmcPointVec_t StEEmcPointMaker::buildSmdPoints( Int_t sector,
126  StEEmcSmdClusterVec_t &u,
127  StEEmcSmdClusterVec_t &v )
128 {
129 
130  StEEmcPointVec_t points;
131 
132  for ( UInt_t i=0; i<u.size(); i++ )
133  {
134 
135  StEEmcSmdCluster uc=u[i];
136  Float_t xu=uc.mean();
137 
138  for ( UInt_t j=0;j<v.size(); j++ )
139  {
140 
141 
142  StEEmcSmdCluster vc=v[j];
143  Float_t xv=vc.mean();
144 
146  TVector3 direct = mEEsmd->getIntersection( sector, xu, xv );
147 
148  Int_t sec,sub,eta;
149  if ( !mEEtow->getTower(direct,sec,sub,eta) )
150  {
151  continue;
152  }
153  else
154  {
156  }
157 
160  if ( sector != sec )
161  continue;
162 
163 
164 
168  Bool_t good = false;
169  if ( mEEanalysis->tower(sec,sub,eta).energy() > 0. )
170  good=true;
171  else if ( mEEanalysis->tower(sec,sub,eta).fail() )
172  {
173  for ( Int_t layer=1;layer<=3;layer++ )
174  {
175  if ( mEEanalysis->tower(sec,sub,eta,layer).energy() > 0. )
176  good=true;
177  }
178  }
179 
180  if ( good ) {
181 
182  StEEmcPoint p;
183  p.cluster( uc, 0 );
184  p.cluster( vc, 1 );
185  p.energy( uc.energy() + vc.energy() );
186  p.tower( mEEanalysis->tower(sec,sub,eta) );
187  TVector3 position=mEEsmd->getIntersection(sector,uc.mean(),vc.mean());
188  p.position(position);
189  points.push_back(p);
190 
191  }
192 
193  }
194 
195  }
196 
197  return points;
198 }
199 
200 
201 
202 
203 
204 
205 
206 
207 
208 
209 
210 
211 
212 
213 
214 
215 
216 // ----------------------------------------------------------------------------
217 Bool_t StEEmcPointMaker::findPoints( Int_t sector,
218  StEEmcSmdClusterVec_t uclusters,
219  StEEmcSmdClusterVec_t vclusters,
220  StEEmcPointVec_t &points )
221 {
222 
223 
225  StEEmcPointVec_t mypoints;
226 
228  std::sort( uclusters.begin(), uclusters.end(), inner );
229  std::sort( vclusters.begin(), vclusters.end(), inner );
230 
232  StEEmcPointVec_t smdpoints = buildSmdPoints( sector, uclusters, vclusters );
233 
234 
236  if ( smdpoints.size() < 1 ) return false;
237 
239  std::sort( smdpoints.begin(), smdpoints.end(), chiSquare );
240 
241 
245  std::map< Int_t, std::vector<Int_t> > u2p, v2p;
246  for ( UInt_t i=0; i<smdpoints.size(); i++ )
247  {
248  u2p[ smdpoints[i].cluster(0).key() ].push_back( i );
249  v2p[ smdpoints[i].cluster(1).key() ].push_back( i );
250  }
251 
257  StEEmcSmdClusterVec_t::iterator uiter=uclusters.begin();
258  StEEmcSmdClusterVec_t::iterator viter=vclusters.begin();
259 
260 
261 
262  // ----------<<<<<<<<< stage one >>>>>>>>>-------------
263 
264 
265 
266  while ( uiter<uclusters.end() || viter<vclusters.end() )
267  {
268 
270  StEEmcSmdCluster ucl;ucl.key(-1);
271  StEEmcSmdCluster vcl;vcl.key(-1);
272  if ( uiter<uclusters.end() ) ucl=(*uiter);
273  if ( viter<vclusters.end() ) vcl=(*viter);
274  Int_t iUV=-1;
275  if ( ucl.key()<0 )
276  iUV=1;
277  else if ( vcl.key()<0 )
278  iUV=0;
279  else if ( (*uiter).mean() < (*viter).mean() )
280  iUV=0;
281  else
282  iUV=1;
283 
285  StEEmcSmdCluster &cluster=(iUV==0)?ucl:vcl;
286 
288  std::vector<Int_t> matched=(iUV==0)?
289  u2p[cluster.key()]:
290  v2p[cluster.key()];
291 
292 
295  if ( matched.size()==0 || matched.size() >1 )
296  {
297  if ( iUV==0 ) uiter++;
298  if ( iUV==1 ) viter++;
299  continue;
300  }
301 
307 
310  StEEmcPoint p=smdpoints[matched.back()];
311 
313  mypoints.push_back(p);
314 
315  if ( iUV==0 ) uiter++;
316  if ( iUV==1 ) viter++;
317 
318  }
319 
320 
323  Float_t chisq=9.0E9;
324  Int_t imin=-1;
325  for ( UInt_t i=0; i<mypoints.size(); i++ )
326  {
327  Float_t eu=mypoints[i].cluster(0).energy();
328  Float_t ev=mypoints[i].cluster(1).energy();
329  Float_t x2=(eu-ev)*(eu-ev);
330  if ( x2 < chisq ) {
331  imin=(Int_t)i;
332  chisq=x2;
333  }
334  }
335 
336 
343  if ( imin >= 0 ) {
344 
345  StEEmcPoint p=mypoints[imin];
346  removeCluster( uclusters, p.cluster(0).key() );
347  removeCluster( vclusters, p.cluster(1).key() );
348  points.push_back(p);
349  findPoints(sector, uclusters, vclusters, points );
350  return true;
351 
352  }
353 
354 
355 
356 
357 
359 
360 
367  uiter=uclusters.begin();
368  viter=vclusters.begin();
369  while ( uiter!=uclusters.end() || viter!=vclusters.end() )
370  {
371 
373  StEEmcSmdCluster ucl;ucl.key(-1);
374  StEEmcSmdCluster vcl;vcl.key(-1);
375  if ( uiter!=uclusters.end() ) ucl=(*uiter);
376  if ( viter!=vclusters.end() ) vcl=(*viter);
377  Int_t iUV=-1;
378  if ( ucl.key()<0 )
379  iUV=1;
380  else if ( vcl.key()<0 )
381  iUV=0;
382  else if ( (*uiter).mean() < (*viter).mean() )
383  iUV=0;
384  else
385  iUV=1;
386 
388  StEEmcSmdCluster &cluster=(iUV==0)?ucl:vcl;
389 
391  std::vector<Int_t> matched=(iUV==0)?
392  u2p[cluster.key()]:
393  v2p[cluster.key()];
394 
397  if ( matched.size()==0 || matched.size()==1 )
398  {
399  if ( iUV==0 ) uiter++;
400  if ( iUV==1 ) viter++;
401  continue;
402  }
403 
406  StEEmcPoint p=smdpoints[matched.front()];
407 
409  mypoints.push_back(p);
410 
411  if ( iUV==0 ) uiter++;
412  if ( iUV==1 ) viter++;
413 
414  }
415 
416 
419  chisq=9.0E9;
420  imin=-1;
421  for ( UInt_t i=0; i<mypoints.size(); i++ )
422  {
423  Float_t eu=mypoints[i].cluster(0).energy();
424  Float_t ev=mypoints[i].cluster(1).energy();
425  Float_t x2=(eu-ev)*(eu-ev);
426  if ( x2 < chisq ) {
427  imin=(Int_t)i;
428  chisq=x2;
429  }
430  }
431 
432 
433 
440  if ( imin >= 0 ) {
441 
442  StEEmcPoint p=mypoints[imin];
443  removeCluster( uclusters, p.cluster(0).key() );
444  removeCluster( vclusters, p.cluster(1).key() );
445  points.push_back(p);
446  findPoints(sector, uclusters, vclusters, points );
447  return true;
448 
449  }
450 
451 
452  return true;
453 
454 }
455 
456 
457 
458 // ----------------------------------------------------------------------------
459 void StEEmcPointMaker::Clear( Option_t *opts )
460 {
461 
462  mEseen=0.;
463  mPoints.clear();
464 
465 }
466 
467 // ----------------------------------------------------------------------------
468 void StEEmcPointMaker::removeCluster( StEEmcSmdClusterVec_t &clusters, Int_t k )
469 {
470  StEEmcSmdClusterVec_t::iterator iter=clusters.begin();
471  while ( iter != clusters.end() )
472  {
473  if ( (*iter).key() == k ) {
474  clusters.erase(iter);
475  return;
476  }
477  iter++;
478  }
479 }
480 
481 // ----------------------------------------------------------------------------
483 {
484 
490 
491 
492 
493  Int_t nrow=(Int_t)mPoints.size();
494  Int_t ncol=(Int_t)mPoints.size();
495  if ( nrow < 1 ) return;
496 
497  std::vector<Float_t> Ef(nrow,0.);
498 
499  TMatrixF fractions(nrow,ncol);
500 
502  for ( Int_t k=0; k<mEEanalysis->numberOfHitTowers(0); k++ )
503  {
504 
506  StEEmcTower tower=mEEanalysis->hittower(k,0);
507 
509  for ( UInt_t i=0; i< mPoints.size(); i++ )
510  {
511 
514  Float_t fi=fracp2t( mPoints[i], tower );
515  if ( fi<=0. ) continue;
516 
518  Ef[i] += fi * tower.energy();
519 
520  for ( UInt_t j=0; j<mPoints.size(); j++ )
521  {
522 
525  Float_t fj=fracp2t( mPoints[j], tower );
526  if (fi*fj<=0.) continue;
527 
528  fractions[i][j] += fi*fj;
529 
530  }
531 
532  }
533 
534  }
535 
536  fractions.Print();
537 
539  Double_t det = 0.;
540  TMatrixF invert= fractions;
541  invert.Invert(&det);
542 
543  invert.Print();
544 
545  TMatrixF test= fractions * invert;
546 
547  test.Print();
548 
549 
551 
552  std::vector<Float_t> epoints(nrow,0.);
553 
554  for ( Int_t i=0; i<nrow; i++ )
555  {
556  for ( Int_t j=0; j<ncol; j++ )
557  {
558  epoints[i] += invert[i][j] * Ef[j];
559  }
560 
561 
562  }
563 
564 
565 
566 
567 
568 
569 }
570 
571 
572 // ----------------------------------------------------------------------------
574 {
575 
577 
578 
581  if ( !t.isNeighbor( p.tower(0) ) ) return 0.;
582  //if ( !(t.index()==p.tower(0).index()) ) return 0.;
583 
584 
586  Float_t xeta=(Float_t)t.etabin();
587  Float_t xphi=(Float_t)t.phibin();
588  Double_t X[]={xphi,xeta};
589 
591  Float_t xeta0=(Float_t)p.tower(0).etabin();
592  Float_t xphi0=(Float_t)p.tower(0).phibin();
593 
596  Int_t sec,sub,eta;
597  Float_t dphi,deta;
598  if ( !mEEtow->getTower(p.position(), sec,sub,eta, dphi,deta ) ) return 0.;
599  dphi/=2.0;
600  deta/=2.0;
601 
603  Float_t xetap=xeta0+deta;
604  Float_t xphip=xphi0+dphi;
605  Double_t P[]={xphip,xetap,1.0};
606 
607  return eeTowerFunction( X, P );
608 
609 }
610 
611 
612 // ----------------------------------------------------------------------------
614 {
615 
616  Int_t limit=mLimit; // algo quickly converges on energy
617  Int_t count=0;
618 
619  while ( count++ < limit )
620  {
621 
624  Float_t sumw[720];
625  for (Int_t i=0;i<720;i++) sumw[i]=0.;
626 
629  for ( UInt_t i=0;i<mPoints.size();i++ )
630  {
631 
633  StEEmcTower tower=point.tower(0);
634 
635  sumw[ tower.index() ] += point.energy() * fracp2t( point, tower );
636 
638  for ( Int_t i=0; i<tower.numberOfNeighbors(); i++ )
639  {
640  StEEmcTower neighbor=tower.neighbor(i);
641  sumw[ neighbor.index() ] += point.energy() * fracp2t( point, neighbor );
642  }
643 
644  }
645 
646 
650 
651  for ( UInt_t i=0;i<mPoints.size();i++ )
652  {
653 
654  StEEmcPoint &point=mPoints[i]; // note the reference
655  StEEmcTower tower=point.tower(0);
656  Float_t epoint=0.;
657 
658  Float_t frac=0.;
659  if ( !tower.fail() && !tower.stat() && sumw[tower.index()]>0. )
660  frac = point.energy() * fracp2t(point,tower) / sumw[ tower.index() ];
661 
662  epoint += tower.energy() * frac;
663 
666  for ( Int_t i=0; i<tower.numberOfNeighbors(); i++ )
667  {
668  StEEmcTower neighbor=tower.neighbor(i);
669  if ( neighbor.stat() || neighbor.fail() || sumw[neighbor.index()]<=0. ) continue;
670  frac = point.energy() * fracp2t(point,neighbor) / sumw[ neighbor.index() ];
671  epoint += frac * neighbor.energy();
672  }
673 
674 
676  point.energy( epoint );
677 
678  }
679 
680  }
681 
682 
683 
684 
685  std::vector<Bool_t> seen(720,false);
686  for ( UInt_t i=0;i<mPoints.size();i++ )
687  {
688 
689  StEEmcTower tow=mPoints[i].tower(0);
690  if ( !seen[ tow.index() ] ) mEseen+=tow.energy();
691  seen[ tow.index() ] = true;
692 
693  for ( Int_t j=0;j<tow.numberOfNeighbors();j++ )
694  {
695  StEEmcTower nei=tow.neighbor(j);
696  if ( !seen[ nei.index() ] ) mEseen += nei.energy();
697  seen[ nei.index() ] = true;
698  }
699 
700  }
701 
702 
703 
704 
705 
706 
707 }
708 
709 
710 
711 
712 
713 // ----------------------------------------------------------------------------
715 {
716 
717 
718  StEvent *stevent=(StEvent*)GetInputDS("StEvent");
719  if ( !stevent ) {
720  Warning("fillStEvent","called, but no StEvent to be found");
721  return;
722  }
723 
724 
726  for ( UInt_t i=0; i<mPoints.size(); i++ )
727  {
728 
729  StEmcPoint *point=mPoints[i].stemc();
730  stevent->emcCollection()->addEndcapPoint( point );
731 
732  mEtoEE[ point ] = mPoints[i];
733 
734 
735  }
736 
737 }
738 
739 
740 
741 // ----------------------------------------------------------------------------
743 {
744 
745  Float_t emc_sum_points = 0.;
746  Float_t eemc_sum_points = 0.;
747 
748  StEvent *stevent=(StEvent*)GetInputDS("StEvent");
749  if ( !stevent ) {
750  Warning("verifyStEvent","called, but no StEvent to be found");
751  return;
752  }
753 
754  StSPtrVecEmcPoint& emcpts = stevent->emcCollection()->endcapPoints();
755  for ( UInt_t i=0;i<emcpts.size();i++ )
756  {
757 
758  StEmcPoint *p=emcpts[i];
759  assert(p);
760  emc_sum_points += p->energy();
761 
762  }
763 
764  for ( UInt_t i=0;i<mPoints.size();i++ )
765  {
766 
767  StEEmcPoint p=mPoints[i];
768  eemc_sum_points += p.energy();
769 
770  }
771 
772  std::cout << "StEEmcPointMaker point checksum: ";
773  if ( emc_sum_points == eemc_sum_points )
774  {
775  std::cout << "passed";
776  }
777  else
778  std::cout << "FAILED";
779  std::cout << std::endl;
780 
781 
782 }
783 
784 
785 // ----------------------------------------------------------------------------
787 {
788 
790  Int_t npoints[720];
791  for ( Int_t i=0;i<720;i++ ) npoints[i]=0;
792 
793  for ( UInt_t i=0;i<mPoints.size();i++ )
794  npoints[ mPoints[i].tower(0).index() ]++;
795 
797  for ( UInt_t i=0;i<mPoints.size();i++ )
798  {
799 
800  StEEmcTower tower=mPoints[i].tower(0);
801  Int_t nn=tower.numberOfNeighbors();
802 
803  Int_t nrel=npoints[ tower.index() ] - 1; // don't count self
804  assert(nrel>=0); // pbck
805 
806  for ( Int_t j=0;j<nn;j++ )
807  {
808  StEEmcTower t2=tower.neighbor(j);
809  nrel+=npoints[ t2.index() ];
810  }
811 
812  mPoints[i].numberOfRelatives(nrel);
813 
814  }
815 
816 
817 }
818 
819 
820 
821 
822 // ----------------------------------------------------------------------------
824 {
825 
828  Float_t sumw[720];for ( Int_t i=0;i<720;i++ )sumw[i]=0.;
829 
830  for ( UInt_t ipoint=0;ipoint<mPoints.size();ipoint++ )
831  {
832 
833 
834  StEEmcPoint point=mPoints[ipoint];
835  StEEmcTower tower=point.tower(0);
836  sumw[tower.index()]+=point.energy();
837 
838  for ( Int_t itow=0;itow<tower.numberOfNeighbors();itow++ )
839  {
840  StEEmcTower t=tower.neighbor(itow);
841  Int_t index=t.index();
842  sumw[index]+=point.energy();
843  }
844 
845  }
846 
849 
850  for ( UInt_t ipoint=0;ipoint<mPoints.size();ipoint++ )
851  {
852 
853  StEEmcPoint point=mPoints[ipoint]; // note reference
854  StEEmcTower tower = point.tower(0);
855  Float_t epoint = 0.;
856  Int_t index = tower.index();
857  epoint += tower.energy() * point.energy() / sumw[index];
858 
859  for ( Int_t itow=0;itow<tower.numberOfNeighbors();itow++ )
860  {
861  StEEmcTower t=tower.neighbor(itow);
862  index = t.index();
863  epoint += t.energy() * point.energy() / sumw[index];
864  }
865 
866  // point.energy( epoint );
867  mPoints[ipoint].energy(epoint);
868 
869  }
870 
871 
872 }
873 
874 
StEEmcA2EMaker * mEEanalysis
ADC2E.
void cluster(const StEEmcSmdCluster &c, Int_t plane)
Add an smd cluster to this point.
Definition: StEEmcPoint.h:40
void shareEnergySmd()
Divide energy of eemc towers between identified smd points in proportion to the smd energy...
Int_t numberOfNeighbors() const
get the number of neighboring towers
Definition: StEEmcTower.h:54
EEmc ADC –&gt; energy maker.
Base class for representing EEMC points.
Definition: StEEmcPoint.h:24
void stat(unsigned s)
Set a status bit for this element.
Definition: StEEmcElement.h:23
Bool_t isNeighbor(const StEEmcTower &t) const
Definition: StEEmcTower.cxx:98
Float_t mEseen
Energy seen by the algorithm.
StEEmcTower & hittower(Int_t hit, Int_t layer)
void energy(Float_t e, Int_t layer=0)
Set the energy of this point.
Definition: StEEmcPoint.h:34
StEEmcPointVec_t mPoints
All fully reconstructed points.
Bool_t mFillStEvent
Option to fill StEvent.
void shareEnergySimple()
Divide energy of eemc towers between identified smd points (doesn&#39;t work as well as smd algo) ...
void neighbor(StEEmcTower *n)
add a tower to list of neighbors
Definition: StEEmcTower.h:52
void verifyStEvent()
Checks that StEvent is properly saved.
StEEmcPoint point(Int_t ipoint)
Return a specified point.
EEmcSmdMap * mEEmap
Tower to smd map.
StEEmcPointVec_t points()
Return vector of all points found in endcap.
Int_t Init()
Initialize.
void shareEnergy()
Divide energy of eemc towers between identified smd points using fit (doesn&#39;t work) ...
Class for building points from smd clusters.
Int_t numberOfHitTowers(Int_t layer) const
Int_t etabin() const
Returns the etabin of this tower, pre- or postshower element.
Definition: StEEmcTower.h:45
void fail(unsigned f)
Set a fail bit for this element.
Definition: StEEmcElement.h:25
StEEmcPointMaker(const Char_t *name)
std::map< StEmcPoint *, StEEmcPoint > mEtoEE
Map connecting StEEmcPoint to StEmcPoint.
void index(Int_t i)
Definition: StEEmcTower.cxx:76
Base class for representing tower, preshower and postshower elements.
Definition: StEEmcTower.h:11
StEEmcPointVec_t buildSmdPoints(Int_t sector, StEEmcSmdClusterVec_t &u, StEEmcSmdClusterVec_t &v)
build smd points and associations between smd points and clusters
Bool_t findPoints(Int_t sector, StEEmcSmdClusterVec_t u, StEEmcSmdClusterVec_t v, StEEmcPointVec_t &p)
find points in the endcap
void countRelatives()
Determine the number of points which share tower energy with another point.
Int_t key()
Return a unique key assigned by the cluster maker.
EEMC simple geometry.
Definition: Stypes.h:40
A base class for representing clusters of EEMC smd strips.
void removeCluster(StEEmcSmdClusterVec_t &clusters, Int_t key)
Remove a cluster from the list of clusters.
void tower(const StEEmcTower &t, Float_t w=1.)
Add a tower with specified weight to the point.
Definition: StEEmcPoint.h:38
Int_t phibin() const
Returns the phibin of this tower.
Definition: StEEmcTower.h:47
TVector3 getIntersection(Int_t iSec, Int_t iUStrip, Int_t iVStrip, const TVector3 &vertex) const
void position(const TVector3 &p)
Set the position of this point at the SMD plane.
Definition: StEEmcPoint.h:32
EEmcGeomSimple * mEEtow
Tower geometry.
StEEmcTower & tower(Int_t index, Int_t layer=0)
bool getTower(const TVector3 &r, int &sec, int &sub, int &etabin, Float_t &dphi, Float_t &deta) const
void Clear(Option_t *opts="")
Clear old points.
EEmcSmdGeom * mEEsmd
Smd geometry.
Int_t mLimit
How many iterations for the tower energy sharing mode.
void energy(Float_t e)
Set the energy (adc-ped+0.5)/gain for this element.
Definition: StEEmcElement.h:21
Float_t fracp2t(StEEmcPoint &p, StEEmcTower &t)
void fillStEvent()
Fills the StEmcPoint collection.
StEEmcClusterMaker * mEEclusters
Clusters.
Int_t Make()
Build points for this event.
Int_t mEnergyMode
Option for dividing energy.
StEEmcSmdClusterVec_t & smdclusters(Int_t sec, Int_t plane)
Return a vector of smd clusters.
A cluster maker for the EEMC.