FEDRA emulsion software from the OPERA Collaboration
EdbShowAlg_OI Class Reference

#include <EdbShowAlg.h>

Inheritance diagram for EdbShowAlg_OI:
Collaboration diagram for EdbShowAlg_OI:

Public Member Functions

 ClassDef (EdbShowAlg_OI, 1)
 
 EdbShowAlg_OI ()
 
void Execute ()
 
void Finalize ()
 
void Initialize ()
 
virtual ~EdbShowAlg_OI ()
 
- 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)
 

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_OI()

EdbShowAlg_OI::EdbShowAlg_OI ( )
1175{
1176 // Default Constructor
1177 Log(2,"EdbShowAlg_OI::EdbShowAlg_OI()","Default Constructor (OfficialImplementation, as it is in libShower module)");
1178
1179// Reset all:
1180 Set0();
1181
1182 eAlgName="OI";
1183 eAlgValue=4; // see default.par_SHOWREC for labeling (labeling identical with ShowRec program)
1184 eParaN =4; // depends on the algorithm, how much parameters it has got.
1185
1186 eParaValue[0]=700.0;
1187 eParaString[0]="CylinderRadius"; // micrometer
1188 eParaValue[1]=0.025;
1189 eParaString[1]="ConeAngle"; // MILLIradiant
1190 eParaValue[2]=150.0;
1191 eParaString[2]="DeltaRToPreceedingBT"; // micrometer
1192 eParaValue[3]=0.13;
1193 eParaString[3]="DeltaThetaComponentwiseToPreceedingBT"; // radiant
1194
1195// 700 0.025 150 0.13
1196 Log(2,"EdbShowAlg_OI::EdbShowAlg_OI()","Default Constructor (OfficialImplementation, as it is in libShower module)...done.");
1197}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
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_OI()

EdbShowAlg_OI::~EdbShowAlg_OI ( )
virtual
1202{
1203 Log(2,"EdbShowAlg_OI::~EdbShowAlg_OI()","Default Destructor called.");
1204// Default Destructor
1205}

Member Function Documentation

◆ ClassDef()

EdbShowAlg_OI::ClassDef ( EdbShowAlg_OI  ,
 
)

◆ Execute()

void EdbShowAlg_OI::Execute ( )
virtual

WARNING ... Somehow the calculation of "STEP" is not correct, when having data ordered patterns, it might go the wrong direction. Better use GetPatternNext(..,+1) to be sure to get the next pattern with increasing Z positions. To be checked: does it work also on Simulated ORFEO data? To be checked: does it work also on the testbeam data?

newActualPID=ActualPID+STEP; // old, mixes up downstream and upstream directions

Reimplemented from EdbShowAlg.

