FEDRA emulsion software from the OPERA Collaboration
EdbShowAlg Class Reference

#include <EdbShowAlg.h>

Inheritance diagram for EdbShowAlg:
Collaboration diagram for EdbShowAlg:

Public Member Functions

void AddRecoShowerArray (EdbShowerP *shower)
 
 ClassDef (EdbShowAlg, 1)
 
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)
 
 EdbShowAlg ()
 
 EdbShowAlg (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
 
EdbShowerPGetShower (Int_t i) const
 
Double_t GetSpatialDist (EdbSegP *s1, EdbSegP *s2)
 
void Help ()
 
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 SetActualAlgParameterset (Int_t ActualAlgParametersetNr)
 
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 par)
 
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)
 
void UpdateShowerIDs ()
 
void UpdateShowerMetaData ()
 
virtual ~EdbShowAlg ()
 

Protected Member Functions

void Set0 ()
 

Protected Attributes

Int_t eActualAlgParametersetNr
 
TString eAlgName
 
Int_t eAlgValue
 
EdbPVReceAli
 
EdbPVReceAli_Sub
 
Int_t eAli_SubNpat
 
Int_t eAliNpat
 
Int_t eFirstPlate_eAliPID
 
TObjArray * eInBTArray
 
Int_t eInBTArrayN
 
Int_t eLastPlate_eAliPID
 
Int_t eMiddlePlate_eAliPID
 
Int_t eNumberPlate_eAliPID
 
Int_t eParaN
 
TString eParaString [10]
 
Float_t eParaValue [10]
 
EdbShowerPeRecoShower
 
TObjArray * eRecoShowerArray
 
Int_t eRecoShowerArrayN
 
Int_t eUseAliSub
 

Constructor & Destructor Documentation

◆ EdbShowAlg() [1/2]

EdbShowAlg::EdbShowAlg ( )

◆ EdbShowAlg() [2/2]

EdbShowAlg::EdbShowAlg ( TString  AlgName,
Int_t  AlgValue 
)
29{
30// Reset all:
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}
Float_t eParaValue[10]
Definition: EdbShowAlg.h:52
TString eParaString[10]
Definition: EdbShowAlg.h:53
TString eAlgName
Definition: EdbShowAlg.h:49
void Set0()
Definition: EdbShowAlg.cxx:53
Int_t eAlgValue
Definition: EdbShowAlg.h:50
EdbPVRec * eAli_Sub
Definition: EdbShowAlg.h:74

◆ ~EdbShowAlg()

EdbShowAlg::~EdbShowAlg ( )
virtual
45{
46// Default Destructor
47 cout << "EdbShowAlg::~EdbShowAlg()"<<endl;
48}

Member Function Documentation

◆ AddRecoShowerArray()

void EdbShowAlg::AddRecoShowerArray ( EdbShowerP shower)
inline
120 {
121 eRecoShowerArray->Add(shower);
123 }
TObjArray * eRecoShowerArray
Definition: EdbShowAlg.h:70
Int_t eRecoShowerArrayN
Definition: EdbShowAlg.h:71

◆ ClassDef()

EdbShowAlg::ClassDef ( EdbShowAlg  ,
 
)

◆ DeltaR_NoPropagation()

Double_t EdbShowAlg::DeltaR_NoPropagation ( EdbSegP s,
EdbSegP stest 
)
397{
398// SAME function as DeltaR_WithoutPropagation !!!
399 return DeltaR_WithoutPropagation(s,stest);
400}
Double_t DeltaR_WithoutPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowAlg.cxx:404
s
Definition: check_shower.C:55

◆ DeltaR_WithoutPropagation()

Double_t EdbShowAlg::DeltaR_WithoutPropagation ( EdbSegP s,
EdbSegP stest 
)
405{
406 return TMath::Sqrt((s->X()-stest->X())*(s->X()-stest->X())+(s->Y()-stest->Y())*(s->Y()-stest->Y()));
407}
Float_t X() const
Definition: EdbSegP.h:173
Float_t Y() const
Definition: EdbSegP.h:174

◆ DeltaR_WithPropagation()

