FEDRA emulsion software from the OPERA Collaboration
EdbShowerAlg Class Reference

Root Class Definition for my Objects. More...

#include <EdbShowerAlg.h>

Inheritance diagram for EdbShowerAlg:
Collaboration diagram for EdbShowerAlg:

Public Member Functions

void AddRecoShowerArray (TObjArray *RecoShowerArray)
 
EdbVertexCalcVertex (TObjArray *segments)
 
 ClassDef (EdbShowerAlg, 1)
 
void Convert_EdbPVRec_To_InBTArray ()
 
Double_t DeltaR_NoPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaR_WithoutPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaR_WithPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaTheta (EdbSegP *s1, EdbSegP *s2)
 
Double_t DeltaThetaComponentwise (EdbSegP *s1, EdbSegP *s2)
 
Double_t DeltaThetaSingleAngles (EdbSegP *s1, EdbSegP *s2)
 
 EdbShowerAlg ()
 
 EdbShowerAlg (TString AlgName, Int_t AlgValue)
 
virtual void Execute ()
 
virtual void Finalize ()
 
TString GetAlgName () const
 
Int_t GetAlgValue () const
 
Double_t GetMinimumDist (EdbSegP *seg1, EdbSegP *seg2)
 
TObjArray * GetRecoShowerArray () const
 
Int_t GetRecoShowerArrayN () const
 
Double_t GetSpatialDist (EdbSegP *s1, EdbSegP *s2)
 
virtual void Initialize ()
 
Bool_t IsInConeTube (EdbSegP *sTest, EdbSegP *sStart, Double_t CylinderRadius, Double_t ConeAngle)
 
void Print ()
 
void PrintAll ()
 
void PrintMore ()
 
void PrintParameters ()
 
void PrintParametersShort ()
 
void PrintRecoShowerArray ()
 
void SetActualAlgParameterset (Int_t ActualAlgParametersetNr)
 
void SetDebug ()
 
void SetEdbPVRec (EdbPVRec *Ali)
 
void SetEdbPVRecPIDNumbers (Int_t FirstPlate_eAliPID, Int_t LastPlate_eAliPID, Int_t MiddlePlate_eAliPID, Int_t NumberPlate_eAliPID)
 
void SetInBTArray (TObjArray *InBTArray)
 
void SetParameter (Int_t parNr, Float_t parvalue)
 
void SetParameters (Float_t *par)
 
void SetRecoShowerArray (TObjArray *RecoShowerArray)
 
void SetRecoShowerArrayN (Int_t RecoShowerArrayN)
 
void SetUseAliSub (Bool_t UseAliSub)
 
void Transform_eAli (EdbSegP *InitiatorBT, Float_t ExtractSize)
 
virtual ~EdbShowerAlg ()
 virtual constructor due to inherited class More...
 

Protected Member Functions

void Set0 ()
 

Protected Attributes

Int_t eActualAlgParametersetNr
 Used when more sets of same algorithm. More...
 
TString eAlgName
 
Int_t eAlgValue
 
EdbPVReceAli
 
EdbPVReceAli_Sub
 
Int_t eAli_SubNpat
 
Int_t eAliNpat
 
Bool_t eDebug
 
Int_t eFirstPlate_eAliPID
 
TObjArray * eInBTArray
 
Int_t eInBTArrayN
 
Int_t eLastPlate_eAliPID
 
Int_t eMiddlePlate_eAliPID
 
Int_t eNumberPlate_eAliPID
 
TString eParaString [10]
 
Float_t eParaValue [10]
 
EdbTrackPeRecoShower
 
TObjArray * eRecoShowerArray
 
Int_t eRecoShowerArrayN
 
Int_t eUseAliSub
 

Detailed Description

Root Class Definition for my Objects.

Constructor & Destructor Documentation

◆ EdbShowerAlg() [1/2]

EdbShowerAlg::EdbShowerAlg ( )

◆ EdbShowerAlg() [2/2]

EdbShowerAlg::EdbShowerAlg ( TString  AlgName,
Int_t  AlgValue 
)

Reset all:

29{
31 Set0();
32
33 eAlgName=AlgName;
34 eAlgValue=AlgValue;
35 for (int i=0; i<10; ++i) {
36 eParaValue[i]=-99999.0;
37 eParaString[i]="UNSPECIFIED";
38 }
39 eAli_Sub=0;
40}
TString eParaString[10]
Definition: EdbShowerAlg.h:42
Int_t eAlgValue
Definition: EdbShowerAlg.h:40
void Set0()
Definition: EdbShowerAlg.cxx:53
TString eAlgName
Definition: EdbShowerAlg.h:39
Float_t eParaValue[10]
Definition: EdbShowerAlg.h:41
EdbPVRec * eAli_Sub
Definition: EdbShowerAlg.h:63

◆ ~EdbShowerAlg()

EdbShowerAlg::~EdbShowerAlg ( )
virtual

virtual constructor due to inherited class

Default Destructor

45{
47 cout << "EdbShowerAlg::~EdbShowerAlg()"<<endl;
48}

Member Function Documentation

◆ AddRecoShowerArray()

void EdbShowerAlg::AddRecoShowerArray ( TObjArray *  RecoShowerArray)

Add an array to the existing RecoShower array.
For retrieving objects and handling them: be carefull that
you add objects of the same class.

