FEDRA emulsion software from the OPERA Collaboration
EdbShowAlg_CA Class Reference

#include <EdbShowAlg.h>

Inheritance diagram for EdbShowAlg_CA:
Collaboration diagram for EdbShowAlg_CA:

Public Member Functions

 ClassDef (EdbShowAlg_CA, 1)
 
 EdbShowAlg_CA ()
 
void Execute ()
 
void Finalize ()
 
void Initialize ()
 
virtual ~EdbShowAlg_CA ()
 
- Public Member Functions inherited from EdbShowAlg
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 ()
 

Private Member Functions

Bool_t FindPrecedingBTs (EdbSegP *s, EdbSegP *InBT, EdbPVRec *gAli, EdbShowerP *shower)
 
Bool_t GetConeTubeDistanceToInBT (EdbSegP *sa, EdbSegP *s, Double_t CylinderRadius, Double_t ConeAngle)
 

Additional Inherited Members

- Protected Member Functions inherited from EdbShowAlg
void Set0 ()
 
- Protected Attributes inherited from EdbShowAlg
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_CA()

EdbShowAlg_CA::EdbShowAlg_CA ( )
794{
795// Default Constructor
796 cout << "EdbShowAlg_CA::EdbShowAlg_CA() Default Constructor"<<endl;
797
798// Reset all:
799 Set0();
800
801 eAlgName="CA";
802 eAlgValue=2; // see default.par_SHOWREC for labelling
803 eParaN =4; // depends on the algorithm, how much parameters it has got.
804
805 eParaValue[0]=700;
806 eParaString[0]="CylinderRadius"; // micrometer
807 eParaValue[1]=0.025;
808 eParaString[1]="ConeAngle"; // radiant
809 eParaValue[2]=110;
810 eParaString[2]="DeltaRToPreceedingBT"; // micrometer
811 eParaValue[3]=0.05;
812 eParaString[3]="DeltaThetaToPreceedingBT"; // radiant
813
814// 700 0.025 110 0.05
815}
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
Int_t eParaN
Definition: EdbShowAlg.h:51

◆ ~EdbShowAlg_CA()

EdbShowAlg_CA::~EdbShowAlg_CA ( )
virtual
820{
821 cout << "EdbShowAlg_CA::~EdbShowAlg_CA()"<<endl;
822// Default Destructor
823}

Member Function Documentation

◆ ClassDef()

EdbShowAlg_CA::ClassDef ( EdbShowAlg_CA  ,
 
)

◆ Execute()

void EdbShowAlg_CA::Execute ( )
virtual

Reimplemented from EdbShowAlg.