Double_t EdbShowAlg::DeltaR_WithPropagation ( EdbSegP s,
EdbSegP stest 
)
412{
413 if (s->Z()==stest->Z()) return TMath::Sqrt((s->X()-stest->X())*(s->X()-stest->X())+(s->Y()-stest->Y())*(s->Y()-stest->Y()));
414 Double_t zorig;
415 Double_t dR;
416 zorig=s->Z();
417 s->PropagateTo(stest->Z());
418 dR=TMath::Sqrt( (s->X()-stest->X())*(s->X()-stest->X())+(s->Y()-stest->Y())*(s->Y()-stest->Y()) );
419 s->PropagateTo(zorig);
420 return dR;
421}
Float_t Z() const
Definition: EdbSegP.h:153

◆ DeltaTheta()

Double_t EdbShowAlg::DeltaTheta ( EdbSegP s1,
EdbSegP s2 
)
426{
427// Be aware that this DeltaTheta function returns the abs() difference between the
428// ABSOLUTE values of dTheta!!! (not componentwise!
429 Double_t tx1,tx2,ty1,ty2;
430 tx1=s1->TX();
431 tx2=s2->TX();
432 ty1=s1->TY();
433 ty2=s2->TY();
434 Double_t dt= TMath::Abs(TMath::Sqrt(tx1*tx1+ty1*ty1) - TMath::Sqrt(tx2*tx2+ty2*ty2));
435 return dt;
436}
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 EdbShowAlg::DeltaThetaComponentwise ( EdbSegP s1,
EdbSegP s2 
)
442{
443// Be aware that this DeltaTheta function returns the difference between the
444// component values of dTheta!!!
445// Acutally this function should be the normal way to calculate dTheta correctly...
446 Double_t tx1,tx2,ty1,ty2;
447 tx1=s1->TX();
448 tx2=s2->TX();
449 ty1=s1->TY();
450 ty2=s2->TY();
451 Double_t dt= TMath::Sqrt( (tx1-tx2)*(tx1-tx2) + (ty1-ty2)*(ty1-ty2) );
452 return dt;
453}

◆ DeltaThetaSingleAngles()

Double_t EdbShowAlg::DeltaThetaSingleAngles ( EdbSegP s1,
EdbSegP s2 
)
456{
457// SAME function as DeltaThetaComponentwise !!!
459}
Double_t DeltaThetaComponentwise(EdbSegP *s1, EdbSegP *s2)
Definition: EdbShowAlg.cxx:441

◆ Execute()

void EdbShowAlg::Execute ( )
virtual

Reimplemented in EdbShowAlg_SA, EdbShowAlg_CA, EdbShowAlg_OI, EdbShowAlg_RC, EdbShowAlg_BW, EdbShowAlg_GS, EdbShowAlg_NN, and EdbShowAlg_N3.

324{
325 cout << "EdbShowAlg::Execute()-----------------------------------------------------------------------" << endl;
326 return;
327}

◆ Finalize()

void EdbShowAlg::Finalize ( )
virtual

◆ GetAlgName()

TString EdbShowAlg::GetAlgName ( ) const
inline
133 {
134 return eAlgName;
135 }

◆ GetAlgValue()

Int_t EdbShowAlg::GetAlgValue ( ) const
inline
130 {
131 return eAlgValue;
132 }

◆ GetMinimumDist()