109{
113
114 Log(3,"EdbShowerAlg::AddRecoShowerArray","AddRecoShowerArray...");
115
116 if (RecoShowerArray==NULL || RecoShowerArray->GetEntries()==0) return;
117
118 // TObject* obj = RecoShowerArray->At(0);
119 // if (obj->ClassName()!="EdbTrackP") { cout << "WARNING AddRecoShowerArray ClassNames dont match!"<< endl; cout << obj->ClassName() << endl; return; } // libShower
120 // if (obj->ClassName()!="EdbShowerP") { cout << "WARNING AddRecoShowerArray ClassNames dont match!"<< endl; return; } // libShowRec
121
122 for (int i=0; i<RecoShowerArray->GetEntries(); ++i) {
123 EdbTrackP* obj = (EdbTrackP*)RecoShowerArray->At(i); // libShower
124 //EdbShowerP* obj = (EdbShowerP*)RecoShowerArray->At(0); // libShowRec
125 eRecoShowerArray->Add(obj);
127 }
128
129 return;
130}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
TObjArray * RecoShowerArray
Definition: Shower_E_FromShowerRoot.C:12
TObjArray * eRecoShowerArray
Definition: EdbShowerAlg.h:59
Int_t eRecoShowerArrayN
Definition: EdbShowerAlg.h:60
Definition: EdbPattern.h:113
#define NULL
Definition: nidaqmx.h:84

◆ CalcVertex()

EdbVertex * EdbShowerAlg::CalcVertex ( TObjArray *  segments)

Exactly same implementation as in EdbEDAUtil.C
but I do not want to link in hardcoded all eda libraries.
calc vertex point with given segments. just topological calculation.
VTA is currently not set.
in case of Single-stop. make a vertex 650 micron upstream ob base.

353 {
354
360
361
362 double xx,yy,zz,Det,Ax,Ay,Az,Bx,By,Bz,Cx,Cy,Cz,Dx,Dy,Dz;
363 Ax=0.;
364 Ay=0.;
365 Az=0.;
366 Bx=0.;
367 By=0.;
368 Bz=0.;
369 Cx=0.;
370 Cy=0.;
371 Cz=0.;
372 Dx=0.;
373 Dy=0.;
374 Dz=0.;
375
376 if (segments->GetEntries()==1) {
377 // in case of Single-stop. make a vertex 650 micron upstream ob base.
378 EdbSegP *s = new EdbSegP(*((EdbSegP *) segments->At(0)));
379 s->PropagateTo(s->Z()-650);
380 xx=s->X();
381 yy=s->Y();
382 zz=s->Z();
383 delete s;
384 }
385 else {
386 for (int i=0; i<segments->GetEntries(); i++) {
387 EdbSegP *s = (EdbSegP *) segments->At(i);
388 double ax = s->TX();
389 double ay = s->TY();
390 double az = 1.0;
391 double x = s->X();
392 double y = s->Y();
393 double z = s->Z();
394 double a = ax*ax+ay*ay+az*az;
395 double c = -ax*x-ay*y-az*z;
396 double b = (ax*ax+ay*ay);
397 // double w = 1.0/a/a; // weight for small angle tracks.
398 double w = 1.0; // no weight
399
400 Ax+=2.0*w/b*( a*(ay*ay+az*az) );
401 Bx+=2.0*w/b*( -a*ax*ay );
402 Cx+=2.0*w/b*( -a*ax*az );
403 Dx+=2.0*w/b*( -(a*x+c*ax)*(ax*ax-a)-(a*y+c*ay)*ax*ay-(a*z+c*az)*az*ax );
404
405 Ay+=2.0*w/b*( -a*ay*ax );
406 By+=2.0*w/b*( a*(az*az+ax*ax) );
407 Cy+=2.0*w/b*( -a*ay*az );
408 Dy+=2.0*w/b*( -(a*y+c*ay)*(ay*ay-a)-(a*z+c*az)*ay*az-(a*x+c*ax)*ax*ay );
409
410 Az+=2.0*w/b*( -a*az*ax );
411 Bz+=2.0*w/b*( -a*az*ay );
412 Cz+=2.0*w/b*( a*b );
413 Dz+=2.0*w/b*( -(a*z+c*az)*(az*az-a)-(a*x+c*ax)*az*ax-(a*y+c*ay)*ay*az );
414
415 }
416
417 Det=fabs( Ax*(By*Cz-Cy*Bz)-Bx*(Ay*Cz-Cy*Az)+Cx*(Ay*Bz-By*Az) );
418 xx=( (By*Cz-Cy*Bz)*Dx-(Bx*Cz-Cx*Bz)*Dy+(Bx*Cy-Cx*By)*Dz)/Det;
419 yy=(-(Ay*Cz-Cy*Az)*Dx+(Ax*Cz-Cx*Az)*Dy-(Ax*Cy-Cx*Ay)*Dz)/Det;
420 zz=( (Ay*Bz-By*Az)*Dx-(Ax*Bz-Bx*Az)*Dy+(Ax*By-Bx*Ay)*Dz)/Det;
421 }
422
423 EdbVertex *v = new EdbVertex();
424 v->SetXYZ(xx,yy,zz);
425
426 return v;
427}
Expr< UnaryOp< Fabs< T >, Expr< A, T, D >, T >, T, D > fabs(const Expr< A, T, D > &rhs)
Definition: UnaryOperators.hh:96
void a()
Definition: check_aligned.C:59
Definition: EdbSegP.h:21
Definition: EdbVertex.h:69
void SetXYZ(float x, float y, float z)
Definition: EdbVertex.h:157
s
Definition: check_shower.C:55
void w(int rid=2, int nviews=2)
Definition: test.C:27

◆ ClassDef()

EdbShowerAlg::ClassDef ( EdbShowerAlg  ,
 
)

◆ Convert_EdbPVRec_To_InBTArray()