1217{
1218 Log(2,"EdbShowAlg_OI::Execute()","DOING MAIN SHOWER RECONSTRUCTION HERE");
1219
1220 EdbSegP* InBT;
1221 EdbSegP* Segment;
1222 EdbShowerP* RecoShower;
1223
1224 Bool_t StillToLoop=kTRUE;
1225 Int_t ActualPID;
1226 Int_t newActualPID;
1227 Int_t STEP=-1;
1228 Int_t NLoopedPattern=0;
1229
1230 Int_t newZ_next_ascending_Z=0;
1231 Int_t newPID_next_ascending_Z=0;
1232
1233
1242// cout << "eFirstPlate_eAliPID = " << eFirstPlate_eAliPID << endl;
1243// cout << "eLastPlate_eAliPID = " << eLastPlate_eAliPID << endl;
1244 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg_OI::Execute--- STEP for patternloop direction = " << STEP << endl;
1245
1246
1247//--- Loop over InBTs:
1248 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg_OI::Execute Loop over InBTs:" << endl;
1249
1250
1251// Since eInBTArray is filled in ascending ordering by zpositon
1252// We use the descending loop to begin with BT with lowest z first.
1253
1254 for (Int_t i=eInBTArrayN-1; i>=0; --i) {
1255
1256// CounterOutPut, only at gEDBDEBUGLEVEL==2
1257 if (gEDBDEBUGLEVEL==2) if ((i%100)==0) cout << eInBTArrayN <<" Initiator Basetracks (InBT) in total, still to do: " << Form("%06d",i) << "\r\r\r\r\r\r" << flush;
1258
1259//-----------------------------------
1260// 1) Make local_gAli with cut parameters:
1261//-----------------------------------
1262
1263// Create new EdbShowerP Object for storage;
1264// See EdbShowerP why I have to call the Constructor as "unique" ideficable value
1265 RecoShower = new EdbShowerP(i,eAlgValue,eActualAlgParametersetNr);
1266
1267// Get InitiatorBT from eInBTArray
1268 InBT=(EdbSegP*)eInBTArray->At(i);
1269 if (gEDBDEBUGLEVEL>2) InBT->PrintNice();
1270 Float_t X0=InBT->X();
1271 Float_t Y0=InBT->Y();
1272
1273// Add InBT to RecoShower:
1274// This has to be done, since by definition the first BT in the RecoShower is the InBT.
1275// Otherwise, also the definition of shower axis and transversal profiles is wrong!
1276 RecoShower -> AddSegment(InBT);
1277 if (gEDBDEBUGLEVEL>4) cout << "Segment (InBT) " << Segment << " was added to RecoShower." << endl;
1278
1279// Transform (make size smaller, extract only events having same MC) the eAli object:
1280 Transform_eAli(InBT);
1281 if (gEDBDEBUGLEVEL>3) eAli_Sub->Print();
1282
1283//-----------------------------------
1284// 2) Loop over (whole) eAli, check BT for Cuts
1285// eAli_Sub
1286//-----------------------------------
1287 ActualPID= InBT->PID() ;
1288 newActualPID= InBT->PID() ;
1289
1290 while (StillToLoop) {
1291 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_OI::Execute--- --- Doing patterloop " << ActualPID << " for patterns Z position=" << eAli_Sub->GetPattern(ActualPID)->Z() << endl;
1292
1293 for (Int_t btloop_cnt=0; btloop_cnt<eAli_Sub->GetPattern(ActualPID)->GetN(); ++btloop_cnt) {
1294
1295 Segment = (EdbSegP*)eAli_Sub->GetPattern(ActualPID)->GetSegment(btloop_cnt);
1296 if (gEDBDEBUGLEVEL>4) Segment->PrintNice();
1297
1298 // Quick Continue Statements:
1299 if ( Segment->MCEvt() > 0 ) if ( Segment->MCEvt()!=InBT->MCEvt() ) continue; // MCEvtNr (>0) or BgMCNr (-999)
1300 if ( Abs(Segment->X()-X0) > 7000 ) continue;
1301 if ( Abs(Segment->Y()-Y0) > 7000 ) continue;
1302 // Now apply cut conditions: OI OfficialImplementation Alg --------------------
1303 if ( !IsInConeTube(Segment, InBT, eParaValue[0], eParaValue[1]) ) continue;
1304 if ( !FindPrecedingBTs(Segment, InBT, eAli_Sub, RecoShower)) continue;
1305 // end of cut conditions: OI OfficialImplementation Alg --------------------
1306
1307
1308// If we arrive here, Basetrack Segment has passed criteria
1309// and is then added to the RecoShower:
1310// Check if its not the InBT which is already added:
1311// Use compare function of EdbSegP, which is safer instead of comparison of adresses
1312 // (in eAli_Sub, Adresses can be new and then it doesnt work)
1313 if ( Segment->Compare(InBT) == 0) {
1314// cout << "Segment->X()==InBT->X()&&Segment->Y()==InBT->Y() segments have same coord, so dont add" << endl;
1315// cout << "check the compare-function" << endl;
1316// cout << "Segment->Compare(InBT) = " << Segment->Compare(InBT) << endl;
1317// int EdbSegP::Compare( const TObject *obj ) const
1318 ; // is InBT, do nothing;
1319 }
1320 else {
1321 RecoShower -> AddSegment(Segment);
1322 }
1323 if (gEDBDEBUGLEVEL>4) cout << "Segment " << Segment << " was added at &Segment : " << &Segment << endl;
1324
1325 } // of btloop_cnt
1326
1327
1328 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_OI::Execute--- --- ActualPID= " << ActualPID << " at Z= " << eAli_Sub->GetPattern(ActualPID)->Z() << " done. Reconstructed shower has up to now: " << RecoShower->N() << " Segments." << endl;
1329
1330// // // // // // // // // // // // // // // // // XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
1331 // cout << " DETERMINE NEXT PATTERN in ASCENDING (INCREASING) Z DIRECTION now: " << endl;
1332 EdbPattern* currentPat = eAli_Sub->GetPattern(ActualPID);
1333 EdbPattern* newPat_next_ascending_Z = eAli_Sub->GetPatternNext(currentPat->Z(),1); // in ascending Z
1334
1335 // cout << "What if there is no next ascending Z pattern? Check it here...... " << endl;
1336 // cout << "newPat_next_ascending_Z = " << newPat_next_ascending_Z << endl;
1337
1338 if ( NULL == newPat_next_ascending_Z ) {
1339 if (gEDBDEBUGLEVEL>2) {
1340 cout << "newPat_next_ascending_Z is NULL pointer ! So it is nonexisting! (possibly last plate); " << endl;
1341 cout << "EdbShowAlg_OI::Execute--- ---Stop Loop since: newPat_next_ascending_Z is NULL pointer "<<endl;
1342 }
1343 StillToLoop=kFALSE;
1344 }
1345 else {
1346 newZ_next_ascending_Z = newPat_next_ascending_Z->Z();
1347 newPID_next_ascending_Z = newPat_next_ascending_Z->PID();
1348 }
1349
1350 // cout << "Set new values for the next ascending pattern: newZ_next_ascending_Z = " << newZ_next_ascending_Z << endl;
1351 // cout << "Set new values for the next ascending pattern: newPID_next_ascending_Z = " << newPID_next_ascending_Z << endl;
1352 // cout << " NEXT PATTERN in ASCENDING (INCREASING) Z DIRECTION has PID()= " << newPID_next_ascending_Z << " and a z-Position of Z= " << newZ_next_ascending_Z << endl;
1353// // // // // // // // // // // // // // // // // XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
1354
1355//------------
1357 newActualPID=newPID_next_ascending_Z;
1358 ++NLoopedPattern;
1359
1360 if (gEDBDEBUGLEVEL>3) {
1361 cout << "EdbShowAlg_OI::Execute--- --- newActualPID= " << newActualPID << endl;
1362 cout << "EdbShowAlg_OI::Execute--- --- newActualZ= " << eAli_Sub->GetPattern(newActualPID)->Z() << endl;
1363 cout << "EdbShowAlg_OI::Execute--- --- NLoopedPattern= " << NLoopedPattern << endl;
1364 cout << "EdbShowAlg_OI::Execute--- --- eNumberPlate_eAliPID= " << eNumberPlate_eAliPID << endl;
1365 cout << "EdbShowAlg_OI::Execute--- --- StillToLoop= " << StillToLoop << endl;
1366 cout << "EdbShowAlg_OI::Execute--- --- We have to decide, if to loop on for the next pattern... (out of the " << eAli_Sub->Npatterns() << " patterns)." << endl;
1367 }
1368
1369 /* NOT NECESSARY ANYMORE SINCE WE DO THE OUT OF BOUNDS CHECK IN ANOTHER WAY .....
1370 * TO BE REMOVED ...
1371 // This if holds in the case of STEP== +1
1372 if (STEP==1) {
1373 if (newActualPID>eLastPlate_eAliPID) StillToLoop=kFALSE;
1374 if (newActualPID>eLastPlate_eAliPID) cout << "EdbShowAlg_OI::Execute--- ---Stop Loop since: newActualPID>eLastPlate_eAliPID"<<endl;
1375 }
1376 // This if holds in the case of STEP== -1
1377 if (STEP==-1) {
1378 if (newActualPID<eLastPlate_eAliPID) StillToLoop=kFALSE;
1379 if (newActualPID<eLastPlate_eAliPID && gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_OI::Execute--- ---Stop Loop since: newActualPID<eLastPlate_eAliPID"<<endl;
1380 }
1381 */
1382
1383 // What is the constraint? 0 <= newActualPID <= Npatterns-1, these patterns are valid....
1384 if (newActualPID<0)
1385 { StillToLoop=kFALSE;
1386 // cout << CStop Loop since: newActualPID<0"<<endl;
1387 }
1388 if (newActualPID>=eNumberPlate_eAliPID-1) {
1389 StillToLoop=kFALSE;
1390 // cout << "EdbShowAlg_OI::Execute--- ---Stop Loop since: newActualPID>=eNumberPlate_eAliPID-1 [0 <= newActualPID <= Npatterns-1]"<<endl;
1391 }
1392
1393// This if holds general, since eNumberPlate_eAliPID is not dependent of the structure of the gAli subject:
1394 if (NLoopedPattern>eNumberPlate_eAliPID) StillToLoop=kFALSE;
1395 if (NLoopedPattern>eNumberPlate_eAliPID && gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_OI::Execute--- ---Stop Loop since: NLoopedPattern>eNumberPlate_eAliPID"<<endl;
1396
1397 ActualPID=newActualPID;
1398 } // of // while (StillToLoop)
1399
1400
1401// Obligatory when Shower Reconstruction is finished:
1402// Set Numbers, plates etc. correct
1403 RecoShower ->Update();
1404
1405 // Debung Print Information:
1406 if (gEDBDEBUGLEVEL>2) RecoShower ->PrintBasics();
1407 if (gEDBDEBUGLEVEL>3) RecoShower ->PrintNice();
1408 if (gEDBDEBUGLEVEL>3) RecoShower ->PrintSegments();
1409
1410
1411// if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_OI::Execute--- Before adding to array delete the histograms...by finalize() of shower."<<endl;
1412// RecoShower ->Finalize();
1413// if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_OI::Execute--- Finalize done..."<<endl;
1414
1415
1416// Add Shower Object to Shower Reco Array.
1417// Not, if its empty or containing only one BT:
1418// Note, this function add also eRecoShowerArrayN directly.
1419 if (RecoShower->N()>1) {
1420 AddRecoShowerArray(RecoShower);
1421 }
1422 else {
1423 // Delete the created shower object. Saves _lots_ of memory and speeds up reconstruction for large #InBTs!
1424 // TODO Also Include this in the other RecoAlgs ...
1425 Log(3,"EdbShowAlg_OI::Execute()","InBT # %d RecoShower->N()<=1: Dont add shower to RecoShowerArray. Delete this shower for memory safing.");
1426 delete RecoShower;
1427 RecoShower=NULL;
1428 }
1429
1430
1431
1432// Set back loop values:
1433 StillToLoop=kTRUE;
1434 NLoopedPattern=0;
1435 } // of // for (Int_t i=eInBTArrayN-1; i>=0; --i)
1436 // Finished Loop over all InBTs here.
1437
1438
1439
1440 // Set ID of the shower, otherwise all reco showers have ID=0:
1441 // TO BE DONE: auch in den anderen reco-Algorithmen die ID() nachträglich setzen!
1443
1444 // Write the MetaData of the Algorithm Values into the showers, that
1445 // information on reconstruction parameters can be retrieved
1446 // (usefull when more parametersets are reconstructed)
1447 // TO BE DONE::: Also put this in the other Execute Functions TODO
1448 cout << "================================================================================="<<endl;
1450 cout << "================================================================================="<<endl;
1451
1452 if (gEDBDEBUGLEVEL==2) cout << endl << "Done reco!" << endl;
1453 cout << "EdbShowAlg_OI::eRecoShowerArray = " << eRecoShowerArray << endl;
1454 cout << "EdbShowAlg_OI::eRecoShowerArrayN = " << eRecoShowerArrayN << endl;
1455
1456 Log(2,"EdbShowAlg_OI::Execute()","DOING MAIN SHOWER RECONSTRUCTION HERE...done.");
1457 return;
1458}
brick X0
Definition: RecDispMC.C:112
Definition: EdbPattern.h:273
int PID() const
Definition: EdbPattern.h:320
Int_t Npatterns() const
Definition: EdbPattern.h:366
void Print() const
Definition: EdbPattern.cxx:1693
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
EdbPattern * GetPatternNext(float z, int dir) const
Definition: EdbPattern.cxx:1857
Definition: EdbSegP.h:21
Float_t X() const
Definition: EdbSegP.h:173
Int_t Compare(const TObject *obj) const
Definition: EdbSegP.cxx:420
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:1479
void UpdateShowerIDs()
Definition: EdbShowAlg.cxx:479
void Transform_eAli(EdbSegP *InitiatorBT, Float_t ExtractSize)
Definition: EdbShowAlg.cxx:107
void AddRecoShowerArray(EdbShowerP *shower)
Definition: EdbShowAlg.h:120
Bool_t IsInConeTube(EdbSegP *sTest, EdbSegP *sStart, Double_t CylinderRadius, Double_t ConeAngle)
Definition: EdbShowAlg.cxx:248
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
Int_t eNumberPlate_eAliPID
Definition: EdbShowAlg.h:68
Int_t eRecoShowerArrayN
Definition: EdbShowAlg.h:71
EdbPVRec * eAli_Sub
Definition: EdbShowAlg.h:74
void UpdateShowerMetaData()
Definition: EdbShowAlg.cxx:492
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_OI::Finalize ( )
virtual

Reimplemented from EdbShowAlg.

1464{
1465 return;
1466}

◆ FindPrecedingBTs()

Bool_t EdbShowAlg_OI::FindPrecedingBTs ( EdbSegP s,
EdbSegP InBT,
EdbPVRec gAli,
EdbShowerP shower 
)
private
1480{
1481 if (gEDBDEBUGLEVEL>3) {
1482 cout << "EdbShowAlg_OI::FindPrecedingBTs Find BTs to be connected with the Test BT:" << endl;
1483 cout << "s_TestBT X,Y,Z " << s->X() << " " << s->Y() << " " << s->Z() << " versus InBT X,Y,Z " << InBT->X() << " " << InBT->Y() << " " << InBT->Z() << endl;
1484 }
1485
1486 EdbSegP* s_TestBT;
1487 Int_t nentries=shower->N();
1488 Double_t dZ;
1489
1490// In case the shower is empty we do not have to search for a preceeding BT:
1491 if (nentries==0) return kTRUE;
1492
1493// We do not test BaseTracks which are before the initiator BaseTrack!
1494 if (s->Z()<InBT->Z()) {
1495 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_OI::FindPrecedingBTs(): s->Z()<InBT->Z()..: We do not test BaseTracks which are before the initiator BaseTrack! return kFALSE;" << endl;
1496 return kFALSE;
1497 }
1498
1499// For the very first Z position we do not test.
1500// If testBT has Preceeders, only if it it has a BT around (case for e+e- coming from gammma):
1501// Take 50microns and 80mrad in (dR/dT) around.
1502// This does not affect the normal results, but helps for
1503// events which may have a second BT close to InBT (like in e+e-)
1504 if (s->Z()==InBT->Z()) {
1505 if (DeltaThetaComponentwise(s, InBT) < 0.08 && DeltaR_WithoutPropagation(s, InBT) < 50.0 ) return kTRUE;
1506 }
1507
1508 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_OI::FindPrecedingBTs(): Testing " << nentries << " entries of the shower." << endl;
1509
1510 for (Int_t i=nentries-1; i>=0; --i) {
1511 s_TestBT = (EdbSegP*)( shower->GetSegment(i) );
1512
1513 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_OI::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;
1514 if (gEDBDEBUGLEVEL>3) s_TestBT->PrintNice();
1515 if (gEDBDEBUGLEVEL>3) s->PrintNice();
1516
1517
1518 dZ=TMath::Abs(s_TestBT->Z()-s->Z());
1519 if (dZ<30) continue; // Exclude the case of same Zpositions...
1521 if (dZ>(3.5*(1300.0+30.0))) continue; // Exclude the case of more than 3 plates before... // modified
1522 // often, experimentally observed, dZ approx 1350 ... so the upper cut really does not
1523 // look for 3 plates, but only for two ...: seems to work now. (24.Oct.2016)
1524
1525
1526 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_OI::FindPrecedingBTs(): Checking dT,dR and dZ for i: " << i << " " << DeltaThetaComponentwise(s, s_TestBT) << " " << DeltaR_WithPropagation(s, s_TestBT) << " "<<dZ << endl;
1527
1528 if (DeltaThetaComponentwise(s, s_TestBT) > eParaValue[3] ) continue;
1529 if (DeltaR_WithPropagation(s, s_TestBT) > eParaValue[2]) continue;
1530
1531 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_OI::FindPrecedingBTs(): Checking dT,dR and dZ for i: " << i << " " << DeltaThetaComponentwise(s, s_TestBT) << " " << DeltaR_WithPropagation(s, s_TestBT) << " "<<dZ << " ok!"<<endl;
1532 return kTRUE;
1533 }
1534//---------------------------------------------
1535 return kFALSE;
1536}
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 DeltaThetaComponentwise(EdbSegP *s1, EdbSegP *s2)
Definition: EdbShowAlg.cxx:441
Double_t DeltaR_WithPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowAlg.cxx:411
EdbSegP * GetSegment(int i) const
Definition: EdbShowerP.h:435
int nentries
Definition: check_shower.C:40
s
Definition: check_shower.C:55

◆ Initialize()

void EdbShowAlg_OI::Initialize ( )
virtual

Reimplemented from EdbShowAlg.

1211{
1212 return;
1213}

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