Double_t EdbShowAlg::GetMinimumDist ( EdbSegP seg1,
EdbSegP seg2 
)
1944{
1945
1946// calculate minimum distance of 2 lines.
1947// use the data of (the selected object)->X(), Y(), Z(), TX(), TY().
1948// means, if the selected object == segment, use the data of the segment.
1949// or it == track, the // use the fitted data.
1950// original code from Tomoko Ariga
1951// double calc_dmin(EdbSegP *seg1, EdbSegP *seg2, double *dminz = NULL){
1952
1953 double x1,y1,z1,ax1,ay1;
1954 double x2,y2,z2,ax2,ay2;
1955 double s1,s2,s1bunsi,s1bunbo,s2bunsi,s2bunbo;
1956 double p1x,p1y,p1z,p2x,p2y,p2z,p1p2;
1957
1958 if (seg1->ID()==seg2->ID()&&seg1->PID()==seg2->PID()) return 0.0;
1959
1960 x1 = seg1->X();
1961 y1 = seg1->Y();
1962 z1 = seg1->Z();
1963 ax1= seg1->TX();
1964 ay1= seg1->TY();
1965
1966 x2 = seg2->X();
1967 y2 = seg2->Y();
1968 z2 = seg2->Z();
1969 ax2= seg2->TX();
1970 ay2= seg2->TY();
1971
1972 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);
1973 s1bunbo=(ax1*ax1+ay1*ay1+1)*(ax2*ax2+ay2*ay2+1) - (ax1*ax2+ay1*ay2+1)*(ax1*ax2+ay1*ay2+1);
1974 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);
1975 s2bunbo=(ax1*ax1+ay1*ay1+1)*(ax2*ax2+ay2*ay2+1) - (ax1*ax2+ay1*ay2+1)*(ax1*ax2+ay1*ay2+1);
1976 s1=s1bunsi/s1bunbo;
1977 s2=s2bunsi/s2bunbo;
1978 p1x=x1+s1*ax1;
1979 p1y=y1+s1*ay1;
1980 p1z=z1+s1*1;
1981 p2x=x2+s2*ax2;
1982 p2y=y2+s2*ay2;
1983 p2z=z2+s2*1;
1984 p1p2=sqrt( (p1x-p2x)*(p1x-p2x)+(p1y-p2y)*(p1y-p2y)+(p1z-p2z)*(p1z-p2z) );
1985
1986 return p1p2;
1987}
Int_t ID() const
Definition: EdbSegP.h:147
Int_t PID() const
Definition: EdbSegP.h:148

◆ GetRecoShowerArray()

TObjArray * EdbShowAlg::GetRecoShowerArray ( ) const
inline
139 {
140 return eRecoShowerArray;
141 }

◆ GetRecoShowerArrayN()

Int_t EdbShowAlg::GetRecoShowerArrayN ( ) const
inline
136 {
137 return eRecoShowerArrayN;
138 }

◆ GetShower()

EdbShowerP * EdbShowAlg::GetShower ( Int_t  i) const
inline
142 {
143 EdbShowerP* shower = (EdbShowerP*)eRecoShowerArray->At(i);
144 return shower;
145 }
Definition: EdbShowerP.h:28

◆ GetSpatialDist()

Double_t EdbShowAlg::GetSpatialDist ( EdbSegP s1,
EdbSegP s2 
)
464{
465// Mainly Z values should dominate... since the are at the order of 10k microns and x,y of 1k microns
466 Double_t x1,x2,y1,y2,z1,z2;
467 x1=s1->X();
468 x2=s2->X();
469 y1=s1->Y();
470 y2=s2->Y();
471 z1=s1->Z();
472 z2=s2->Z();
473 Double_t dist= TMath::Sqrt( (x1-x2)*(x1-x2) + (y1-y2)*(y1-y2) + (z1-z2)*(z1-z2) );
474//cout << "dist = " << dist << endl;
475 return dist;
476}

◆ Help()

void EdbShowAlg::Help ( )
514{
515 cout << "---------------------------------------------------------" << endl;
516 cout << "EdbShowAlg::Help()" << endl << endl;
517 cout << " What is this?" << endl;
518 cout << " Primary Reconstruction Class for the Showers" << endl;
519 cout << " Commonly, it takes as basic data input an EdbPVRec object (gAli)" << endl;
520 cout << endl;
521 cout << " What can we do with this?" << endl;
522 cout << " First of all, you can choose different types of algos:" << endl;
523 cout << " ShowAlg_OI (default)" << endl;
524 cout << " ShowAlg_AS" << endl;
525 cout << " ShowAlg_RC" << endl;
526 cout << " ShowAlg_TC" << endl;
527 cout << " ShowAlg_NN" << endl;
528 cout << " ShowAlg_RC (*** to be implemented ***)" << endl;
529 cout << " ShowAlg_BW (*** to be implemented ***)" << endl;
530 cout << " ShowAlg_CL (*** to be skipped ***)" << endl;
531 cout << " ShowAlg_GS (*** to be implemented ***)" << endl;
532 cout << " ...." << endl;
533 cout << " Then you set you input via: " << endl;
534 cout << " ->SetEdbPVRec(gAli) : " << endl;
535 cout << " Then you execute the algorithm via: " << endl;
536 cout << " ->Execute() : " << endl;
537 cout << " Finally you get the TObjArray of reconstructed showers (EdbShowerP*) via: " << endl;
538 cout << " ->GetRecoShowerArray() : " << endl;
539 cout << endl;
540 cout << " Technical details, information about the class:" << endl;
541 cout << " ...." << endl;
542 cout << " ...." << endl;
543 cout << endl;
544 cout << "---------------------------------------------------------" << endl;
545 return;
546}