void EdbShowerAlg::Convert_EdbPVRec_To_InBTArray ( )

This function takes all basetracks of the EdbPVRec object,
and fills the InBTArray with those.
This function is useful, if one wants to start from all
possible basetracks.

1117{
1122
1123 if (eInBTArray!=NULL) {
1124 eInBTArray->Clear();
1125 cout << "EdbShowerAlg::Convert_EdbPVRec_To_InBTArray() eInBTArray->Clear() done." << endl;
1126 }
1127 if (eInBTArray==NULL) {
1128 // (Should never happen.)
1129 eInBTArray = new TObjArray();
1130 cout << " eInBTArray = new TObjArray()" << endl;
1131 }
1132
1133 cout << "EdbShowerAlg::Convert_EdbPVRec_To_InBTArray() eInBTArray at address "<< eInBTArray <<endl;
1134
1135
1136 eInBTArray->Clear();
1137 Int_t npat=eAli->Npatterns();
1138 EdbSegP* Segment2;
1139 EdbPattern* pat_one;
1140 for (Int_t j=0; j<npat; ++j) {
1141 pat_one=(EdbPattern*)eAli->GetPattern(j);
1142 for (Int_t i=0; i< pat_one->N(); i++) {
1143 Segment2 = (EdbSegP*)pat_one->GetSegment(i);
1144 eInBTArray->Add(Segment2);
1145 }
1146 }
1147 eInBTArrayN=eInBTArray->GetEntries();
1148
1149 cout << "Convert_EdbPVRec_To_InBTArray eInBTArray->GetEntries() ";
1150 cout << eInBTArray->GetEntries() << endl;
1151 cout << "EdbShowerAlg::Convert_EdbPVRec_To_InBTArray()...done."<<endl;
1152 return;
1153}
Int_t npat
Definition: Xi2HatStartScript.C:33
Definition: EdbPattern.h:273
Int_t Npatterns() const
Definition: EdbPattern.h:366
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
EdbPVRec * eAli
Definition: EdbShowerAlg.h:48
Int_t eInBTArrayN
Definition: EdbShowerAlg.h:52
TObjArray * eInBTArray
Definition: EdbShowerAlg.h:51

◆ DeltaR_NoPropagation()

Double_t EdbShowerAlg::DeltaR_NoPropagation ( EdbSegP s,
EdbSegP stest 
)

Calculate the distance between two basetracks by
sqrt((deltaX)^2+(deltaY)^2)
Do NOT propagate tracks onto same Z, in case they are different.
SAME function as DeltaR_WithoutPropagation !!!

557{
562 return DeltaR_WithoutPropagation(s,stest);
563}
Double_t DeltaR_WithoutPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowerAlg.cxx:567

◆ DeltaR_WithoutPropagation()

Double_t EdbShowerAlg::DeltaR_WithoutPropagation ( EdbSegP s,
EdbSegP stest 
)

Calculate the distance between two basetracks by
sqrt((deltaX)^2+(deltaY)^2)
Do NOT propagate tracks onto same Z, in case they are different.

568{
572 return TMath::Sqrt((s->X()-stest->X())*(s->X()-stest->X())+(s->Y()-stest->Y())*(s->Y()-stest->Y()));
573}
Float_t X() const
Definition: EdbSegP.h:173
Float_t Y() const
Definition: EdbSegP.h:174

◆ DeltaR_WithPropagation()

Double_t EdbShowerAlg::DeltaR_WithPropagation ( EdbSegP s,
EdbSegP stest 
)

Calculate the distance between two basetracks by
sqrt((deltaX)^2+(deltaY)^2)
DO propagate tracks onto same Z, in case they are different.
(propagation is done from stest onto Z-position of s)

578{
583 if (s->Z()==stest->Z()) return TMath::Sqrt((s->X()-stest->X())*(s->X()-stest->X())+(s->Y()-stest->Y())*(s->Y()-stest->Y()));
584 Double_t zorig;
585 Double_t dR;
586 zorig=s->Z();
587 s->PropagateTo(stest->Z());
588 dR=TMath::Sqrt( (s->X()-stest->X())*(s->X()-stest->X())+(s->Y()-stest->Y())*(s->Y()-stest->Y()) );
589 s->PropagateTo(zorig);
590 return dR;
591}
Float_t Z() const
Definition: EdbSegP.h:153

◆ DeltaTheta()

Double_t EdbShowerAlg::DeltaTheta ( EdbSegP s1,
EdbSegP s2 
)

