NCFS-Pack
A generic (astro)particle physics analysis framework
 
Loading...
Searching...
No Matches
NcJet.cxx
Go to the documentation of this file.
1
31
33
114
115#include "NcJet.h"
116#include "Riostream.h"
117
118ClassImp(NcJet); // Class implementation to enable ROOT I/O
119
122{
130
131 Init();
132 Reset();
133 SetNtinit();
134}
135
137{
143
144 fTracks=0;
145 fNtinit=0;
146 fTrackCopy=0;
147 fRef=0;
148 fSelected=0;
149 fEscale=1;
150}
151
152NcJet::NcJet(Int_t n) : TNamed(),Nc4Vector()
153{
160
161 Init();
162 Reset();
163 if (n > 0)
164 {
165 SetNtinit(n);
166 }
167 else
168 {
169 cout << endl;
170 cout << " *NcJet* Initial max. number of tracks entered : " << n << endl;
171 cout << " This is invalid. Default initial maximum will be used." << endl;
172 cout << endl;
173 SetNtinit();
174 }
175}
176
178{
184
185 if (fTracks)
186 {
187 delete fTracks;
188 fTracks=0;
189 }
190 if (fRef)
191 {
192 delete fRef;
193 fRef=0;
194 }
195 if (fSelected)
196 {
197 delete fSelected;
198 fSelected=0;
199 }
200}
201
202void NcJet::SetOwner(Bool_t own)
203{
226
227 Int_t mode=1;
228 if (!own) mode=0;
229 if (fTracks) fTracks->SetOwner(own);
230 fTrackCopy=mode;
231}
232
233NcJet::NcJet(const NcJet& j) : TNamed(j),Nc4Vector(j)
234{
240
242 fNtmax=j.fNtmax;
243 fQ=j.fQ;
245 fNtrk=j.fNtrk;
248 if (j.fRef) fRef=new NcPositionObj(*(j.fRef));
249
250 fSelected=0;
251
252 fTracks=0;
253 if (fNtrk)
254 {
255 fTracks=new TObjArray(fNtmax);
256 if (fTrackCopy) fTracks->SetOwner();
257 }
258
259 for (Int_t i=1; i<=fNtrk; i++)
260 {
261 NcTrack* tx=j.GetTrack(i);
262 if (fTrackCopy)
263 {
264 fTracks->Add(tx->Clone());
265 }
266 else
267 {
268 fTracks->Add(tx);
269 }
270 }
271}
272
273void NcJet::SetNtinit(Int_t n)
274{
280
281 fNtinit=n;
282 fNtmax=n;
283
284 if (fTracks)
285 {
286 delete fTracks;
287 fTracks=0;
288 }
289 if (fRef)
290 {
291 delete fRef;
292 fRef=0;
293 }
294}
295
297{
305
306 fNtrk=0;
307 fQ=0;
308 fUserId=0;
309 SetZero();
310 if (fNtinit > 0) SetNtinit(fNtinit);
311}
312
314{
339
340
341 AddTrack(t,1);
342}
343
344void NcJet::AddTrack(NcTrack& t,Int_t copy)
345{
367
368 if (!fTracks)
369 {
370 fTracks=new TObjArray(fNtmax);
371 if (fTrackCopy) fTracks->SetOwner();
372 }
373 else if (!fTrackCopy || !copy) // Check if this track is already present
374 {
375 for (Int_t i=0; i<fNtrk; i++)
376 {
377 NcTrack* tx=(NcTrack*)fTracks->At(i);
378 if (tx == &t) return;
379 }
380 }
381
382 if (fNtrk == fNtmax) // Check if maximum track number is reached
383 {
385 fTracks->Expand(fNtmax);
386 }
387
388 // Add the track to this jet
389 fNtrk++;
390 if (fTrackCopy && copy)
391 {
392 fTracks->Add(t.Clone());
393 }
394 else
395 {
396 fTracks->Add(&t);
397 }
398
399 fQ+=t.GetCharge();
400
401 Nc4Vector p4=(Nc4Vector)t;
402 Float_t tscale=t.GetEscale();
403 if ((tscale/fEscale > 1.1) || (fEscale/tscale > 1.1)) p4=p4*(tscale/fEscale);
404 (*this)+=p4;
405
406}
407
408void NcJet::Data(TString f,TString u)
409{
422
423 const char* name=GetName();
424 const char* title=GetTitle();
425
426 cout << " *NcJet::Data*";
427 if (strlen(name)) cout << " Name : " << GetName();
428 if (strlen(title)) cout << " Title : " << GetTitle();
429 cout << endl;
430 cout << " Id : " << fUserId << " Invmass : " << GetInvmass() << " Charge : " << fQ
431 << " Momentum : " << GetMomentum() << " Energy scale : " << fEscale << " GeV" << endl;
432
433 ShowTracks(0);
434
435 if (fV.HasVector()) Nc4Vector::Data(f,u);
436}
437
438void NcJet::List(TString f,TString u,TObjArray* tracks)
439{
458
459 Data(f,u); // Information of the current jet
460 if (fRef) { cout << " Ref-point :"; fRef->Data(f,u); }
461
462 Int_t ntk=0;
463 if (tracks)
464 {
465 ntk=tracks->GetEntries();
466 }
467 else
468 {
469 ntk=GetNtracks();
470 }
471
472 if (!ntk)
473 {
474 cout << " *NcJet::List* No tracks are present." << endl;
475 return;
476 }
477
478 // The tracks of this jet
479 TObject* obj;
480 NcTrack* tx;
481 for (Int_t it=1; it<=ntk; it++)
482 {
483 if (tracks)
484 {
485 obj=tracks->At(it-1);
486 if (!obj) continue;
487 if (!obj->InheritsFrom("NcTrack")) continue;
488 tx=(NcTrack*)obj;
489 }
490 else
491 {
492 tx=GetTrack(it);
493 if (!tx) continue;
494 }
495 cout << " ---Track no. " << it << endl;
496 cout << " ";
497 tx->Data(f,u);
498 }
499}
500
501void NcJet::ListAll(TString f,TString u,TObjArray* tracks)
502{
521
522 Data(f,u); // Information of the current jet
523 if (fRef) { cout << " Ref-point :"; fRef->Data(f,u); }
524
525 Int_t ntk=0;
526 if (tracks)
527 {
528 ntk=tracks->GetEntries();
529 }
530 else
531 {
532 ntk=GetNtracks();
533 }
534
535 if (!ntk)
536 {
537 cout << " *NcJet::ListAll* No tracks are present." << endl;
538 return;
539 }
540
541 // The tracks of this jet
542 TObject* obj;
543 NcTrack* tx;
544 for (Int_t it=1; it<=ntk; it++)
545 {
546 if (tracks)
547 {
548 obj=tracks->At(it-1);
549 if (!obj) continue;
550 if (!obj->InheritsFrom("NcTrack")) continue;
551 tx=(NcTrack*)obj;
552 }
553 else
554 {
555 tx=GetTrack(it);
556 if (!tx) continue;
557 }
558 cout << " ---Track no. " << it << endl;
559 cout << " ";
560 tx->ListAll(f,u);
561 }
562}
563
564Int_t NcJet::GetNtracks(Int_t idmode,Int_t chmode,Int_t pcode)
565{
584
585 Int_t n=0;
586 if (idmode==0 && chmode==2 && pcode==0)
587 {
588 return fNtrk;
589 }
590 else
591 {
592 TObjArray tracks;
593 GetTracks(idmode,chmode,pcode,&tracks);
594 n=tracks.GetEntries();
595 return n;
596 }
597}
598
599Int_t NcJet::GetNtracks(TString name,Int_t mode)
600{
619
620 Int_t n=0;
621
622 TObjArray tracks;
623 GetTracks(name,mode,&tracks);
624 n=tracks.GetEntries();
625 return n;
626}
627
628Double_t NcJet::GetEnergy(Float_t scale)
629{
642
643 Double_t E=GetScalar();
644 if (E>0)
645 {
646 if (scale>0)
647 {
648 E*=fEscale/scale;
649 fDresult*=fEscale/scale;
650 }
651 return E;
652 }
653 else
654 {
655 return 0;
656 }
657}
658
659Double_t NcJet::GetMomentum(Float_t scale)
660{
673
674 Double_t norm=fV.GetNorm();
675 fDresult=fV.GetResultError();
676 if (scale>0)
677 {
678 norm*=fEscale/scale;
679 fDresult*=fEscale/scale;
680 }
681 return norm;
682}
683
684Nc3Vector NcJet::Get3Momentum(Float_t scale) const
685{
697
699 if (scale>0) p*=fEscale/scale;
700 return p;
701}
702
703Double_t NcJet::GetInvmass(Float_t scale)
704{
717
718 Double_t inv=Dot(*this);
719 Double_t dinv=GetResultError();
720 Double_t dm=0;
721 if (inv >= 0)
722 {
723 Double_t m=sqrt(inv);
724 if (m) dm=dinv/(2.*m);
725 if (scale>0)
726 {
727 m*=fEscale/scale;
728 dm*=fEscale/scale;
729 }
730 fDresult=dm;
731 return m;
732 }
733 else
734 {
735 fDresult=dm;
736 return 0;
737 }
738}
739
740Float_t NcJet::GetCharge() const
741{
747
748 return fQ;
749}
750
752{
758
759 if (!fTracks) return 0;
760
761 if (i<=0 || i>fNtrk)
762 {
763 cout << " *NcJet*::GetTrack* Invalid argument i : " << i
764 << " Ntrk = " << fNtrk << endl;
765 return 0;
766 }
767 else
768 {
769 return (NcTrack*)fTracks->At(i-1);
770 }
771}
772
774{
780
781 if (!fTracks) return 0;
782
783 NcTrack* tx=0;
784 for (Int_t i=0; i<fNtrk; i++)
785 {
786 tx=(NcTrack*)fTracks->At(i);
787 if (id == tx->GetId()) return tx;
788 }
789 return 0; // No matching id found
790}
791
792TObjArray* NcJet::GetTracks(Int_t idmode,Int_t chmode,Int_t pcode,TObjArray* tracks)
793{
843
844 if (tracks)
845 {
846 tracks->Clear();
847 }
848 else
849 {
850 if (fSelected)
851 {
852 fSelected->Clear();
853 }
854 else
855 {
856 fSelected=new TObjArray();
857 }
858 }
859
860 if (!fTracks)
861 {
862 if (tracks)
863 {
864 return 0;
865 }
866 else
867 {
868 return fSelected;
869 }
870 }
871
872 NcTrack* tx=0;
873 Int_t code=0;
874 Int_t id=0;
875 Float_t q=0;
876 for (Int_t i=0; i<fNtrk; i++)
877 {
878 tx=(NcTrack*)fTracks->At(i);
879 if (!tx) continue;
880
881 code=tx->GetParticleCode();
882 if (pcode && abs(pcode)!=abs(code)) continue;
883
884 id=tx->GetId();
885 if (idmode==-1 && id>=0) continue;
886 if (idmode==1 && id<=0) continue;
887
888 q=tx->GetCharge();
889 if (chmode==-1 && q>=0) continue;
890 if (chmode==0 && fabs(q)>1e-10) continue;
891 if (chmode==1 && q<=0) continue;
892 if (chmode==3 && fabs(q)<1e-10) continue;
893
894 if (tracks)
895 {
896 tracks->Add(tx);
897 }
898 else
899 {
900 fSelected->Add(tx);
901 }
902 }
903
904 if (tracks)
905 {
906 return 0;
907 }
908 else
909 {
910 return fSelected;
911 }
912}
913
914TObjArray* NcJet::GetTracks(TString name,Int_t mode,TObjArray* tracks)
915{
940
941 if (tracks)
942 {
943 tracks->Clear();
944 }
945 else
946 {
947 if (fSelected)
948 {
949 fSelected->Clear();
950 }
951 else
952 {
953 fSelected=new TObjArray();
954 }
955 }
956
957 if (!fTracks)
958 {
959 if (tracks)
960 {
961 return 0;
962 }
963 else
964 {
965 return fSelected;
966 }
967 }
968
969 NcTrack* tx=0;
970 TString s;
971 Int_t isel=0;
972 for (Int_t i=0; i<fNtrk; i++)
973 {
974 tx=(NcTrack*)fTracks->At(i);
975 if (!tx) continue;
976
977 s=tx->GetName();
978 isel=0;
979 if (name=="*") isel=1;
980 if (!mode && s==name) isel=1;
981 if (mode==1 && s.Contains(name.Data())) isel=1;
982
983 if (!isel) continue;
984
985 if (tracks)
986 {
987 tracks->Add(tx);
988 }
989 else
990 {
991 fSelected->Add(tx);
992 }
993 }
994
995 if (tracks)
996 {
997 return 0;
998 }
999 else
1000 {
1001 return fSelected;
1002 }
1003}
1004
1006{
1013
1014 if (!t || !fTracks) return;
1015
1016 RemoveTrack(t,1);
1017}
1018
1019void NcJet::RemoveTracks(TString name,Int_t mode)
1020{
1038
1039 if (!fTracks) return;
1040
1041 TObjArray arr;
1042 GetTracks(name,mode,&arr);
1043
1044 Int_t ntk=arr.GetEntries();
1045 if (!ntk) return;
1046
1047 NcTrack* tx=0;
1048 for (Int_t i=0; i<ntk; i++)
1049 {
1050 tx=(NcTrack*)arr.At(i);
1051 if (tx) RemoveTrack(tx,0);
1052 }
1053 fTracks->Compress();
1054 fNtrk=fTracks->GetEntries();
1055}
1056
1057void NcJet::RemoveTracks(Int_t idmode,Int_t chmode,Int_t pcode)
1058{
1067
1068 if (!fTracks) return;
1069
1070 TObjArray arr;
1071 GetTracks(idmode,chmode,pcode,&arr);
1072
1073 Int_t ntk=arr.GetEntries();
1074 if (!ntk) return;
1075
1076 NcTrack* tx=0;
1077 for (Int_t i=0; i<ntk; i++)
1078 {
1079 tx=(NcTrack*)arr.At(i);
1080 if (tx) RemoveTrack(tx,0);
1081 }
1082 fTracks->Compress();
1083 fNtrk=fTracks->GetEntries();
1084}
1085
1087{
1094
1095 if (!told || !tnew || !fTracks) return;
1096
1097 Int_t index=fTracks->IndexOf(told);
1098
1099 if (index<0 || index>=fTracks->GetEntries()) return;
1100
1101 RemoveTrack(told,0);
1102
1103 if (fTrackCopy)
1104 {
1105 fTracks->AddAt(tnew->Clone(),index);
1106 }
1107 else
1108 {
1109 fTracks->AddAt(tnew,index);
1110 }
1111
1112 fQ+=tnew->GetCharge();
1113
1114 Nc4Vector p4=(Nc4Vector)(*tnew);
1115 Float_t tscale=tnew->GetEscale();
1116 if ((tscale/fEscale > 1.1) || (fEscale/tscale > 1.1)) p4=p4*(tscale/fEscale);
1117 (*this)+=p4;
1118}
1119
1120void NcJet::RemoveTrack(NcTrack* t,Int_t compress)
1121{
1135
1136 if (!t || !fTracks) return;
1137
1138 TObject* obj=fTracks->Remove(t);
1139 if (obj)
1140 {
1141 fQ-=t->GetCharge();
1142 Nc4Vector p4=(Nc4Vector)(*t);
1143 Float_t tscale=t->GetEscale();
1144 if ((tscale/fEscale > 1.1) || (fEscale/tscale > 1.1)) p4=p4*(tscale/fEscale);
1145 (*this)-=p4;
1146 if (compress)
1147 {
1148 fTracks->Compress();
1149 fNtrk=fTracks->GetEntries();
1150 }
1151 if (fTracks->IsOwner())
1152 {
1153 delete t;
1154 t=0;
1155 }
1156 }
1157}
1158
1159void NcJet::ShowTracks(Int_t mode,TString f,TString u,TObjArray* tracks)
1160{
1181
1182 Int_t ntk=0;
1183 if (tracks)
1184 {
1185 ntk=tracks->GetEntries();
1186 }
1187 else
1188 {
1189 ntk=GetNtracks();
1190 }
1191
1192 if (!ntk)
1193 {
1194 cout << " No tracks are present." << endl;
1195 return;
1196 }
1197
1198 if (!mode)
1199 {
1200 cout << " There are " << ntk << " tracks available." << endl;
1201 return;
1202 }
1203
1204 if (mode==1)
1205 {
1206 cout << " The following " << ntk << " tracks are available :" << endl;
1207 TObject* obj=0;
1208 NcTrack* tx=0;
1209 for (Int_t i=1; i<=ntk; i++)
1210 {
1211 if (tracks)
1212 {
1213 obj=tracks->At(i-1);
1214 if (!obj) continue;
1215 if (!obj->InheritsFrom("NcTrack")) continue;
1216 tx=(NcTrack*)obj;
1217 }
1218 else
1219 {
1220 tx=GetTrack(i);
1221 if (!tx) continue;
1222 }
1223 const char* name=tx->GetName();
1224 const char* title=tx->GetTitle();
1225 cout << " Track : " << i;
1226 cout << " Id : " << tx->GetId();
1227 cout << " Q : " << tx->GetCharge() << " m : " << tx->GetMass() << " p : " << tx->GetMomentum();
1228 if (strlen(name)) cout << " Name : " << name;
1229 if (strlen(title)) cout << " Title : " << title;
1230 cout << endl;
1231 }
1232 }
1233
1234 if (mode==2) List(f,u,tracks);
1235
1236 if (mode==3) ListAll(f,u,tracks);
1237}
1238
1239Double_t NcJet::GetPt(Float_t scale)
1240{
1253
1254 Nc3Vector v;
1255 v=GetVecTrans();
1256 Double_t norm=v.GetNorm();
1258 if (scale>0)
1259 {
1260 norm*=fEscale/scale;
1261 fDresult*=fEscale/scale;
1262 }
1263
1264 return norm;
1265}
1266
1267Double_t NcJet::GetPl(Float_t scale)
1268{
1282
1283 Nc3Vector v;
1284 v=GetVecLong();
1285
1286 Double_t pl=v.GetNorm();
1288
1289 Double_t a[3];
1290 v.GetVector(a,"sph");
1291 if (cos(a[1])<0) pl=-pl;
1292 if (scale>0)
1293 {
1294 pl*=fEscale/scale;
1295 fDresult*=fEscale/scale;
1296 }
1297
1298 return pl;
1299}
1300
1301Double_t NcJet::GetEt(Float_t scale)
1302{
1315
1316 Double_t et=GetScaTrans();
1317 if (scale>0)
1318 {
1319 et*=fEscale/scale;
1320 fDresult*=fEscale/scale;
1321 }
1322
1323 return et;
1324}
1325
1326Double_t NcJet::GetEl(Float_t scale)
1327{
1341
1342 Double_t el=GetScaLong();
1343 if (scale>0)
1344 {
1345 el*=fEscale/scale;
1346 fDresult*=fEscale/scale;
1347 }
1348
1349 return el;
1350}
1351
1352Double_t NcJet::GetMt(Float_t scale)
1353{
1366
1367 Double_t pt=GetPt();
1368 Double_t dpt=GetResultError();
1369 Double_t m=GetInvmass();
1370 Double_t dm=GetResultError();
1371
1372 Double_t mt=sqrt(pt*pt+m*m);
1373 Double_t dmt2=0;
1374 if (mt) dmt2=(pow((pt*dpt),2)+pow((m*dm),2))/(mt*mt);
1375
1376 fDresult=sqrt(dmt2);
1377 if (scale>0)
1378 {
1379 mt*=fEscale/scale;
1380 fDresult*=fEscale/scale;
1381 }
1382 return mt;
1383}
1384
1386{
1396
1397 Double_t e=GetEnergy();
1398 Double_t de=GetResultError();
1399 Double_t pl=GetPl();
1400 Double_t dpl=GetResultError();
1401 Double_t sum=e+pl;
1402 Double_t dif=e-pl;
1403
1404 Double_t y=9999,dy2=0;
1405 if (sum && dif) y=0.5*log(sum/dif);
1406
1407 if (sum*dif) dy2=(1./(sum*dif))*(pow((pl*de),2)+pow((e*dpl),2));
1408
1409 fDresult=sqrt(dy2);
1410 return y;
1411}
1412
1414{
1427
1428 if (!fTracks)
1429 {
1430 if (j==0 || j==1)
1431 {
1432 fTrackCopy=j;
1433 }
1434 else
1435 {
1436 cout << "*NcJet::SetTrackCopy* Invalid argument : " << j << endl;
1437 }
1438 }
1439 else
1440 {
1441 cout << "*NcJet::SetTrackCopy* Storage already contained tracks."
1442 << " ==> TrackCopy mode not changed." << endl;
1443 }
1444}
1445
1447{
1455
1456 return fTrackCopy;
1457}
1458
1459void NcJet::SetId(Int_t id)
1460{
1466
1467 fUserId=id;
1468}
1469
1470Int_t NcJet::GetId() const
1471{
1477
1478 return fUserId;
1479}
1480
1482{
1496
1497 if (fRef) delete fRef;
1498 fRef=new NcPositionObj(p);
1499}
1500
1502{
1516
1517 return fRef;
1518}
1519
1520TObjArray* NcJet::SortTracks(Int_t mode,TObjArray* tracks,TObjArray* ordered)
1521{
1562
1563 if (ordered) ordered->Clear();
1564
1565 TObjArray* atracks=tracks;
1566 if (!atracks) atracks=fTracks;
1567
1568 Int_t ntracks=0;
1569 if (atracks) ntracks=atracks->GetEntries();
1570
1571 if (!mode || abs(mode)>13 || !tracks || !ntracks) return 0;
1572
1573 if (ordered)
1574 {
1575 ordered->Expand(ntracks);
1576 }
1577 else
1578 {
1579 if (fSelected) delete fSelected;
1580 fSelected=new TObjArray(ntracks);
1581 }
1582
1583 Double_t val1,val2; // Values of the observable to be tested upon
1584 TObjArray* arr=0; // Generic sorted array pointer to be used in the logic below
1585
1586 if (ordered)
1587 {
1588 arr=ordered;
1589 }
1590 else
1591 {
1592 arr=fSelected;
1593 }
1594
1595 TObject* obj=0;
1596 NcTrack* tx=0;
1597 Int_t nord=0;
1598 for (Int_t i=0; i<ntracks; i++) // Loop over all tracks of the array
1599 {
1600 obj=atracks->At(i);
1601 if (!obj) continue;
1602 if (!obj->InheritsFrom("NcTrack")) continue;
1603 tx=(NcTrack*)obj;
1604
1605 if (nord == 0) // store the first track at the first ordered position
1606 {
1607 nord++;
1608 arr->AddAt(tx,nord-1);
1609 continue;
1610 }
1611
1612 for (Int_t j=0; j<=nord; j++) // put track in the right ordered position
1613 {
1614 if (j == nord) // track has smallest (mode<0) or largest (mode>0) observable value seen so far
1615 {
1616 nord++;
1617 arr->AddAt(tx,j); // add track at the end
1618 break; // go for next track
1619 }
1620
1621 val1=0;
1622 val2=0;
1623
1624 switch (abs(mode))
1625 {
1626 case 1:
1627 val1=tx->GetNsignals();
1628 val2=((NcTrack*)arr->At(j))->GetNsignals();
1629 break;
1630 case 2:
1631 val1=tx->GetEnergy(1);
1632 val2=((NcTrack*)arr->At(j))->GetEnergy(1);
1633 break;
1634 case 3:
1635 val1=tx->GetMomentum(1);
1636 val2=((NcTrack*)arr->At(j))->GetMomentum(1);
1637 break;
1638 case 4:
1639 val1=tx->GetMass(1);
1640 val2=((NcTrack*)arr->At(j))->GetMass(1);
1641 break;
1642 case 5:
1643 val1=tx->GetPt(1);
1644 val2=((NcTrack*)arr->At(j))->GetPt(1);
1645 break;
1646 case 6:
1647 val1=tx->GetPl(1);
1648 val2=((NcTrack*)arr->At(j))->GetPl(1);
1649 break;
1650 case 7:
1651 val1=tx->GetEt(1);
1652 val2=((NcTrack*)arr->At(j))->GetEt(1);
1653 break;
1654 case 8:
1655 val1=tx->GetEl(1);
1656 val2=((NcTrack*)arr->At(j))->GetEl(1);
1657 break;
1658 case 9:
1659 val1=tx->GetMt(1);
1660 val2=((NcTrack*)arr->At(j))->GetMt(1);
1661 break;
1662 case 10:
1663 val1=tx->GetRapidity();
1664 val2=((NcTrack*)arr->At(j))->GetRapidity();
1665 break;
1666 case 11:
1667 val1=tx->GetPseudoRapidity();
1668 val2=((NcTrack*)arr->At(j))->GetPseudoRapidity();
1669 break;
1670 case 12:
1671 val1=tx->GetCharge();
1672 val2=((NcTrack*)arr->At(j))->GetCharge();
1673 break;
1674 case 13:
1675 val1=tx->GetProb();
1676 val2=((NcTrack*)arr->At(j))->GetProb();
1677 break;
1678 }
1679
1680 if (mode<0 && val1 <= val2) continue;
1681 if (mode>0 && val1 >= val2) continue;
1682
1683 nord++;
1684 for (Int_t k=nord-1; k>j; k--) // create empty position
1685 {
1686 arr->AddAt(arr->At(k-1),k);
1687 }
1688 arr->AddAt(tx,j); // put track at empty position
1689 break; // go for next track
1690 }
1691 }
1692
1693 if (ordered)
1694 {
1695 return 0;
1696 }
1697 else
1698 {
1699 return fSelected;
1700 }
1701}
1702
1703Double_t NcJet::GetDistance(NcPosition* p,Float_t scale)
1704{
1725
1726 Double_t dist=-1.;
1727 fDresult=0.;
1728
1729 if (!p) return dist;
1730
1731 // Obtain a defined position on this jet
1732 NcPosition* rx=fRef;
1733
1734 if (!rx) return dist;
1735
1737
1738 if (pj.GetNorm() <= 0.) return dist;
1739
1740 NcTrack tj;
1741 tj.Set3Momentum(pj);
1742 tj.SetReferencePoint(*rx);
1743 dist=tj.GetDistance(p,scale);
1745 return dist;
1746}
1747
1748Double_t NcJet::GetDistance(NcTrack* t,Float_t scale)
1749{
1770
1771 Double_t dist=-1.;
1772 fDresult=0.;
1773
1774 if (!t) return dist;
1775
1776 // Obtain a defined position on this jet
1777 NcPosition* rx=fRef;
1778
1779 if (!rx) return dist;
1780
1782
1783 if (pj.GetNorm() <= 0.) return dist;
1784
1785 NcTrack tj;
1786 tj.Set3Momentum(pj);
1787 tj.SetReferencePoint(*rx);
1788 dist=tj.GetDistance(t,scale);
1790 return dist;
1791}
1792
1793Double_t NcJet::GetDistance(NcJet* j,Float_t scale)
1794{
1817
1818 Double_t dist=-1.;
1819 fDresult=0.;
1820
1821 if (!j) return dist;
1822
1823 // Obtain a defined position on jet j
1825
1826 if (!rx) return dist;
1827
1828 Nc3Vector pj=j->Get3Momentum();
1829
1830 if (pj.GetNorm() <= 0.) return dist;
1831
1832 NcTrack tj;
1833 tj.Set3Momentum(pj);
1834 tj.SetReferencePoint(*rx);
1835 dist=GetDistance(tj,scale);
1836 return dist;
1837}
1838
1839Int_t NcJet::GetNsignals(TString classname,Int_t par) const
1840{
1856
1857 if (fNtrk<1) return 0;
1858
1859 TObjArray arr;
1860
1861 Int_t n=0;
1862 NcTrack* tx=0;
1863 Int_t exists=0;
1864 for (Int_t i=1; i<=fNtrk; i++)
1865 {
1866 tx=GetTrack(i);
1867 for (Int_t j=1; j<=tx->GetNsignals(); j++)
1868 {
1869 NcSignal* sx=tx->GetSignal(j);
1870 if (!sx) continue;
1871 exists=0;
1872 for (Int_t k=0; k<arr.GetEntries(); k++)
1873 {
1874 if (sx==(NcSignal*)arr.At(k))
1875 {
1876 exists=1;
1877 break;
1878 }
1879 }
1880 if (!exists)
1881 {
1882 if ((par==0 || par==2) && sx->InheritsFrom(classname.Data()))
1883 {
1884 arr.Add(sx);
1885 continue;
1886 }
1887 NcDevice* parent=sx->GetDevice();
1888 if (!parent) continue;
1889 if ((par==1 || par==2) && parent->InheritsFrom(classname.Data())) arr.Add(sx);
1890 }
1891 }
1892 }
1893 n=arr.GetEntries();
1894 return n;
1895}
1896
1897TObjArray* NcJet::GetSignals(TString classname,Int_t par,TObjArray* signals)
1898{
1924
1925 if (signals)
1926 {
1927 signals->Clear();
1928 }
1929 else
1930 {
1931 if (fSelected)
1932 {
1933 fSelected->Clear();
1934 }
1935 else
1936 {
1937 fSelected=new TObjArray();
1938 }
1939 }
1940
1941 if (fNtrk<1)
1942 {
1943 if (signals)
1944 {
1945 return 0;
1946 }
1947 else
1948 {
1949 return fSelected;
1950 }
1951 }
1952
1953 // Generic array pointer to be used in the logic below
1954 TObjArray* arr=0;
1955 if (signals)
1956 {
1957 arr=signals;
1958 }
1959 else
1960 {
1961 arr=fSelected;
1962 }
1963
1964 NcTrack* tx=0;
1965 Int_t exists=0;
1966 for (Int_t i=1; i<=fNtrk; i++)
1967 {
1968 tx=GetTrack(i);
1969 for (Int_t j=1; j<=tx->GetNsignals(); j++)
1970 {
1971 NcSignal* sx=tx->GetSignal(j);
1972 if (!sx) continue;
1973 exists=0;
1974 for (Int_t k=0; k<arr->GetEntries(); k++)
1975 {
1976 if (sx==(NcSignal*)arr->At(k))
1977 {
1978 exists=1;
1979 break;
1980 }
1981 }
1982 if (!exists)
1983 {
1984 if ((par==0 || par==2) && sx->InheritsFrom(classname.Data()))
1985 {
1986 arr->Add(sx);
1987 continue;
1988 }
1989 NcDevice* parent=sx->GetDevice();
1990 if (!parent) continue;
1991 if ((par==1 || par==2) && parent->InheritsFrom(classname.Data())) arr->Add(sx);
1992 }
1993 }
1994 }
1995
1996 if (signals)
1997 {
1998 return 0;
1999 }
2000 else
2001 {
2002 return fSelected;
2003 }
2004}
2005
2006void NcJet::ShowSignals(TString classname,Int_t par,Int_t mode,TString f,TString u)
2007{
2025
2026 TObjArray hits;
2027 GetSignals(classname.Data(),par,&hits);
2028
2029 Int_t nhits=hits.GetEntries();
2030
2031 cout << " *NcJet::ShowSignals* There are " << nhits
2032 << " signals recorded for (device) class " << classname.Data() << endl;
2033
2034 if (!nhits || !mode) return;
2035
2036 NcPosition r;
2037 for (Int_t i=0; i<nhits; i++)
2038 {
2039 NcSignal* sx=(NcSignal*)hits.At(i);
2040 if (sx) sx->Data(f,u);
2041 if (mode==2)
2042 {
2043 NcDevice* dev=(NcDevice*)sx->GetDevice();
2044 if (!dev) continue;
2045 r=dev->GetPosition();
2046 cout << " Device Position";
2047 r.Data(f,u);
2048 }
2049 }
2050}
2051
2052Double_t NcJet::GetSignalValue(TString classname,TString varname,Int_t mode,Int_t par)
2053{
2070
2071 TObjArray hits;
2072 GetSignals(classname.Data(),par,&hits);
2073
2074 Int_t nhits=hits.GetEntries();
2075
2076 if (!nhits) return 0;
2077
2078 Double_t val=0;
2079 for (Int_t i=0; i<nhits; i++)
2080 {
2081 NcSignal* sx=(NcSignal*)hits.At(i);
2082 if (!sx) continue;
2083 val+=sx->GetSignal(varname,mode);
2084 }
2085
2086 return val;
2087}
2088
2089void NcJet::SetEscale(Float_t scale)
2090{
2102
2103 if (scale>0)
2104 {
2105 fEscale=scale;
2106 }
2107 else
2108 {
2109 cout << " *NcJet::SetEscale* Invalid scale value : " << scale << endl;
2110 }
2111}
2112
2113Float_t NcJet::GetEscale() const
2114{
2124
2125 return fEscale;
2126}
2127
2128TObject* NcJet::Clone(const char* name) const
2129{
2141
2142 NcJet* jet=new NcJet(*this);
2143 if (name)
2144 {
2145 if (strlen(name)) jet->SetName(name);
2146 }
2147 return jet;
2148}
2149
ClassImp(NcJet)
Handling of 3-vectors in various reference frames.
Definition Nc3Vector.h:15
Double_t GetResultError() const
void GetVector(Double_t *v, TString f, TString u="rad") const
Double_t GetNorm()
Double_t GetPseudoRapidity()
Double_t GetScaTrans()
virtual void SetZero()
Double_t Dot(Nc4Vector &q)
Double_t GetScaLong()
Nc3Vector GetVecLong() const
Double_t GetResultError() const
virtual void Data(TString f="car", TString u="rad")
Double_t GetScalar()
Nc3Vector GetVecTrans() const
Nc3Vector Get3Vector() const
Nc3Vector fV
Definition Nc4Vector.h:74
Double32_t fDresult
! The error on the scalar result of an operation (e.g. dotproduct)
Definition Nc4Vector.h:77
Signal (Hit) handling of a generic device.
Definition NcDevice.h:14
Creation and investigation of a jet of particle tracks.
Definition NcJet.h:18
Double_t GetInvmass(Float_t scale=-1)
Definition NcJet.cxx:703
Double_t GetEnergy(Float_t scale=-1)
Definition NcJet.cxx:628
void SetReferencePoint(NcPosition &p)
Definition NcJet.cxx:1481
void ShowTracks(Int_t mode=1, TString f="car", TString u="rad", TObjArray *tracks=0)
Definition NcJet.cxx:1159
Double_t GetEl(Float_t scale=-1)
Definition NcJet.cxx:1326
NcPosition * GetReferencePoint()
Definition NcJet.cxx:1501
virtual void List(TString f="car", TString u="rad", TObjArray *tracks=0)
Definition NcJet.cxx:438
virtual ~NcJet()
Definition NcJet.cxx:177
TObjArray * GetTracks(Int_t idmode=0, Int_t chmode=2, Int_t pcode=0, TObjArray *tracks=0)
Definition NcJet.cxx:792
Int_t fNtrk
Definition NcJet.h:83
Double_t GetPl(Float_t scale=-1)
Definition NcJet.cxx:1267
Int_t fNtinit
Definition NcJet.h:80
Int_t GetNtracks(Int_t idmode=0, Int_t chmode=2, Int_t pcode=0)
Definition NcJet.cxx:564
void RemoveTracks(Int_t idmode=0, Int_t chmode=2, Int_t pcode=0)
Definition NcJet.cxx:1057
Int_t GetTrackCopy() const
Definition NcJet.cxx:1446
NcTrack * GetTrack(Int_t i) const
Definition NcJet.cxx:751
Double_t GetDistance(NcPosition *p, Float_t scale=-1)
Definition NcJet.cxx:1703
Double_t GetRapidity()
Definition NcJet.cxx:1385
NcTrack * GetIdTrack(Int_t id) const
Definition NcJet.cxx:773
void RemoveTrack(NcTrack *t)
Definition NcJet.cxx:1005
TObjArray * fSelected
! Temp. array to hold user selected or ordered objects
Definition NcJet.h:88
NcPositionObj * fRef
Definition NcJet.h:87
Double_t GetMomentum(Float_t scale=-1)
Definition NcJet.cxx:659
Int_t GetId() const
Definition NcJet.cxx:1470
Double_t GetSignalValue(TString classname, TString varname, Int_t mode=0, Int_t par=2)
Definition NcJet.cxx:2052
NcJet()
Definition NcJet.cxx:121
Double_t GetPt(Float_t scale=-1)
Definition NcJet.cxx:1239
void Init()
Definition NcJet.cxx:136
Int_t fNtmax
Definition NcJet.h:81
Int_t fTrackCopy
Definition NcJet.h:85
Double_t GetEt(Float_t scale=-1)
Definition NcJet.cxx:1301
void SetEscale(Float_t scale)
Definition NcJet.cxx:2089
Int_t fUserId
Definition NcJet.h:86
TObjArray * GetSignals(TString classname, Int_t par=0, TObjArray *signals=0)
Definition NcJet.cxx:1897
virtual void Reset()
Definition NcJet.cxx:296
Float_t fQ
Definition NcJet.h:82
void ShowSignals(TString classname, Int_t par=0, Int_t mode=1, TString f="car", TString u="rad")
Definition NcJet.cxx:2006
virtual void SetOwner(Bool_t own=kTRUE)
Definition NcJet.cxx:202
TObjArray * fTracks
Definition NcJet.h:84
virtual TObject * Clone(const char *name="") const
Definition NcJet.cxx:2128
TObjArray * SortTracks(Int_t mode=-1, TObjArray *tracks=0, TObjArray *ordered=0)
Definition NcJet.cxx:1520
void SetId(Int_t id)
Definition NcJet.cxx:1459
Float_t GetCharge() const
Definition NcJet.cxx:740
virtual void Data(TString f="car", TString u="rad")
Definition NcJet.cxx:408
void SetNtinit(Int_t n=2)
Definition NcJet.cxx:273
void SetTrackCopy(Int_t j)
Definition NcJet.cxx:1413
Float_t fEscale
Definition NcJet.h:89
Double_t GetMt(Float_t scale=-1)
Definition NcJet.cxx:1352
virtual void ListAll(TString f="car", TString u="rad", TObjArray *tracks=0)
Definition NcJet.cxx:501
Int_t GetNsignals(TString classname="TObject", Int_t par=0) const
Definition NcJet.cxx:1839
void ReplaceTrack(NcTrack *told, NcTrack *tnew)
Definition NcJet.cxx:1086
Float_t GetEscale() const
Definition NcJet.cxx:2113
Nc3Vector Get3Momentum(Float_t scale=-1) const
Definition NcJet.cxx:684
void AddTrack(NcTrack &t)
Definition NcJet.cxx:313
Handling of positions (with timestamps) in various reference frames.
Definition NcPosition.h:18
void GetPosition(Double_t *r, TString f, TString u="rad", Float_t s=-1) const
virtual void Data(TString f="car", TString u="rad") const
Handling of positions in various reference frames.
Generic handling of (extrapolated) detector signals.
Definition NcSignal.h:23
NcDevice * GetDevice() const
virtual Float_t GetSignal(Int_t j=1, Int_t mode=0) const
Definition NcSignal.cxx:651
virtual void Data(TString f="car", TString u="rad") const
Definition NcSignal.cxx:959
Handling of the attributes of a reconstructed particle track.
Definition NcTrack.h:19
Double_t GetEt(Float_t scale=-1)
Definition NcTrack.cxx:1586
virtual void ListAll(TString f="car", TString u="rad")
Definition NcTrack.cxx:541
Double_t GetRapidity()
Definition NcTrack.cxx:1670
void SetReferencePoint(NcPosition &p)
Definition NcTrack.cxx:1469
virtual TObject * Clone(const char *name="") const
Definition NcTrack.cxx:2167
Int_t GetId() const
Definition NcTrack.cxx:1783
Double_t GetPl(Float_t scale=-1)
Definition NcTrack.cxx:1552
NcSignal * GetSignal(Int_t j) const
Definition NcTrack.cxx:1017
Int_t GetParticleCode() const
Definition NcTrack.cxx:1867
Int_t GetNsignals() const
Definition NcTrack.cxx:970
Double_t GetDistance(NcPosition *p, Float_t scale=-1)
Definition NcTrack.cxx:2011
Double_t GetEl(Float_t scale=-1)
Definition NcTrack.cxx:1611
virtual void Data(TString f="car", TString u="rad")
Definition NcTrack.cxx:455
void Set3Momentum(Nc3Vector &p)
Definition NcTrack.cxx:401
Float_t GetProb() const
Definition NcTrack.cxx:1911
Double_t GetEnergy(Float_t scale=-1)
Definition NcTrack.cxx:721
Double_t GetMass(Float_t scale=-1)
Definition NcTrack.cxx:671
Double_t GetMt(Float_t scale=-1)
Definition NcTrack.cxx:1637
Float_t GetEscale() const
Definition NcTrack.cxx:1841
Float_t GetCharge() const
Definition NcTrack.cxx:710
Double_t GetPt(Float_t scale=-1)
Definition NcTrack.cxx:1524
Double_t GetMomentum(Float_t scale=-1)
Definition NcTrack.cxx:627