◆ Initialize()

void EdbShowAlg::Initialize ( )
virtual

◆ IsInConeTube()

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

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

249{
250// General Function which returns Bool if the Testing BaeTrack is in a cone defined
251// by the StartingBaseTrack. In case of starting same Z position, a distance cut of
252// 20microns is assumed....
253// In case of TestingSegment==StartingSegment this function should correctly return kTRUE also...
254 //if (gEDBDEBUGLEVEL>3) cout << "Bool_t EdbShowAlg::IsInConeTube() Test Segment " << TestingSegment << " vs. Starting Segment " << StartingSegment << endl;
255
256// We reject any TestingSegment segments which have lower Z than the StartingSegment .
257 if (StartingSegment->Z()>TestingSegment->Z() ) return kFALSE;
258
259 TVector3 x1(StartingSegment->X(),StartingSegment->Y(),StartingSegment->Z());
260 TVector3 x2(TestingSegment->X(),TestingSegment->Y(),TestingSegment->Z());
261 TVector3 direction_x1(StartingSegment->TX()*1300,StartingSegment->TY()*1300,1300);
262
263// u1 is the difference vector of the position!
264 TVector3 u1=x2-x1;
265
266 Double_t direction_x1_norm= direction_x1.Mag();
267 Double_t cosangle= (direction_x1*u1)/(u1.Mag()*direction_x1_norm);
268 Double_t angle = TMath::ACos(cosangle);
269
270// This is the old version of angle calculation. It does not give the same results as in ROOT
271// when use TVector3.Angle(&TVector3). // For this IsInConeTube() we use therefore the ROOT calculation.
272 angle=u1.Angle(direction_x1);
273
274// For the case where the two basetracks have same z position
275// the angle is about 90 degree so it makes no sense to calculate it...
276// therefore we set it artificially to zero:
277 if (StartingSegment->Z()==TestingSegment->Z() ) {
278 if (gEDBDEBUGLEVEL>3) cout << "same Z position of TestingSegment and StartingSegment, Set angle artificially to zero" << endl;
279 angle=0.0;
280// Check here for dR manually:
281//cout << DeltaR_WithoutPropagation(StartingSegment,TestingSegment) << endl;
282//StartingSegment->PrintNice();
283//TestingSegment->PrintNice();
284// cout << StartingSegment->Flag() << " " << TestingSegment->Flag() << endl;
285// cout << StartingSegment->P() << " " << TestingSegment->P() << endl;
286
287// Check for position distance for 20microns if
288// Testing Segment is in same Z as StartingSegment
289 if (gEDBDEBUGLEVEL>3) cout << "Check for position distance for 20microns if Testing Segment is in same Z as StartingSegment" << endl;
290 if (gEDBDEBUGLEVEL>3) cout << "DeltaR_WithoutPropagation(StartingSegment,TestingSegment) = "<< DeltaR_WithoutPropagation(StartingSegment,TestingSegment) << endl;
291 if (DeltaR_WithoutPropagation(StartingSegment,TestingSegment)<20) return kTRUE;
292 if (DeltaR_WithoutPropagation(StartingSegment,TestingSegment)>=20) return kFALSE;
293 }
294
296 if (gEDBDEBUGLEVEL>3) cout << "Check if AngleVector now within the ConeAngleVector (<"<< ConeAngle<<"): " << angle << endl;
297 if (angle>ConeAngle) {
298 return kFALSE;
299 }
300
302 Double_t TubeDistance = 1.0/direction_x1_norm * ( (x2-x1).Cross(direction_x1) ).Mag();
303
304 if (gEDBDEBUGLEVEL>3) cout << "Check if TestingSegment is now within the Tube (<"<< CylinderRadius<<"): " << TubeDistance << endl;
305
306 if (TubeDistance>CylinderRadius) {
307 return kFALSE;
308 }
309
310 return kTRUE;
311}
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ Print()