Be aware that this DeltaTheta function returns the abs() difference between the
ABSOLUTE values of dTheta!!! (not componentwise!

596{
599 Double_t tx1,tx2,ty1,ty2;
600 tx1=s1->TX();
601 tx2=s2->TX();
602 ty1=s1->TY();
603 ty2=s2->TY();
604 Double_t dt= TMath::Abs(TMath::Sqrt(tx1*tx1+ty1*ty1) - TMath::Sqrt(tx2*tx2+ty2*ty2));
605 return dt;
606}
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ DeltaThetaComponentwise()

Double_t EdbShowerAlg::DeltaThetaComponentwise ( EdbSegP s1,
EdbSegP s2 
)

Be aware that this DeltaTheta function returns the difference between the component values of dTheta!!!
Acutally this function should be the normal way to calculate dTheta correctly...

612{
616 Double_t tx1,tx2,ty1,ty2;
617 tx1=s1->TX();
618 tx2=s2->TX();
619 ty1=s1->TY();
620 ty2=s2->TY();
621 Double_t dt= TMath::Sqrt( (tx1-tx2)*(tx1-tx2) + (ty1-ty2)*(ty1-ty2) );
622 return dt;
623}

◆ DeltaThetaSingleAngles()

Double_t EdbShowerAlg::DeltaThetaSingleAngles ( EdbSegP s1,
EdbSegP s2 
)

SAME function as DeltaThetaComponentwise !!!

626{
629}
Double_t DeltaThetaComponentwise(EdbSegP *s1, EdbSegP *s2)
Definition: EdbShowerAlg.cxx:611

◆ Execute()

void EdbShowerAlg::Execute ( )
virtual

Reimplemented in EdbShowerAlg_GS.

440{
441 cout << "EdbShowerAlg::Execute() -----------------------------------------------------------" << endl;
442
443 cout << "EdbShowerAlg::Execute() This function should be implemented in the inheriting classes!!" << endl;
444 return;
445}

◆ Finalize()

void EdbShowerAlg::Finalize ( )
virtual

Reimplemented in EdbShowerAlg_GS.

450{
451
452 return;
453}

◆ GetAlgName()

TString EdbShowerAlg::GetAlgName ( ) const
inline
126 {
127 return eAlgName;
128 }

◆ GetAlgValue()

Int_t EdbShowerAlg::GetAlgValue ( ) const
inline
123 {
124 return eAlgValue;
125 }

◆ GetMinimumDist()

Double_t EdbShowerAlg::GetMinimumDist ( EdbSegP seg1,
EdbSegP seg2 
)

Calculate minimum distance of 2 lines.
Use the data of (the selected object)->X(), Y(), Z(), TX(), TY().
means, if the selected object == segment, use the data of the segment. or it == track, the use the fitted data.
Original code from Tomoko Ariga
(double calc_dmin(EdbSegP *seg1, EdbSegP *seg2, double *dminz = NULL))

652{
658
659 double x1,y1,z1,ax1,ay1;
660 double x2,y2,z2,ax2,ay2;
661 double s1,s2,s1bunsi,s1bunbo,s2bunsi,s2bunbo;
662 double p1x,p1y,p1z,p2x,p2y,p2z,p1p2;
663
664 if (seg1->ID()==seg2->ID()&&seg1->PID()==seg2->PID()) return 0.0;
665
666 x1 = seg1->X();
667 y1 = seg1->Y();
668 z1 = seg1->Z();
669 ax1= seg1->TX();
670 ay1= seg1->TY();
671
672 x2 = seg2->X();
673 y2 = seg2->Y();
674 z2 = seg2->Z();
675 ax2= seg2->TX();
676 ay2= seg2->TY();
677
678 s1bunsi=(ax2*ax2+ay2*ay2+1)*(ax1*(x2-x1)+ay1*(y2-y1)+z2-z1) - (ax1*ax2+ay1*ay2+1)*(ax2*(x2-x1)+ay2*(y2-y1)+z2-z1);
679 s1bunbo=(ax1*ax1+ay1*ay1+1)*(ax2*ax2+ay2*ay2+1) - (ax1*ax2+ay1*ay2+1)*(ax1*ax2+ay1*ay2+1);
680 s2bunsi=(ax1*ax2+ay1*ay2+1)*(ax1*(x2-x1)+ay1*(y2-y1)+z2-z1) - (ax1*ax1+ay1*ay1+1)*(ax2*(x2-x1)+ay2*(y2-y1)+z2-z1);
681 s2bunbo=(ax1*ax1+ay1*ay1+1)*(ax2*ax2+ay2*ay2+1) - (ax1*ax2+ay1*ay2+1)*(ax1*ax2+ay1*ay2+1);
682 s1=s1bunsi/s1bunbo;
683 s2=s2bunsi/s2bunbo;
684 p1x=x1+s1*ax1;
685 p1y=y1+s1*ay1;
686 p1z=z1+s1*1;
687 p2x=x2+s2*ax2;
688 p2y=y2+s2*ay2;
689 p2z=z2+s2*1;
690 p1p2=sqrt( (p1x-p2x)*(p1x-p2x)+(p1y-p2y)*(p1y-p2y)+(p1z-p2z)*(p1z-p2z) );
691
692 return p1p2;
693}
Int_t ID() const
Definition: EdbSegP.h:147
Int_t PID() const
Definition: EdbSegP.h:148

◆ GetRecoShowerArray()

TObjArray * EdbShowerAlg::GetRecoShowerArray ( ) const
inline
132 {
133 return eRecoShowerArray;
134 }

◆ GetRecoShowerArrayN()

Int_t EdbShowerAlg::GetRecoShowerArrayN ( ) const
inline
129 {
130 return eRecoShowerArrayN;
131 }

◆ GetSpatialDist()

Double_t EdbShowerAlg::GetSpatialDist ( EdbSegP s1,
EdbSegP s2 
)

Mainly Z values should dominate...
since the are at the order of 10k microns and x,y of 1k microns

634{
637 Double_t x1,x2,y1,y2,z1,z2;
638 x1=s1->X();
639 x2=s2->X();
640 y1=s1->Y();
641 y2=s2->Y();
642 z1=s1->Z();
643 z2=s2->Z();
644 Double_t dist= TMath::Sqrt( (x1-x2)*(x1-x2) + (y1-y2)*(y1-y2) + (z1-z2)*(z1-z2) );
645 //cout << "dist = " << dist << endl;
646 return dist;
647}

◆ Initialize()

void EdbShowerAlg::Initialize ( )
virtual

Reimplemented in EdbShowerAlg_GS.

433{
434
435 return;
436}

◆ IsInConeTube()

Bool_t EdbShowerAlg::IsInConeTube ( EdbSegP sTest,
EdbSegP sStart,
Double_t  CylinderRadius,
Double_t  ConeAngle 
)

General Function which returns Bool if the Testing BaeTrack is in a cone defined
by the StartingBaseTrack. In case of starting same Z position, a distance cut of
20microns is assumed....
In case of TestingSegment==StartingSegment this function should correctly return kTRUE also...

Outside if angle greater than ConeAngle (to be fulfilled for Cone and Tube in both cases)

if angle smaller than ConeAngle, then you can differ between Tuberadius and CylinderRadius

286{
291 if (gEDBDEBUGLEVEL>3) cout << "Bool_t EdbShowerAlg::IsInConeTube() Test Segment " << TestingSegment << " vs. Starting Segment " << StartingSegment << endl;
292
293 // We reject any TestingSegment segments which have lower Z than the StartingSegment .
294 if (StartingSegment->Z()>TestingSegment->Z() ) return kFALSE;
295
296 TVector3 x1(StartingSegment->X(),StartingSegment->Y(),StartingSegment->Z());
297 TVector3 x2(TestingSegment->X(),TestingSegment->Y(),TestingSegment->Z());
298 TVector3 direction_x1(StartingSegment->TX()*1300,StartingSegment->TY()*1300,1300);
299
300 // u1 is the difference vector of the position!
301 TVector3 u1=x2-x1;
302
303 Double_t direction_x1_norm= direction_x1.Mag();
304 Double_t cosangle= (direction_x1*u1)/(u1.Mag()*direction_x1_norm);
305 Double_t angle = TMath::ACos(cosangle);
306
307 // This is the old version of angle calculation. It does not give the same results as in ROOT
308 // when use TVector3.Angle(&TVector3). // For this IsInConeTube() we use therefore the ROOT calculation.
309 angle=u1.Angle(direction_x1);
310
311 // For the case where the two basetracks have same z position
312 // the angle is about 90 degree so it makes no sense to calculate it...
313 // therefore we set it artificially to zero:
314 if (StartingSegment->Z()==TestingSegment->Z() ) {
315 if (gEDBDEBUGLEVEL>3) cout << "same Z position of TestingSegment and StartingSegment, Set angle artificially to zero" << endl;
316 angle=0.0;
317 // Check here for dR manually:
318 //cout << DeltaR_WithoutPropagation(StartingSegment,TestingSegment) << endl;
319 //StartingSegment->PrintNice();
320 //TestingSegment->PrintNice();
321// cout << StartingSegment->Flag() << " " << TestingSegment->Flag() << endl;
322// cout << StartingSegment->P() << " " << TestingSegment->P() << endl;
323
324 // Check for position distance for 20microns if
325 // Testing Segment is in same Z as StartingSegment
326 if (gEDBDEBUGLEVEL>3) cout << "Check for position distance for 20microns if Testing Segment is in same Z as StartingSegment" << endl;
327 if (gEDBDEBUGLEVEL>3) cout << "DeltaR_WithoutPropagation(StartingSegment,TestingSegment) = "<< DeltaR_WithoutPropagation(StartingSegment,TestingSegment) << endl;
328 if (DeltaR_WithoutPropagation(StartingSegment,TestingSegment)<20) return kTRUE;
329 if (DeltaR_WithoutPropagation(StartingSegment,TestingSegment)>=20) return kFALSE;
330 }
331
333 if (gEDBDEBUGLEVEL>3) cout << "Check if AngleVector now within the ConeAngleVector (<"<< ConeAngle<<"): " << angle << endl;
334 if (angle>ConeAngle) {
335 return kFALSE;
336 }
337
339 Double_t TubeDistance = 1.0/direction_x1_norm * ( (x2-x1).Cross(direction_x1) ).Mag();
340
341 if (gEDBDEBUGLEVEL>3) cout << "Check if TestingSegment is now within the Tube (<"<< CylinderRadius<<"): " << TubeDistance << endl;
342
343 if (TubeDistance>CylinderRadius) {
344 return kFALSE;
345 }
346
347 return kTRUE;
348}
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ Print()

void EdbShowerAlg::Print ( )
458{
459 cout << "EdbShowerAlg::Print()" << endl;
460 cout << eAlgName << " ; AlgValue= " << eAlgName << " ." << endl;
461 for (int i=0; i<10; i++) cout << eParaString[i] << "=" << eParaValue[i] << "; ";
462 cout << endl;
463 cout << "PID numbers: ..." << eFirstPlate_eAliPID << " " << eLastPlate_eAliPID << " " << eMiddlePlate_eAliPID << " " << eNumberPlate_eAliPID << " " << endl;
464 cout << "UseAliSub= " << eUseAliSub << " ." << endl;
465 cout << "EdbShowerAlg::Print()...done." << endl;
466 return;
467}
Int_t eNumberPlate_eAliPID
Definition: EdbShowerAlg.h:57
Int_t eMiddlePlate_eAliPID
Definition: EdbShowerAlg.h:56
Int_t eLastPlate_eAliPID
Definition: EdbShowerAlg.h:55
Int_t eUseAliSub
Definition: EdbShowerAlg.h:67
Int_t eFirstPlate_eAliPID
Definition: EdbShowerAlg.h:54

◆ PrintAll()

void EdbShowerAlg::PrintAll ( )
542{
543 cout << "------------------------------------------------------------" << endl;
544 cout << "EdbShowerAlg::PrintAll()" << endl;
545 cout << "------------------------" << endl;
546 Print();
547 PrintMore();
548 cout << "------------------------------------------------------------" << endl;
549 return;
550}
void PrintMore()
Definition: EdbShowerAlg.cxx:494
void Print()
Definition: EdbShowerAlg.cxx:457

◆ PrintMore()

void EdbShowerAlg::PrintMore ( )
495{
496 cout << "EdbShowerAlg::PrintMore()" << endl;
497 cout << "eInBTArray->GetEntries() = " << eInBTArray->GetEntries() << " ." << endl;
498 cout << "eRecoShowerArray->GetEntries() = " << eRecoShowerArray->GetEntries() << " ." << endl;
499 cout << "EdbShowerAlg::PrintMore()...done." << endl;
500 return;
501}

◆ PrintParameters()

void EdbShowerAlg::PrintParameters ( )
472{
473 cout << "EdbShowerAlg::PrintParameters()" << endl;
474 cout << eAlgName<< " :" << endl;
475 for (int i=0; i<5; i++) cout << setw(6) << eParaString[i];
476 cout << endl;
477 for (int i=0; i<5; i++) cout << setw(6) << eParaValue[i];
478 cout << " ."<<endl;
479 return;
480}

◆ PrintParametersShort()

void EdbShowerAlg::PrintParametersShort ( )
485{
486 cout << eAlgName<< " :";
487 for (int i=0; i<5; i++) cout << setw(6) << eParaValue[i];
488 cout << " ."<<endl;
489 return;
490}

◆ PrintRecoShowerArray()

void EdbShowerAlg::PrintRecoShowerArray ( )
506{
507 cout << "EdbShowerAlg::PrintRecoShowerArray()..." << endl;
508 cout << "EdbShowerAlg::PrintRecoShowerArray() eRecoShowerArray->GetEntries() = " << eRecoShowerArray->GetEntries() << " ." << endl;
509
510 if (eRecoShowerArrayN<1) return;
511 printf("i X Y Z NBT \n");
512
513 Bool_t caseMC=kFALSE;
514
515// case a) MC event info there:
517 EdbSegP* seg= tr->GetSegment(0);
518
519 if (seg->MCEvt()>0) caseMC =kTRUE;
520
521 for (Int_t i=0; i<eRecoShowerArray->GetEntries(); ++i) {
522 if (caseMC==kTRUE) {
523
525 printf("%06d %06f %06f %06f %06d \n", i, tr->X(), tr->Y(), tr->Z(), tr->N());
526 }
527 else {
529 printf("%06d %06f %06f %06f %06d \n", i, tr->X(), tr->Y(), tr->Z(), tr->N());
530 }
531 }
532
533// case a) MC event info there:
534 cout << "EdbShowerAlg::PrintRecoShowerArray()...done." << endl;
535 return;
536}
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
Int_t MCEvt() const
Definition: EdbSegP.h:145

◆ Set0()

void EdbShowerAlg::Set0 ( )
protected

Reset internal variable to default values.

54{
55
57
58 cout << "EdbShowerAlg::Set0()"<<endl;
59
60 eAlgName="UNSPECIFIED";
61 eAlgValue=-999999;
62 for (int i=0; i<10; ++i) {
63 eParaValue[i]=-99999.0;
64 eParaString[i]="UNSPECIFIED";
65 }
66 eAli_Sub=0;
67
68 // do not use use small eAli Object by default:
69 // (this solution is memory safe....)
70 eUseAliSub=0;
71
73
78}
Int_t eActualAlgParametersetNr
Used when more sets of same algorithm.
Definition: EdbShowerAlg.h:43

◆ SetActualAlgParameterset()

void EdbShowerAlg::SetActualAlgParameterset ( Int_t  ActualAlgParametersetNr)
inline
119 {
120 eActualAlgParametersetNr=ActualAlgParametersetNr;
121 }

◆ SetDebug()

void EdbShowerAlg::SetDebug ( )
inline
170 {
171 eDebug=kTRUE;
172 }
Bool_t eDebug
Definition: EdbShowerAlg.h:80

◆ SetEdbPVRec()

void EdbShowerAlg::SetEdbPVRec ( EdbPVRec Ali)
inline
95 {
96 eAli = Ali;
98 }
Int_t eAliNpat
Definition: EdbShowerAlg.h:49

◆ SetEdbPVRecPIDNumbers()

void EdbShowerAlg::SetEdbPVRecPIDNumbers ( Int_t  FirstPlate_eAliPID,
Int_t  LastPlate_eAliPID,
Int_t  MiddlePlate_eAliPID,
Int_t  NumberPlate_eAliPID 
)
inline
103 {
104 eFirstPlate_eAliPID=FirstPlate_eAliPID;
105 eLastPlate_eAliPID=LastPlate_eAliPID;
106 eMiddlePlate_eAliPID=MiddlePlate_eAliPID;
107 eNumberPlate_eAliPID=NumberPlate_eAliPID;
108 }

◆ SetInBTArray()

void EdbShowerAlg::SetInBTArray ( TObjArray *  InBTArray)
inline
99 {
100 eInBTArray = InBTArray;
101 eInBTArrayN=eInBTArray->GetEntries();
102 }

◆ SetParameter()

void EdbShowerAlg::SetParameter ( Int_t  parNr,
Float_t  parvalue 
)
83{
84 if (parNr>9) {
85 cout << "EdbShowerAlg::SetParameter() WARNING parNr>9 . Do nothing."<<endl;
86 return;
87 }
88 eParaValue[parNr]=parvalue;
89 cout << "EdbShowerAlg::SetParameter()...done."<<endl;
90}

◆ SetParameters()

void EdbShowerAlg::SetParameters ( Float_t *  par)
95{
96 cout << "EdbShowerAlg::SetParameters()..."<<endl;
97 // SetParameters
98 for (int i=0; i<10; ++i) {
99 eParaValue[i]=par[i];
100 cout << "EdbShowerAlg::SetParameters()...Parameter " << i << " set to: " << eParaValue[i] <<endl;
101 }
102 cout << "EdbShowerAlg::SetParameters()...done."<<endl;
103}

◆ SetRecoShowerArray()

void EdbShowerAlg::SetRecoShowerArray ( TObjArray *  RecoShowerArray)
inline

◆ SetRecoShowerArrayN()

void EdbShowerAlg::SetRecoShowerArrayN ( Int_t  RecoShowerArrayN)
inline
113 {
114 eRecoShowerArrayN=RecoShowerArrayN;
115 }

◆ SetUseAliSub()

void EdbShowerAlg::SetUseAliSub ( Bool_t  UseAliSub)
inline
136 {
137 eUseAliSub=UseAliSub;
138 }

◆ Transform_eAli()

void EdbShowerAlg::Transform_eAli ( EdbSegP InitiatorBT,
Float_t  ExtractSize = 1500 
)

-----------------------------------------------------------------------------------—
Transform eAli to eAli_Sub:
the lenght of eAli_Sub is not changed.
Only XY-size (and MC) cuts are applied.

-----------------------------------------------------------------------------------— \n

— Whereas in ShowRec.cpp the treebranch file is written directly after each
— BT reconstruction, it was not a problem when the eAliSub was deleted each time.
— But now situation is different, since the BTs of the showers have the adresses from
— the eAli_Sub !each! so if one eAli_Sub is deleted the EdbShowerP object has lost its
— BT adresses.
— So if we do not delete the eAli_Sub, reconstrtuction is fast (compared to eAli) but
— memory increases VERY fast.
— If we do use eAli, then memory consumption will not increase fast, but
— reconstruction is slow.

— Possible workarounds:
— For few InBTs (EDA use for data...): use eAliSub and dont delete it.
— For many InBTs (parameterstudies): use eAli—.
— Use always eAliSub, search BT correspond in eAli, add BT from eAli....May take also long time...


— SEVERE warning: IF gAli is in wrong order it can be that no showers
— SEVERE warning: are reconstructed since most implemented algorithms
— SEVERE warning: rely on the order that plate 2 comes directly behind

— SEVERE warning: the InBT. THIS IS STILL ON BUGFIXING LIST!!! \n

DEBUG if (eAli_Sub) { delete eAli_Sub;eAli_Sub=0;} // original, but keeps not adresses of segment in eAli.

do nothing now... let it live... delete eAli_Sub;eAli_Sub=0;}