835{
836 cout << "EdbShowAlg_CA::Execute()" << endl;
837 cout << "EdbShowAlg_CA::Execute DOING MAIN SHOWER RECONSTRUCTION HERE" << endl;
838
839 EdbSegP* InBT;
840 EdbSegP* Segment;
841 EdbShowerP* RecoShower=0;
842
843 Bool_t StillToLoop=kTRUE;
844 Int_t ActualPID;
845 Int_t newActualPID;
846 Int_t STEP=-1;
847 Int_t NLoopedPattern=0;
848
849 if (gEDBDEBUGLEVEL>3) {
850 cout << "======================================================"<<endl;
851 cout << " eFirstPlate_eAliPID = " << eFirstPlate_eAliPID << endl;
852 cout << " eLastPlate_eAliPID = " << eLastPlate_eAliPID << endl;
853 cout << " eNumberPlate_eAliPID = " << eNumberPlate_eAliPID << endl;
854 cout << " eAli->GetPattern(eFirstPlate_eAliPID)->Z() = " << eAli->GetPattern(eFirstPlate_eAliPID)->Z() << endl;
855 cout << " eAli->GetPattern(eLastPlate_eAliPID)->Z() = " << eAli->GetPattern(eLastPlate_eAliPID)->Z() << endl;
856 cout << " eInBTArrayN = " << eInBTArrayN << endl;
857 cout << "======================================================"<<endl;
858 }
859
861 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- STEP for patternloop direction = " << STEP << endl;
862
863
864//--- Loop over InBTs:
865
866// Since eInBTArray is filled in ascending ordering by zpositon
867// We use the descending loop to begin with BT with lowest z first.
868 for (Int_t i=eInBTArrayN-1; i>=0; --i) {
869
870// CounterOutPut
871 if (gEDBDEBUGLEVEL==2) if ((i%50)==0) cout << eInBTArrayN <<" InBT in total, still to do:"<<Form("%4d",i)<< "\r\r\r\r"<<flush;
872 if (gEDBDEBUGLEVEL>=4) cout << "Doing i=" << i << endl;
873//-----------------------------------
874// 1) Make local_gAli with cut parameters:
875//-----------------------------------
876
877// Create new EdbShowerP Object for storage;
878// See EdbShowerP why I have to call the Constructor as "unique" ideficable value
879 // cout << "Create new EdbShowerP Object for storage."<<endl;
880 // cout << i << endl;
881 // cout << eAlgValue << endl;
882 // cout << eActualAlgParametersetNr << endl;
884 // cout << "Create new EdbShowerP Object for storage;...done."<<endl;
885
886// Get InitiatorBT from eInBTArray
887 InBT=(EdbSegP*)eInBTArray->At(i);
888 if (gEDBDEBUGLEVEL>2) InBT->PrintNice();
889 Float_t X0=InBT->X();
890 Float_t Y0=InBT->Y();
891
892// Add InBT to RecoShower:
893// This has to be done, since by definition the first BT in the RecoShower is the InBT.
894// Otherwise, also the definition of shower axis and transversal profiles is wrong!
895 RecoShower -> AddSegment(InBT);
896 if (gEDBDEBUGLEVEL>4) cout << "Segment (InBT) " << Segment << " was added to RecoShower." << endl;
897
898
899// Transform (make size smaller, extract only events having same MC) the eAli object:
900 Transform_eAli(InBT);
901 if (gEDBDEBUGLEVEL>3) eAli_Sub->Print();
902
903 if (gEDBDEBUGLEVEL>2) cout << " eAli_Sub->GetPattern(eFirstPlate_eAliPID)->Z() = " << eAli_Sub->GetPattern(eFirstPlate_eAliPID)->Z() << endl;
904 if (gEDBDEBUGLEVEL>2) cout << " eAli_Sub->GetPattern(eLastPlate_eAliPID)->Z() = " << eAli_Sub->GetPattern(eLastPlate_eAliPID)->Z() << endl;
905
906
907//-----------------------------------
908// 2) Loop over (whole) eAli, check BT for Cuts
909// eAli_Sub
910//-----------------------------------
911 ActualPID= InBT->PID() ;
912 newActualPID= InBT->PID() ;
913
914 while (StillToLoop) {
915 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- --- Doing patterloop id= " << ActualPID << " for pattern with Z position= " << eAli_Sub->GetPattern(ActualPID)->Z() << endl;
916
917 for (Int_t btloop_cnt=0; btloop_cnt<eAli_Sub->GetPattern(ActualPID)->GetN(); ++btloop_cnt) {
918
919 Segment = (EdbSegP*)eAli_Sub->GetPattern(ActualPID)->GetSegment(btloop_cnt);
920 if (gEDBDEBUGLEVEL>4) Segment->PrintNice();
921
922
923// Now apply cut conditions: CA ConeAdvanced Alg --------------------
924 if ( Segment->MCEvt() > 0 ) if ( Segment->MCEvt()!=InBT->MCEvt() ) continue; // MCEvtNr (>0) or BgMCNr (-999)
925 if ( Abs(Segment->X()-X0) > 7000 ) continue;
926 if ( Abs(Segment->Y()-Y0) > 7000 ) continue;
927 if ( !GetConeTubeDistanceToInBT(Segment, InBT, eParaValue[0], eParaValue[1]) ) continue;
928// if ( !FindPrecedingBTs(Segment, InBT, eAli, RecoShower)) continue;
929 if ( !FindPrecedingBTs(Segment, InBT, eAli_Sub, RecoShower)) continue;
930// end of cut conditions: CA ConeAdvanced Alg --------------------
931
932// If we arrive here, Basetrack Segment has passed criteria
933// and is then added to the RecoShower:
934// Check if its not the InBT which is already added:
935 if (Segment->X()==InBT->X()&&Segment->Y()==InBT->Y()) {
936 ; // is InBT, do nothing;
937 }
938 else {
939 RecoShower -> AddSegment(Segment);
940 }
941 if (gEDBDEBUGLEVEL>4) cout << "Segment " << Segment << " was added at &Segment : " << &Segment << endl;
942
943
944 } // of btloop_cnt
945
946 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- --- ActualPID= " << newActualPID << " done. Reconstructed shower has up to now: " << RecoShower->N() << " Segments." << endl;
947
948//------------
949 newActualPID=ActualPID+STEP;
950 ++NLoopedPattern;
951
952 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- --- newActualPID= " << newActualPID << endl;
953 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- --- NLoopedPattern= " << NLoopedPattern << endl;
954 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- --- eNumberPlate_eAliPID= " << eNumberPlate_eAliPID << endl;
955 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- --- StillToLoop= " << StillToLoop << endl;
956
957// This if holds in the case of STEP== +1
958 if (STEP==1) {
959 if (newActualPID>eLastPlate_eAliPID) StillToLoop=kFALSE;
960 if (newActualPID>eLastPlate_eAliPID) cout << "EdbShowAlg_CA::Execute--- ---Stop Loop since: newActualPID>eLastPlate_eAliPID"<<endl;
961 }
962// This if holds in the case of STEP== -1
963 if (STEP==-1) {
964 if (newActualPID<eLastPlate_eAliPID) StillToLoop=kFALSE;
965 if (newActualPID<eLastPlate_eAliPID && gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- ---Stop Loop since: newActualPID<eLastPlate_eAliPID"<<endl;
966 }
967// This if holds general, since eNumberPlate_eAliPID is not dependent of the structure of the gAli subject:
968 if (NLoopedPattern>eNumberPlate_eAliPID) StillToLoop=kFALSE;
969 if (NLoopedPattern>eNumberPlate_eAliPID && gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- ---Stop Loop since: NLoopedPattern>eNumberPlate_eAliPID"<<endl;
970
971 ActualPID=newActualPID;
972 } // of // while (StillToLoop)
973
974// Obligatory when Shower Reconstruction is finished!
975 RecoShower->Update();
976
977
978 if (gEDBDEBUGLEVEL>3) RecoShower ->PrintBasics();
979 if (gEDBDEBUGLEVEL>3) RecoShower ->PrintNice();
980 if (gEDBDEBUGLEVEL>3) RecoShower ->PrintSegments();
981
982// if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- Before adding to array delete the histograms...by finalize() of shower."<<endl;
983// RecoShower ->Finalize();
984// if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_CA::Execute--- Finalize done..."<<endl;
985
986
987// Add Shower Object to Shower Reco Array.
988// Not, if its empty or containing only one BT:
989 if (RecoShower->N()>1) {
990 AddRecoShowerArray(RecoShower);
991 }
992 else {
993 // Delete the created shower object. Saves _lots_ of memory and speeds up reconstruction for large #InBTs!
994 // TODO Also Include this in the other RecoAlgs ...
995 Log(3,"EdbShowAlg_CA::Execute()","InBT # %d RecoShower->N()<=1: Dont add shower to RecoShowerArray. Delete this shower for memory safing.");
996 delete RecoShower;
997 RecoShower=NULL;
998 }
999
1000//cout << " eRecoShowerArray->GetEntries() " << eRecoShowerArray->GetEntries() << endl;
1001
1002
1003// Set back loop values:
1004 StillToLoop=kTRUE;
1005 NLoopedPattern=0;
1006 } // of // for (Int_t i=eInBTArrayN-1; i>=0; --i) {
1007
1008
1009// Set new value for eRecoShowerArrayN (may now be < eInBTArrayN).
1011
1012
1013 cout << "EdbShowAlg_CA::eRecoShowerArray() Entries: " << eRecoShowerArray->GetEntries() << endl;
1014 cout << "EdbShowAlg_CA::Execute()...done." << endl;
1015 return;
1016}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
brick X0
Definition: RecDispMC.C:112
void Print() const
Definition: EdbPattern.cxx:1693
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
Definition: EdbSegP.h:21
Float_t X() const
Definition: EdbSegP.h:173
Float_t Y() const
Definition: EdbSegP.h:174
void PrintNice() const
Definition: EdbSegP.cxx:392
Int_t PID() const
Definition: EdbSegP.h:148
Int_t MCEvt() const
Definition: EdbSegP.h:145
Float_t Z() const
Definition: EdbPattern.h:84
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
Bool_t FindPrecedingBTs(EdbSegP *s, EdbSegP *InBT, EdbPVRec *gAli, EdbShowerP *shower)
Definition: EdbShowAlg.cxx:1088
Bool_t GetConeTubeDistanceToInBT(EdbSegP *sa, EdbSegP *s, Double_t CylinderRadius, Double_t ConeAngle)
Definition: EdbShowAlg.cxx:1027
void Transform_eAli(EdbSegP *InitiatorBT, Float_t ExtractSize)
Definition: EdbShowAlg.cxx:107
void AddRecoShowerArray(EdbShowerP *shower)
Definition: EdbShowAlg.h:120
void SetRecoShowerArrayN(Int_t RecoShowerArrayN)
Definition: EdbShowAlg.h:117
TObjArray * eRecoShowerArray
Definition: EdbShowAlg.h:70
Int_t eFirstPlate_eAliPID
Definition: EdbShowAlg.h:65
Int_t eLastPlate_eAliPID
Definition: EdbShowAlg.h:66
Int_t eInBTArrayN
Definition: EdbShowAlg.h:63
TObjArray * eInBTArray
Definition: EdbShowAlg.h:62
EdbPVRec * eAli
Definition: EdbShowAlg.h:59
Int_t eNumberPlate_eAliPID
Definition: EdbShowAlg.h:68
EdbPVRec * eAli_Sub
Definition: EdbShowAlg.h:74
Int_t eActualAlgParametersetNr
Definition: EdbShowAlg.h:54
Definition: EdbShowerP.h:28
Int_t N() const
Definition: EdbShowerP.h:412
void PrintNice()
Definition: EdbShowerP.cxx:2339
void PrintSegments()
Definition: EdbShowerP.cxx:2367
void PrintBasics()
Definition: EdbShowerP.cxx:2360
void Update()
Definition: EdbShowerP.cxx:975
gEDBDEBUGLEVEL
Definition: energy.C:7
#define NULL
Definition: nidaqmx.h:84

◆ Finalize()

void EdbShowAlg_CA::Finalize ( )
virtual

Reimplemented from EdbShowAlg.

1022{
1023 return;
1024}

◆ FindPrecedingBTs()

Bool_t EdbShowAlg_CA::FindPrecedingBTs ( EdbSegP s,
EdbSegP InBT,
EdbPVRec gAli,
EdbShowerP shower 
)
private
1089{
1090 EdbSegP* s_TestBT;
1091// Int_t nentries=shower->GetEntries();
1092 Int_t nentries=shower->N();
1093 Double_t dZ;
1094
1095// In case the shower is empty we do not have to search for a preceeding BT:
1096 if (nentries==0) return kTRUE;
1097
1098// We do not test BaseTracks which are before the initiator BaseTrack!
1099 if (s->Z()<InBT->Z()) {
1100 cout << "--- --- EdbShowAlg_CA::FindPrecedingBTs(): s->Z()<InBT->Z()..: We do not test BaseTracks which are before the initiator BaseTrack! return kFALSE;" << endl;
1101 return kFALSE;
1102 }
1103
1104// For the very first Z position we do not test
1105// if testBT has Preceeders, only if it it has a BT around (case for e+e- coming from gammma):
1106// Take 50microns and 80mrad in (dR/dT) around.
1107// This does not affect the normal results, but helps for
1108// events which may have a second BT close to InBT (like in e+e-)
1109 if (s->Z()==InBT->Z()) {
1110//cout << DeltaTheta(s, InBT) << " " << DeltaR_WithPropagation(s, InBT) << endl;
1111 if (DeltaTheta(s, InBT) < 0.08 && DeltaR_WithoutPropagation(s, InBT) < 50.0 ) return kTRUE;
1112 }
1113
1114 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_CA::FindPrecedingBTs(): Testing " << nentries << " entries of the shower." << endl;
1115
1116 for (Int_t i=nentries-1; i>=0; --i) {
1117 s_TestBT = (EdbSegP*)( shower->GetSegment(i) );
1118
1119 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_CA::FindPrecedingBTs(): Do s_TestBT->ID() s->ID() s_TestBT->MCEvt() s_TestBT->Z() s->Z() "<< s_TestBT->ID() << " " << s->ID() << " " << s_TestBT->MCEvt() <<" " << s_TestBT->Z() << " " << s->Z() <<endl;
1120 if (gEDBDEBUGLEVEL>3) s_TestBT->PrintNice();
1121 if (gEDBDEBUGLEVEL>3) s->PrintNice();
1122
1123
1124 dZ=TMath::Abs(s_TestBT->Z()-s->Z());
1125 if (dZ<30) continue; // Exclude the case of same Zpositions...
1127 if (dZ>(3.5*(1300.0+30.0))) continue; // Exclude the case of more than 3 plates before... // modified
1128 // often, experimentally observed, dZ approx 1350 ... so the upper cut really does not
1129 // look for 3 plates, but only for two ...: seems to work now. (24.Oct.2016)
1130
1131 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_CA::FindPrecedingBTs(): Checking dT,dR and dZ for i: " << i << " " << DeltaTheta(s, s_TestBT) << " " << DeltaR_WithPropagation(s, s_TestBT) << " "<<dZ << endl;
1132
1133 if (DeltaTheta(s, s_TestBT) > eParaValue[3] ) continue;
1134 if (DeltaR_WithPropagation(s, s_TestBT) > eParaValue[2]) continue;
1135
1136 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_CA::FindPrecedingBTs(): Checking dT,dR and dZ for i: " << i << " " << DeltaTheta(s, s_TestBT) << " " << DeltaR_WithPropagation(s, s_TestBT) << " "<<dZ << " ok!"<<endl;
1137 return kTRUE;
1138 }
1139//---------------------------------------------
1140 return kFALSE;
1141}
Int_t ID() const
Definition: EdbSegP.h:147
Float_t Z() const
Definition: EdbSegP.h:153
Double_t DeltaR_WithoutPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowAlg.cxx:404
Double_t DeltaR_WithPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowAlg.cxx:411
Double_t DeltaTheta(EdbSegP *s1, EdbSegP *s2)
Definition: EdbShowAlg.cxx:425
EdbSegP * GetSegment(int i) const
Definition: EdbShowerP.h:435
int nentries
Definition: check_shower.C:40
s
Definition: check_shower.C:55

◆ GetConeTubeDistanceToInBT()

Bool_t EdbShowAlg_CA::GetConeTubeDistanceToInBT ( EdbSegP sa,
EdbSegP s,
Double_t  CylinderRadius,
Double_t  ConeAngle 
)
private

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

1028{
1029// This is an updated version since the TMath::ACos(cosangle) does not really refer to a coneangle:
1030// We rather use from now on the ROOT implementation of TVector3::Angle()
1031
1032 TVector3 x1(InBT->X(),InBT->Y(),InBT->Z());
1033 TVector3 x2(sa->X(),sa->Y(),sa->Z());
1034 TVector3 direction_x1(InBT->TX()*1300,InBT->TY()*1300,1300);
1035 TVector3 direction_x2(sa->TX()*1300,sa->TY()*1300,1300);
1036 TVector3 u1=x2-x1;
1037
1038 Double_t direction_x1_norm= direction_x1.Mag();
1039 Double_t cosangle= (direction_x1*u1)/(u1.Mag()*direction_x1_norm);
1040 Double_t angle = TMath::ACos(cosangle); // NO THIS IS NOT THE CONE ANGLE!! (this is rather an old relict)
1041
1042 TVector3 direction_1(InBT->TX()*1300,InBT->TY()*1300,1300);
1043 TVector3 direction_2(sa->TX()*1300,sa->TY()*1300,1300);
1044
1045 /*
1046 // DEBUG STUFF ... TO BE REMOVED IF EVERYTHING IS CLEAR
1047 cout <<"------"<<endl;
1048 cout << "angle = " << angle << endl;
1049 cout << "angle V2 = " << x1.Angle(x2) << endl;
1050 cout << "angle V3 = " << direction_x1.Angle(direction_x2) << endl;
1051 cout << "angle V4 = " << direction_1.Angle(direction_2) << " " << direction_1.Angle(direction_2)/3.14*180.0 << endl;
1052
1053 cout << "angle V5 direction_x1.Angle(u1 = " << direction_x1.Angle(u1) << endl;
1054 cout << "angle V5 u1.Angle(direction_x1 = " << u1.Angle(direction_x1) << endl;
1055 */
1056 angle=u1.Angle(direction_x1);
1057
1058// For the case where the two basetracks have same z position
1059// the angle is about 90 degree so it makes no sense to calculate it...
1060// therefore we set it artificially to zero:
1061 if (InBT->Z()==sa->Z() ) {
1062 angle=0.0;
1063//if (gEDBDEBUGLEVEL>3) //cout << "same z position, set angle artificially to zero" << endl;
1064 }
1065
1067 if (angle>ConeAngle) {
1068 return kFALSE;
1069 }
1070
1072 Double_t TubeDistance = 1.0/direction_x1_norm * ( (x2-x1).Cross(direction_x1) ).Mag();
1073// cout << "CylinderRadius= " << CylinderRadius << " TubeDistance= " << TubeDistance << endl;
1074 if (TubeDistance>CylinderRadius) {
1075 return kFALSE;
1076 }
1077
1078 return kTRUE;
1079}
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176

◆ Initialize()

void EdbShowAlg_CA::Initialize ( )
virtual

Reimplemented from EdbShowAlg.

829{
830 return;
831}

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