void EdbShowAlg::Print ( )
339{
340 Log(2,"EdbShowAlg::Print()","Print relevant settings for shower reconstruction algorithm");
341
342 Log(2,"EdbShowAlg::Print()","Algorithm name abbreviation: %s", eAlgName.Data() );
343 Log(2,"EdbShowAlg::Print()","Algorithm number (internal): %2d", eAlgValue );
344
345 Log(2,"EdbShowAlg::Print()","Algorithm parameters N: %2d", eParaN );
346 for (int i=0; i<eParaN; i++) Log(2,"EdbShowAlg::Print()","Algorithm parameter: number, name, value: %d, %40s, %.04f", i, eParaString[i].Data(), eParaValue[i] );
347
348// cout << "Plates: " << endl;
349// cout << eFirstPlate_eAliPID << " " << eLastPlate_eAliPID << " " << eMiddlePlate_eAliPID << " " << eNumberPlate_eAliPID << " " << endl;
350// cout << "UseAliSub= " << eUseAliSub << " ." << endl;
351// cout << "EdbShowAlg::Print()...done." << endl;
352 Log(2,"EdbShowAlg::Print()","Print relevant settings for shower reconstruction algorithm...done.");
353 return;
354}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
Int_t eParaN
Definition: EdbShowAlg.h:51
@ Data
Definition: tlg2couples.C:53

◆ PrintAll()

void EdbShowAlg::PrintAll ( )
388{
389 return;
390}

◆ PrintMore()

void EdbShowAlg::PrintMore ( )
377{
378 cout << "EdbShowAlg::PrintMore()" << endl;
379 cout << "eInBTArray->GetEntries();= " << eInBTArray->GetEntries() << " ." << endl;
380 cout << "EdbShowAlg::PrintMore()...done." << endl;
381 return;
382}
TObjArray * eInBTArray
Definition: EdbShowAlg.h:62

◆ PrintParameters()

void EdbShowAlg::PrintParameters ( )
359{
360 Print();
361 return;
362}
void Print()
Definition: EdbShowAlg.cxx:338

◆ PrintParametersShort()

void EdbShowAlg::PrintParametersShort ( )
367{
368 cout << eAlgName<< " :";
369 for (int i=0; i<5; i++) cout << setw(6) << eParaValue[i];
370 cout << "."<<endl;
371 return;
372}

◆ Set0()

void EdbShowAlg::Set0 ( )
protected
54{
55// Reset all internal variables.
56 eParaN =10; // depends on the algorithm, how much parameters it has got. to be set in the constructor of the specific class...
57 eAlgName="UNSPECIFIED";
58 eAlgValue=-999999;
59 for (int i=0; i<10; ++i) {
60 eParaValue[i]=-99999.0;
61 eParaString[i]="UNSPECIFIED";
62 }
63 eAli_Sub=0;
64
65// Do not use use small eAli Object by default:
66// (this solution is memory safe....)
67 eUseAliSub=0;
68
70
75}
Int_t eInBTArrayN
Definition: EdbShowAlg.h:63
Int_t eUseAliSub
Definition: EdbShowAlg.h:78
Int_t eActualAlgParametersetNr
Definition: EdbShowAlg.h:54
#define NULL
Definition: nidaqmx.h:84

◆ SetActualAlgParameterset()

void EdbShowAlg::SetActualAlgParameterset ( Int_t  ActualAlgParametersetNr)
inline
126 {
127 eActualAlgParametersetNr=ActualAlgParametersetNr;
128 }

◆ SetEdbPVRec()

void EdbShowAlg::SetEdbPVRec ( EdbPVRec Ali)
inline
100 {
101 eAli = Ali;
103 }
Int_t Npatterns() const
Definition: EdbPattern.h:366
Int_t eAliNpat
Definition: EdbShowAlg.h:60
EdbPVRec * eAli
Definition: EdbShowAlg.h:59

◆ SetEdbPVRecPIDNumbers()

void EdbShowAlg::SetEdbPVRecPIDNumbers ( Int_t  FirstPlate_eAliPID,
Int_t  LastPlate_eAliPID,
Int_t  MiddlePlate_eAliPID,
Int_t  NumberPlate_eAliPID 
)
inline
108 {
109 eFirstPlate_eAliPID=FirstPlate_eAliPID;
110 eLastPlate_eAliPID=LastPlate_eAliPID;
111 eMiddlePlate_eAliPID=MiddlePlate_eAliPID;
112 eNumberPlate_eAliPID=NumberPlate_eAliPID;
113 }
Int_t eFirstPlate_eAliPID
Definition: EdbShowAlg.h:65
Int_t eLastPlate_eAliPID
Definition: EdbShowAlg.h:66
Int_t eMiddlePlate_eAliPID
Definition: EdbShowAlg.h:67
Int_t eNumberPlate_eAliPID
Definition: EdbShowAlg.h:68