if (eAli_Sub) { delete eAli_Sub;eAli_Sub=0;} Try not to delete it maybe then it works.....

138{
166
167
168 /*
169 cout << "void EdbShowerAlg::Transform_eAli() SEVERE WARNING: IF gAli is in wrong order it can be that no showers " << endl;
170 cout << "void EdbShowerAlg::Transform_eAli() SEVERE WARNING: are reconstructed since most implemented algorithms " << endl;
171 cout << "void EdbShowerAlg::Transform_eAli() SEVERE WARNING: rely on the order that plate 2 comes directly behind " << endl;
172 cout << "void EdbShowerAlg::Transform_eAli() SEVERE WARNING: the InBT. /// DEBUG TODO // DEBUG TODO " << endl;
173 */
174
175
176// gEDBDEBUGLEVEL=3;
177// cout << "void EdbShowerAlg::Transform_eAli() ===== eUseAliSub " << eUseAliSub << endl;
178
179 // IF TO RECREATE THE gALI_SUB in RECONSTRUCTION or if to use gAli global (slowlier but maybe memory leak safe).
180 if (!eUseAliSub)
181 {
182 if (gEDBDEBUGLEVEL>2) cout << "EdbShowerAlg::Transform_eAli UseAliSub==kFALSE No new eAli_Sub created. Use eAli instead. "<<endl;
184 eAli_SubNpat=eAli_Sub->Npatterns(); //number of plates
185 if (gEDBDEBUGLEVEL>3) eAli_Sub->Print();
186 return;
187 }
188
189 Int_t npat;
190 npat = eAli->Npatterns(); //number of plates
191
192 // has to be deleted in some part of the script outside this function...
193 // Dont forget , otherwise memory heap overflow!
195 if (eAli_Sub) {
196 ;
197 }
198// eAli_Sub = new EdbPVRec();
199
200
201 if (eUseAliSub) {
202 if (gEDBDEBUGLEVEL>2) cout << "EdbShowerAlg::Transform_eAli UseAliSub==kTRUE Will now create new eAli_Sub. "<<endl;
203
205 eAli_Sub = new EdbPVRec();
206// cout << "Adress of eAli_Sub " << eAli_Sub << endl;
207// cout << "eAli_Sub->Npatterns(); " << eAli_Sub->Npatterns() << endl;
208// eAli_Sub->Print();
209 }
210
211
212
213// cout << "PROBLE THRANFOMING, CAUSE ADRESES ARE DELETED ALSO...." << endl;
214// cout << "TEMPORAY SOLUTION: comment the delete eAli_Sub ( will lead to large memory consumption when run for long time" << endl;
215
216
217 // Create SubPattern objects
218 EdbSegP* ExtrapolateInitiatorBT=0;
219 ExtrapolateInitiatorBT = (EdbSegP*)InitiatorBT->Clone();
220
221 Int_t InitiatorBTMCEvt=InitiatorBT->MCEvt();
222
223 // Create Variables For ExtractSubpattern boundaries
224 Float_t mini[5];
225 Float_t maxi[5];
226 //Float_t ExtractSize=1000;// now in fucntion header
227 mini[0]=ExtrapolateInitiatorBT->X()-ExtractSize;
228 mini[1]=ExtrapolateInitiatorBT->Y()-ExtractSize;
229 maxi[0]=ExtrapolateInitiatorBT->X()+ExtractSize;
230 maxi[1]=ExtrapolateInitiatorBT->Y()+ExtractSize;
231 mini[2]=-0.7;
232 mini[3]=-0.7;
233 mini[4]=0.0;
234 maxi[2]=0.7;
235 maxi[3]=0.7;
236 maxi[4]=100.0;
237
238 EdbPattern* singlePattern;
239 Float_t ExtrapolateInitiatorBT_zpos_orig=ExtrapolateInitiatorBT->Z();
240
241 // Add the subpatterns in a loop for the plates:
242 // in reverse ordering.due to donwstream behaviour (!):
243 // (Only downstream is supported now...)
244 for (Int_t ii=0; ii<npat; ++ii) {
245
246 Float_t zpos=eAli->GetPattern(ii)->Z();
247
248 ExtrapolateInitiatorBT->PropagateTo(zpos);
249
250 mini[0]=ExtrapolateInitiatorBT->X()-ExtractSize;
251 mini[1]=ExtrapolateInitiatorBT->Y()-ExtractSize;
252 maxi[0]=ExtrapolateInitiatorBT->X()+ExtractSize;
253 maxi[1]=ExtrapolateInitiatorBT->Y()+ExtractSize;
254
255 singlePattern=(EdbPattern*)eAli->GetPattern(ii)->ExtractSubPattern(mini,maxi,InitiatorBTMCEvt);
256 // This sets PID() analogue to (upstream), nut not PID of the BTs !
257 singlePattern-> SetID(eAli->GetPattern(ii)->ID());
258 // This sets PID() analogue to (upstream), nut not PID of the BTs !
259 singlePattern-> SetPID(eAli->GetPattern(ii)->PID());
260
261 eAli_Sub->AddPattern(singlePattern);
262
263 // Propagate back...! (in order not to change the original BT)
264 ExtrapolateInitiatorBT->PropagateTo(ExtrapolateInitiatorBT_zpos_orig);
265 }
266
267 delete ExtrapolateInitiatorBT;
268
269 eAli_SubNpat=eAli_Sub->Npatterns(); //number of plates
270
271
272 if (gEDBDEBUGLEVEL>2) {
273 cout << "EdbShowerAlg::Transform_eAli eAli_Sub created."<<endl;
274 cout << "Adress of eAli_Sub " << eAli_Sub << endl;
275 cout << "eAli_Sub->Npatterns(); " << eAli_Sub->Npatterns() << endl;
276 }
277 return;
278}
Definition: EdbPVRec.h:148
int PID() const
Definition: EdbPattern.h:320
int ID() const
Definition: EdbPattern.h:319
EdbPattern * ExtractSubPattern(float min[5], float max[5], int MCevt=-1)
Definition: EdbPattern.cxx:1470
void AddPattern(EdbPattern *pat)
Definition: EdbPattern.cxx:1707
void Print() const
Definition: EdbPattern.cxx:1693
void PropagateTo(float z)
Definition: EdbSegP.cxx:293
Float_t Z() const
Definition: EdbPattern.h:84
Int_t eAli_SubNpat
Definition: EdbShowerAlg.h:64

Member Data Documentation

◆ eActualAlgParametersetNr

Int_t EdbShowerAlg::eActualAlgParametersetNr
protected

Used when more sets of same algorithm.

◆ eAlgName

TString EdbShowerAlg::eAlgName
protected

◆ eAlgValue

Int_t EdbShowerAlg::eAlgValue
protected

◆ eAli

EdbPVRec* EdbShowerAlg::eAli
protected

◆ eAli_Sub

EdbPVRec* EdbShowerAlg::eAli_Sub
protected

◆ eAli_SubNpat

Int_t EdbShowerAlg::eAli_SubNpat
protected

◆ eAliNpat

Int_t EdbShowerAlg::eAliNpat
protected

◆ eDebug

Bool_t EdbShowerAlg::eDebug
protected

◆ eFirstPlate_eAliPID

Int_t EdbShowerAlg::eFirstPlate_eAliPID
protected

◆ eInBTArray

TObjArray* EdbShowerAlg::eInBTArray
protected

◆ eInBTArrayN

Int_t EdbShowerAlg::eInBTArrayN
protected

◆ eLastPlate_eAliPID

Int_t EdbShowerAlg::eLastPlate_eAliPID
protected

◆ eMiddlePlate_eAliPID

Int_t EdbShowerAlg::eMiddlePlate_eAliPID
protected

◆ eNumberPlate_eAliPID

Int_t EdbShowerAlg::eNumberPlate_eAliPID
protected

◆ eParaString

TString EdbShowerAlg::eParaString[10]
protected

◆ eParaValue

Float_t EdbShowerAlg::eParaValue[10]
protected

◆ eRecoShower

EdbTrackP* EdbShowerAlg::eRecoShower
protected

◆ eRecoShowerArray

TObjArray* EdbShowerAlg::eRecoShowerArray
protected

◆ eRecoShowerArrayN

Int_t EdbShowerAlg::eRecoShowerArrayN
protected

◆ eUseAliSub

Int_t EdbShowerAlg::eUseAliSub
protected

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