◆ SetInBTArray()

void EdbShowAlg::SetInBTArray ( TObjArray *  InBTArray)
inline
104 {
105 eInBTArray = InBTArray;
106 eInBTArrayN=eInBTArray->GetEntries();
107 }

◆ SetParameter()

void EdbShowAlg::SetParameter ( Int_t  parNr,
Float_t  par 
)
91{
92// SetParameter from the given value.
93 if (parNr>9) {
94 cout << "EdbShowAlg::SetParameters() parNr>9 not supported. Dont set par value!"<<endl;
95 }
96 eParaValue[parNr]=par;
97 cout << "EdbShowAlg::SetParameter()...done"<<endl;
98}

◆ SetParameters()

void EdbShowAlg::SetParameters ( Float_t *  par)
80{
81// SetParameters from the given array.
82 for (int i=0; i<10; ++i) {
83 eParaValue[i]=par[i];
84 }
85 cout << "EdbShowAlg::SetParameters()...done"<<endl;
86}

◆ SetRecoShowerArray()

void EdbShowAlg::SetRecoShowerArray ( TObjArray *  RecoShowerArray)
inline
114 {
116 }
TObjArray * RecoShowerArray
Definition: Shower_E_FromShowerRoot.C:12

◆ SetRecoShowerArrayN()

void EdbShowAlg::SetRecoShowerArrayN ( Int_t  RecoShowerArrayN)
inline
117 {
118 eRecoShowerArrayN=RecoShowerArrayN;
119 }

◆ SetUseAliSub()

void EdbShowAlg::SetUseAliSub ( Bool_t  UseAliSub)
inline
147 {
148 eUseAliSub=UseAliSub;
149 }

◆ Transform_eAli()

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

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.....

108{
109
110// Transform eAli to eAli_Sub:
111// the lenght of eAli_Sub is not changed.
112// Only XY-size (and MC) cuts are applied.
113// --------------------------------------------------------------------------------------
114// ---
115// --- Whereas in ShowRec.cpp the treebranch file is written directly after each
116// --- BT reconstruction, it was not a problem when the eAliSub was deleted each time.
117// --- But now situation is different, since the BTs of the showers have the adresses from
118// --- the eAli_Sub !each! so if one eAli_Sub is deleted the EdbShowerP object has lost its
119// --- BT adresses.
120// --- So if we do not delete the eAli_Sub, reconstrtuction is fast (compared to eAli) but
121// --- memory increases VERY fast.
122// --- If we do use eAli, then memory consumption will not increase fast, but
123// --- reconstruction is slow.
124// ---
125// --- Possible workarounds:
126// --- For few InBTs (EDA use for data...): use eAliSub and dont delete it.
127// --- For many InBTs (parameterstudies): use eAli---.
128// --- Use always eAliSub, search BT correspond in eAli, add BT from eAli....May take also long time...
129// ---
130// --- SEVERE wARNING: IF gAli is in wrong order it can be that no showers.
131// --- SEVERE wARNING: are reconstructed since most implemented algorithms
132// --- SEVERE wARNING: rely on the order that plate 2 comes directly behind
133// --- the InBT. on the todo list....
134// ---
135// --------------------------------------------------------------------------------------
136
137
138 Log(3,"EdbShowAlg::Transform_eAli(EdbSegP* InitiatorBT, Float_t ExtractSize=1500)","EdbShowAlg::Transform_eAli");
139
140// IF TO RECREATE THE gALI_SUB in RECONSTRUCTION or if to use gAli global (slowlier but maybe memory leak safe).
141 if (!eUseAliSub)
142 {
143 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg::Transform_eAli UseAliSub==kFALSE No new eAli_Sub created. Use eAli instead. "<<endl;
145 eAli_SubNpat=eAli_Sub->Npatterns(); //number of plates
146 if (gEDBDEBUGLEVEL>3) eAli_Sub->Print();
147 return;
148 }
149
150 Int_t npat;
151 npat = eAli->Npatterns(); //number of plates
152
153// has to be deleted in some part of the script outside this function...
154// Dont forget , otherwise memory heap overflow!
156 if (eAli_Sub) {
157 ;
158 }
159// eAli_Sub = new EdbPVRec();
160
161
162 if (eUseAliSub) {
163 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg::Transform_eAli UseAliSub==kTRUE Will now create new eAli_Sub. "<<endl;
164
166 eAli_Sub = new EdbPVRec();
167// cout << "Adress of eAli_Sub " << eAli_Sub << endl;
168// cout << "eAli_Sub->Npatterns(); " << eAli_Sub->Npatterns() << endl;
169// eAli_Sub->Print();
170 }
171
172
173
174// cout << "PROBLE THRANFOMING, CAUSE ADRESES ARE DELETED ALSO...." << endl;
175// cout << "TEMPORAY SOLUTION: comment the delete eAli_Sub ( will lead to large memory consumption when run for long time" << endl;
176
177
178// Create SubPattern objects
179 EdbSegP* ExtrapolateInitiatorBT=0;
180 ExtrapolateInitiatorBT = (EdbSegP*)InitiatorBT->Clone();
181
182 Int_t InitiatorBTMCEvt=InitiatorBT->MCEvt();
183
184// Create Variables For ExtractSubpattern boundaries
185 Float_t mini[5];
186 Float_t maxi[5];
187//Float_t ExtractSize=1000;// now in fucntion header
188 mini[0]=ExtrapolateInitiatorBT->X()-ExtractSize;
189 mini[1]=ExtrapolateInitiatorBT->Y()-ExtractSize;
190 maxi[0]=ExtrapolateInitiatorBT->X()+ExtractSize;
191 maxi[1]=ExtrapolateInitiatorBT->Y()+ExtractSize;
192 mini[2]=-0.7;
193 mini[3]=-0.7;
194 mini[4]=0.0;
195 maxi[2]=0.7;
196 maxi[3]=0.7;
197 maxi[4]=100.0;
198
199 EdbPattern* singlePattern;
200 Float_t ExtrapolateInitiatorBT_zpos_orig=ExtrapolateInitiatorBT->Z();
201
202// Add the subpatterns in a loop for the plates:
203// in reverse ordering.due to donwstream behaviour (!):
204// (Only downstream is supported now...)
205 for (Int_t ii=0; ii<npat; ++ii) {
206
207 Float_t zpos=eAli->GetPattern(ii)->Z();
208
209 ExtrapolateInitiatorBT->PropagateTo(zpos);
210
211 mini[0]=ExtrapolateInitiatorBT->X()-ExtractSize;
212 mini[1]=ExtrapolateInitiatorBT->Y()-ExtractSize;
213 maxi[0]=ExtrapolateInitiatorBT->X()+ExtractSize;
214 maxi[1]=ExtrapolateInitiatorBT->Y()+ExtractSize;
215
216 singlePattern=(EdbPattern*)eAli->GetPattern(ii)->ExtractSubPattern(mini,maxi,InitiatorBTMCEvt);
217// This sets PID() analogue to (upstream), nut not PID of the BTs !
218 singlePattern-> SetID(eAli->GetPattern(ii)->ID());
219// This sets PID() analogue to (upstream), nut not PID of the BTs !
220 singlePattern-> SetPID(eAli->GetPattern(ii)->PID());
221
222 eAli_Sub->AddPattern(singlePattern);
223
224// Propagate back...! (in order not to change the original BT)
225 ExtrapolateInitiatorBT->PropagateTo(ExtrapolateInitiatorBT_zpos_orig);
226 }
227
228 delete ExtrapolateInitiatorBT;
229
230 eAli_SubNpat=eAli_Sub->Npatterns(); //number of plates
231
232
233 if (gEDBDEBUGLEVEL>2) {
234 cout << "EdbShowAlg::Transform_eAli eAli_Sub created at : "<< eAli_Sub << endl;
235 cout << "EdbShowAlg::Transform_eAli eAli_Sub->Npatterns(): "<< eAli_Sub->Npatterns() << endl;
236 }
237
238 Log(3,"EdbShowAlg::Transform_eAli(EdbSegP* InitiatorBT, Float_t ExtractSize=1500)...done.","EdbShowAlg::Transform_eAli...done.");
239
240 return;
241}
Int_t npat
Definition: Xi2HatStartScript.C:33
Definition: EdbPVRec.h:148
Definition: EdbPattern.h:273
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
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
Definition: EdbSegP.h:21
void PropagateTo(float z)
Definition: EdbSegP.cxx:293
Int_t MCEvt() const
Definition: EdbSegP.h:145
Float_t Z() const
Definition: EdbPattern.h:84
Int_t eAli_SubNpat
Definition: EdbShowAlg.h:75

◆ UpdateShowerIDs()

void EdbShowAlg::UpdateShowerIDs ( )
480{
481 // Shower IDs are numbered in the ordering in the RecoShowerArray
482 if ( NULL == eRecoShowerArray ) return;
483 Int_t idNumber=0;
484 for (int i=0; i<eRecoShowerArrayN; ++i) {
485 EdbShowerP* shower = GetShower(i);
486 shower->SetID(i);
487 }
488 return;
489}
void SetID(int id)
Definition: EdbSegP.h:128
EdbShowerP * GetShower(Int_t i) const
Definition: EdbShowAlg.h:142

◆ UpdateShowerMetaData()

void EdbShowAlg::UpdateShowerMetaData ( )
493{
494 // To know which Algorithm and Parameters were used
495 // for this shower to be reconstructed.
496 if ( NULL == eRecoShowerArray ) return;
497 Int_t idNumber=0;
498 for (int i=0; i<eRecoShowerArrayN; ++i) {
499 EdbShowerP* shower = GetShower(i);
500 shower->SetAlgName(eAlgName);
501 shower->SetAlgValue(eAlgValue);
502 for (Int_t j=0; j<10; ++j) {
503 shower->SetAlgParaValue(j, eParaValue[j]);
504 shower->SetAlgParaString(j, eParaString[j]);
505 }
506 }
507 return;
508}
void SetAlgParaString(Int_t ParaStringNr, TString ParaString)
Definition: EdbShowerP.h:623
void SetAlgValue(Int_t AlgValue)
Definition: EdbShowerP.h:615
void SetAlgName(TString AlgName)
Definition: EdbShowerP.h:612
void SetAlgParaValue(Int_t ParaValueNr, Float_t ParaValue)
Definition: EdbShowerP.h:619

Member Data Documentation

◆ eActualAlgParametersetNr

Int_t EdbShowAlg::eActualAlgParametersetNr
protected

◆ eAlgName

TString EdbShowAlg::eAlgName
protected

◆ eAlgValue

Int_t EdbShowAlg::eAlgValue
protected

◆ eAli

EdbPVRec* EdbShowAlg::eAli
protected

◆ eAli_Sub

EdbPVRec* EdbShowAlg::eAli_Sub
protected

◆ eAli_SubNpat

Int_t EdbShowAlg::eAli_SubNpat
protected

◆ eAliNpat

Int_t EdbShowAlg::eAliNpat
protected

◆ eFirstPlate_eAliPID

Int_t EdbShowAlg::eFirstPlate_eAliPID
protected

◆ eInBTArray

TObjArray* EdbShowAlg::eInBTArray
protected

◆ eInBTArrayN

Int_t EdbShowAlg::eInBTArrayN
protected

◆ eLastPlate_eAliPID

Int_t EdbShowAlg::eLastPlate_eAliPID
protected

◆ eMiddlePlate_eAliPID

Int_t EdbShowAlg::eMiddlePlate_eAliPID
protected

◆ eNumberPlate_eAliPID

Int_t EdbShowAlg::eNumberPlate_eAliPID
protected

◆ eParaN

Int_t EdbShowAlg::eParaN
protected

◆ eParaString

TString EdbShowAlg::eParaString[10]
protected

◆ eParaValue

Float_t EdbShowAlg::eParaValue[10]
protected

◆ eRecoShower

EdbShowerP* EdbShowAlg::eRecoShower
protected

◆ eRecoShowerArray

TObjArray* EdbShowAlg::eRecoShowerArray
protected

◆ eRecoShowerArrayN

Int_t EdbShowAlg::eRecoShowerArrayN
protected

◆ eUseAliSub

Int_t EdbShowAlg::eUseAliSub
protected

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