FEDRA emulsion software from the OPERA Collaboration
EdbDSRec Class Reference

#include <EdbDataStore.h>

Inheritance diagram for EdbDSRec:
Collaboration diagram for EdbDSRec:

Public Member Functions

void Clear (bool hard=false)
 
int DoDecaySearch ()
 
int DoFindBlkSeg (EdbVertex *v, int w0, double ImpMax=50., double RMax=3000, int Dpat=1)
 TODO. More...
 
int DoMomEst ()
 
int DoTracking (bool use_btk=true, int p0=0, int p1=100)
 
int DoTracking0 (bool use_btk=true, int p0=0, int p1=100)
 
int DoVertexing ()
 
 EdbDSRec ()
 
void FillECovPV (EdbPatternsVolume *, EdbScanCond *cnd=0)
 prepare segments' cov matrix More...
 
void FillECovSeg (EdbSegP *seg, EdbScanCond *cnd=0)
 
void FillECovTrks ()
 
void FillErrorsCOV ()
 
 ~EdbDSRec ()
 
- Public Member Functions inherited from EdbDataStore
void AddPattern (EdbPattern *pat)
 
EdbSegPAddSegment (EdbSegP *seg, EdbSegmentCut *cut=0, int Plt0=0, int Plt1=100)
 add methods More...
 
void AddTrack (EdbTrackP *tr)
 
void AddVertex (EdbVertex *v)
 
void Clear (bool hard=false)
 clear methods More...
 
void ClearGeom ()
 
void ClearRaw (bool hard=false)
 
void ClearSeg (bool hard=false)
 
void ClearTracks (bool hard=false)
 
void ClearVTX ()
 
void DoEfficiency (TF1 *eff_seg, TF1 *eff_mtk)
 
void DoSmearing (EdbScanCond *cond_btk, EdbScanCond *cond_mtk=0)
 methods for simulation: More...
 
 EdbDataStore ()
 
EdbLayerFindLayer (int plate, int side=0)
 
EdbTrackPFindLongTrk (int nsmin=8)
 
EdbPatternFindPattern (int plate, int side=0)
 
EdbVertexFindPrimVtx ()
 
EdbTrackPFindTrack (int id)
 find methods More...
 
EdbVertexFindVertex (int id)
 
long Gen_mtk_BG (long NBG, int Plate, int Side, TH2 *pdf_Ang, TH2 *pdf_WT=0)
 
EdbPatternGetPattern (int n, bool btk=true)
 
EdbPatternsVolumeGetPV (bool btk=true)
 
EdbPatternGetRawPat (int n)
 
EdbPatternGetSegPat (int n)
 
EdbTrackPGetTrack (int n)
 get methods More...
 
EdbVertexGetVertex (int n)
 
void LoadMCVertices (TObjArray *vtx)
 restore MC info methods More...
 
void MakePattern (double z, int plate, int side)
 
int Nplt ()
 
int Nt ()
 count methods: More...
 
int Nv ()
 
void PrintBrief ()
 print methods More...
 
void PrintPatterns ()
 
void PrintTracks (int vlev=0)
 
void Restore_PatFromGeom (int np0=0, int np1=1000)
 
void Restore_PIDFromID ()
 
void Restore_SegFromTrx (EdbSegmentCut *cut=0, int Plt0=0, int Plt1=1000)
 
void Restore_TrxFromVtx ()
 
void SavePlateToRaw (const char *fname, int PID, Option_t *option="RECREATE")
 
void SaveToRaw (const char *dir="./", const EdbID &idset="0.0.0.0", Option_t *option="RECREATE", bool doaff=true)
 save methods: More...
 
void SetOwnTracks (bool own=true)
 setown methods: More...
 
void SetOwnTrkSegs ()
 
void SetOwnVertices (bool own=true)
 
void TransferGeometry (EdbDataStore *ds)
 
void TransferTo (EdbDataStore *ds, char level, EdbSegmentCut *cut=0, int FromPlate=0, int ToPlate=57)
 transfer methods More...
 
 ~EdbDataStore ()
 

Public Attributes

EdbScanCond eCond_b
 
EdbScanCond eCond_m
 
EdbMomentumEstimator eMomEst
 
EdbVertexRec eVRec
 
- Public Attributes inherited from EdbDataStore
EdbBrickPeBrick
 
EdbPatternsVolume eRawPV
 geometry More...
 
EdbPatternsVolume eSegPV
 
TObjArray eTracks
 
TObjArray eVTX
 

Additional Inherited Members

- Static Public Member Functions inherited from EdbDataStore
static void TransferSegs (EdbPatternsVolume *pv0, EdbPatternsVolume *pv1, EdbSegmentCut *cut=0, int FromPlate=0, int ToPlate=57)
 

Constructor & Destructor Documentation

◆ EdbDSRec()

EdbDSRec::EdbDSRec ( )

init Vertexing:

init momentum estimator:

466 {
469 eVRec.eDZmax=5000;
470 eVRec.eProbMin=0.001;
471 eVRec.eImpMax=50;
472 eVRec.eUseMom=0;
477 eMomEst.eAlg = 0;
478 eMomEst.eMinEntr = 3;
479}
EdbVertexRec eVRec
Definition: EdbDataStore.h:110
EdbMomentumEstimator eMomEst
Definition: EdbDataStore.h:111
TObjArray eTracks
Definition: EdbDataStore.h:84
int eAlg
select the algorithm for PMS estimation
Definition: EdbMomentumEstimator.h:25
void SetParPMS_Mag()
Definition: EdbMomentumEstimator.cxx:73
int eMinEntr
min number of entries in the cell to accept it for fitting (def=1)
Definition: EdbMomentumEstimator.h:28
Bool_t eUseMom
use or not track momentum for vertex calculations
Definition: EdbVertex.h:181
Int_t eQualityMode
vertex quality estimation method (0:=Prob/(sigVX^2+sigVY^2); 1:= inverse average track-vertex distanc...
Definition: EdbVertex.h:183
Bool_t eUseSegPar
use only the nearest measured segments for vertex fit (as Neuchatel)
Definition: EdbVertex.h:182
Float_t eImpMax
maximal acceptable impact parameter (preliminary check)
Definition: EdbVertex.h:179
Float_t eDZmax
maximum z-gap in the track-vertex group
Definition: EdbVertex.h:177
Float_t eProbMin
minimum acceptable probability for chi2-distance between tracks
Definition: EdbVertex.h:178
TObjArray * eEdbTracks
Definition: EdbVertex.h:204

◆ ~EdbDSRec()

EdbDSRec::~EdbDSRec ( )
inline
93{};

Member Function Documentation

◆ Clear()

void EdbDSRec::Clear ( bool  hard = false)

482 {
484 eVRec.Reset();
485}
void Clear(bool hard=false)
clear methods
Definition: EdbDataStore.h:30
void Reset()
Definition: EdbVertex.cxx:722

◆ DoDecaySearch()

int EdbDSRec::DoDecaySearch ( )

716 {
717 return 0;
718}

◆ DoFindBlkSeg()

int EdbDSRec::DoFindBlkSeg ( EdbVertex v,
int  w0,
double  ImpMax = 50.,
double  RMax = 3000,
int  Dpat = 1 
)

TODO.


segment assigned to track

not black enough

681 {
682 int nn=0;
684 EdbVTA* vta=0;
685 printf("start\n");
686 EdbPattern* pat=0;
687 EdbSegP* seg=0;
688 double dzmax=Dpat*1300;
689 double dx,dy,dz;
690 double imp=0;
691 for(int np=0;np<eSegPV.Npatterns();++np){
692 pat=eSegPV.GetPattern(np);
693 dz=fabs(pat->Z()-v->Z());
694 if(dz<0 || dz>dzmax)continue;
695 for(int ns=0;ns<pat->N();++ns){
696 seg=pat->GetSegment(ns);
697 if(seg->Track()!=-1)continue;
698 if(seg->W()<w0)continue;
699 imp=v->DistSeg(seg);
700 printf("imp =%2.1f\n",imp);
701 if (imp > ImpMax) continue;
702 dx=seg->X()-v->X();
703 dy=seg->Y()-v->Y();
704 vta = new EdbVTA((EdbTrackP *)seg, v);
705 vta->SetZpos(1);
706 vta->SetFlag(1);
707 vta->SetImp(imp);
708 vta->SetDist(sqrt(dx*dx+dy*dy+dz*dz));
709 v->AddVTA(vta);
710 ++nn;
711 }
712 }
713 return nn;
714}
brick dz
Definition: RecDispMC.C:107
Expr< UnaryOp< Fabs< T >, Expr< A, T, D >, T >, T, D > fabs(const Expr< A, T, D > &rhs)
Definition: UnaryOperators.hh:96
EdbPatternsVolume eSegPV
Definition: EdbDataStore.h:83
Definition: EdbPattern.h:273
Int_t Npatterns() const
Definition: EdbPattern.h:366
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
Definition: EdbSegP.h:21
Int_t Track() const
Definition: EdbSegP.h:150
Float_t X() const
Definition: EdbSegP.h:173
Float_t Y() const
Definition: EdbSegP.h:174
Float_t W() const
Definition: EdbSegP.h:151
Float_t Z() const
Definition: EdbPattern.h:84
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
Definition: EdbPattern.h:113
Definition: EdbVertex.h:26
void SetImp(float imp)
Definition: EdbVertex.h:57
void SetDist(float dist)
Definition: EdbVertex.h:58
void SetZpos(int zpos)
Definition: EdbVertex.h:55
void SetFlag(int flag)
Definition: EdbVertex.h:56
Float_t DistSeg(EdbSegP *seg, float X0=0.)
Definition: EdbVertex.cxx:672
void AddVTA(EdbVTA *vta)
Definition: EdbVertex.cxx:355
Float_t X() const
Definition: EdbVertex.h:130
void ClearNeighborhood()
Definition: EdbVertex.cxx:275
Float_t Z() const
Definition: EdbVertex.h:132
Float_t Y() const
Definition: EdbVertex.h:131
float ImpMax
Definition: check_vertex.C:30

◆ DoMomEst()

int EdbDSRec::DoMomEst ( )

647 {
648
649 float p,p0,p1;
650 float tl,tt;
651 EdbTrackP* trk=0;
652 EdbTrackP* clon=0;
653 int Nest=0;
654 for(int nt=0;nt<eTracks.GetEntries();++nt){
655 trk=GetTrack(nt);
656 clon=(EdbTrackP*)trk->Clone();
657 Log(1,"EdbDSRec::DoMomEst",Form("Mom Est: track #%d (%dseg)",nt,trk->N()));
658// if(trk->N()<5){trk->SetP(0); printf("track #%d (%dseg) - skip!\n",nt,trk->N());continue;}
659 p=eMomEst.PMSang(*clon);
660// if(eMomEst.eStatus==-1){printf("track #%d (%dseg) - skip!\n",nt,trk->N());continue;}
661 printf("track #%d (%dseg) mom %f - status %d\n",nt,trk->N(),p,eMomEst.eStatus);
662 p0=eMomEst.ePmin;
663 p1=eMomEst.ePmax;
664/* tl=eMomEst.eGX->GetY()[0];
665 tt=eMomEst.eGY->GetY()[0]; */
666// printf("track #%d (%dseg) mom %f [%2.2f <> %2.2f]\n",nt,trk->N(),p,p0,p1);
667// printf("dTL=%2.4g, dTT=%2.4g\n",tl,tt);
668// if(p<0 && (tl>0.015 || tt>0.015)){
669// printf("BAD!\n");
670// trk->SetP(-p);
671// continue;
672// }
673// if(p<0)p=0;
674 trk->SetP(p);
675 trk->SetPerr(p0,p1);
676 delete clon;
677 }
678 return Nest;
679}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
EdbTrackP * GetTrack(int n)
get methods
Definition: EdbDataStore.h:45
float ePmax
momentum 90% errors range
Definition: EdbMomentumEstimator.h:44
float PMSang(EdbTrackP &tr)
Definition: EdbMomentumEstimator.cxx:255
int eStatus
status of the estimation (-1 nothing was done)
Definition: EdbMomentumEstimator.h:26
float ePmin
Definition: EdbMomentumEstimator.h:44
void SetP(float p)
Definition: EdbSegP.h:133
Int_t N() const
Definition: EdbPattern.h:177
void SetPerr(Float_t perrDown, Float_t perrUp)
Definition: EdbPattern.cxx:1230
p
Definition: testBGReduction_AllMethods.C:8

◆ DoTracking()

int EdbDSRec::DoTracking ( bool  use_btk = true,
int  p0 = 0,
int  p1 = 100 
)

init tracker:

doing tracking using microtracks or basetracks

488 {
490 EdbTrackAssembler tracker;
491 tracker.eCond=use_btk?eCond_b:eCond_m;
492// .SetSigma0(3,3,0.005,0.005);
493 tracker.eDTmax=0.07;
494 tracker.eDRmax=45.;
495// tracker.eDZGapMax=10000;
496 tracker.InitTrZMap( 2400, eBrick->Xmin(), eBrick->Xmax(), 2000, eBrick->Ymin(), eBrick->Ymax(), 30 );
498 int nsegmin=use_btk?2:3;
499 SetOwnTracks(1);
500 EdbPatternsVolume* pv=use_btk?(&eSegPV):(&eRawPV);
501 int npl=pv->Npatterns();
502 EdbPattern *pat=0;
503 EdbLayer *plate=0;
504 for(int ipass=0; ipass<2; ipass++) {
505 printf("\n\n*************** ipass=%d ************\n",ipass);
506 tracker.eCollisionsRate=0;
507 for(int i=0; i<npl; i++) {
508 if(i<p0 || i>p1)continue;
509 pat=pv->GetPattern(i);
510 if(pat->N()==0)continue;
511 Log(3,Form("EdbDSRec::DoTracking(%d,%d,%d)",use_btk,p0,p1),Form("plate=%d at Z=%4.1f side=%d [%d segm]\n",pat->Plate(),pat->Z(),pat->Side(),pat->N()));
512 plate = eBrick->GetPlate(pat->Plate());
513// printf("pattern z = (%2.4f vs %2.4f)\n", pat->Z(),plate->Z());
514 pat->SetSegmentsZ();
515 pat->SetSegmentsPID();
516 pat->Transform( plate->GetAffineXY() );
517 pat->TransformShr( plate->Shr() );
518 pat->TransformA( plate->GetAffineTXTY() );
519 pat->SetSegmentsPlate(pat->Plate());
520// printf("flag0=%d\n",pat->GetSegment(0)->Flag());
521
522 // pat=ds_rec->eRawPV.GetPattern(i);
523
524// printf("pat #%d has %d segments\n",i,pat->N());
525 if(i>0) tracker.ExtrapolateTracksToZ(pat->Z());
526 tracker.FillTrZMap();
527 tracker.AddPattern(*pat);
528// printf("ok\n");
529 }
530 }
531
532 int Ntr=0;
533 int ntr = tracker.Tracks().GetEntries();
534 for( int i=0; i<ntr; i++ ) {
535 EdbTrackP *t = (EdbTrackP*)(tracker.Tracks().At(i));
536 int nseg=t->N();
537// printf("trk %d/%d: nseg=%d (min %d)\n",i,ntr,nseg);
538 if(nseg>=nsegmin) {
539 for(int j=0; j<nseg; j++) {
540 EdbSegP *s = t->GetSegment(j);
541 s->SetErrors();
542 tracker.eCond.FillErrorsCov(s->TX(),s->TY(),s->COV());
543 }
544 t->SetP(1.);
545 t->FitTrackKFS(0);
546 //fit.FitTrackLine(*t);
547 tracker.RecalculateSegmentsProb(*t);
548 t->SetCounters();
549 AddTrack(new EdbTrackP(*t));
550
551 Ntr++;
552 }
553// else
554 //if(t->N()>2) t->PrintNice();
555 }
556 return Ntr;
557}
int Ntr
Definition: RecDispNU.C:88
EdbPlateP * GetPlate(int i)
Definition: EdbBrick.h:53
EdbScanCond eCond_m
Definition: EdbDataStore.h:112
EdbScanCond eCond_b
Definition: EdbDataStore.h:112
void SetOwnTracks(bool own=true)
setown methods:
Definition: EdbDataStore.h:41
EdbPatternsVolume eRawPV
geometry
Definition: EdbDataStore.h:82
void AddTrack(EdbTrackP *tr)
Definition: EdbDataStore.h:62
EdbBrickP * eBrick
Definition: EdbDataStore.h:81
Definition: EdbLayer.h:39
float Ymin() const
Definition: EdbLayer.h:87
float Ymax() const
Definition: EdbLayer.h:88
float Xmax() const
Definition: EdbLayer.h:86
float Xmin() const
Definition: EdbLayer.h:85
void SetSegmentsPID()
Definition: EdbPattern.cxx:1455
Int_t Side() const
Definition: EdbPattern.h:328
Int_t Plate() const
Definition: EdbPattern.h:327
Definition: EdbPattern.h:334
virtual void Transform(const EdbAffine2D *a)
Definition: EdbVirtual.cxx:155
void FillErrorsCov(float tx, float ty, TMatrixD &cov)
Definition: EdbScanCond.cxx:161
void TransformA(const EdbAffine2D *affA)
Definition: EdbPattern.cxx:324
void TransformShr(const float shr)
Definition: EdbPattern.cxx:341
void SetSegmentsZ()
Definition: EdbPattern.cxx:250
void SetSegmentsPlate(int plate)
Definition: EdbPattern.cxx:241
generic class for the tracks assembling from segments
Definition: EdbScanTracking.h:17
TObjArray & Tracks()
Definition: EdbScanTracking.h:75
void FillTrZMap()
Definition: EdbScanTracking.cxx:290
EdbScanCond eCond
Definition: EdbScanTracking.h:40
void RecalculateSegmentsProb(EdbTrackP &t)
Definition: EdbScanTracking.cxx:181
void ExtrapolateTracksToZ(float z, int nsegmin=0)
Definition: EdbScanTracking.cxx:259
void InitTrZMap(const char *str)
Definition: EdbScanTracking.cxx:301
void AddPattern(EdbPattern &p)
Definition: EdbScanTracking.cxx:105
float eDRmax
position acceptance for the fast preselection
Definition: EdbScanTracking.h:34
float eDTmax
angular acceptance for the fast preselection
Definition: EdbScanTracking.h:33
int eCollisionsRate
Definition: EdbScanTracking.h:39
TTree * t
Definition: check_shower.C:4
s
Definition: check_shower.C:55
Int_t plate
Definition: merge_Energy_SytematicSources_Electron.C:1
int nsegmin
Definition: check_vertex.C:23

◆ DoTracking0()

int EdbDSRec::DoTracking0 ( bool  use_btk = true,
int  p0 = 0,
int  p1 = 100 
)

transfer patterns

propagate

add track here

561 {
562 float momentum = 0.3;
563 float mass = 0.139; // particle mass
564 float ProbMinP = 0.001; // minimal probability to accept segment on propagation
565 int nsegmin = 1; // minimal number of segments to propagate this track
566 int ngapmax = 3; // maximal gap for propagation
567
569 EdbPatternsVolume* pv=use_btk?&eSegPV:&eRawPV;
570
572 EdbPattern* pat,*pat1;
573 for(int np=0; np<pv->Npatterns(); ++np){
574 pat=pv->GetPattern(np);
575 if(pat->Plate()<p0 || pat->Plate()>p1)continue;
576 pat1=(EdbPattern*)pat->Clone();
577 pat1->SetSegmentsZ();
578 pat1->SetID(np);
579 pat1->SetSegmentsPID();
580 Log(1,"EdbDSRec",Form("add pattern %d with %d segments",pat1->PID(),pat1->N()));
581 gAli.AddPattern(pat1);
582 }
583
584 gAli.SetScanCond(use_btk?&eCond_b:&eCond_m);
588 Log(2,"EdbDSRec::DoTracking0",Form("ntr = %d\n",ntr));
589
592
593 EdbTrackP *tr=0;
594 ntr = gAli.eTracks->GetEntries();
595
596 for(int i=0; i<ntr; i++) {
597 tr = (EdbTrackP*)(gAli.eTracks->At(i));
598 tr->SetID(i);
599 tr->SetSegmentsTrack();
600 tr->SetFlag(0);
601 }
602
603 int nadd = 0;
604 int nseg=0;
605 int Ntr=0;
606// nadd=gAli.PropagateTracks(0,57,0.001,3,0);
607 for(int i=0; i<ntr; i++) {
608 tr = (EdbTrackP*)(gAli.eTracks->At(i));
609
610 float p=momentum;
611
612 tr->SetErrorP(0.2*0.2*p*p);
613
614 nseg = tr->N();
615 tr->SetP(p);
616 printf("pid=%d\n",tr->GetSegmentFirst()->PID());
617 if(tr->Flag()<0) continue;
618 // tr->PrintNice();
619 nadd += gAli.PropagateTrack( *tr, true, 0.001, 3, 0 );
620 if(tr->Flag()<0) printf("%d flag=%d\n",i,tr->Flag());
621 //if(tr->N() != nseg) printf("%d nseg=%d (%d) \t p = %f\n",i,tr->N(),nseg,tr->P());
623 AddTrack((EdbTrackP*)tr->Clone());
624 Ntr++;
625 }
626 printf("nadd = %d\n",nadd);
627
628 return Ntr;
629}
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
static int LinkTracksWithFlag(EdbPVRec *ali, float p, float probmin, int nsegmin, int maxgap, int flag, float mass=0.1396)
Definition: EdbDataSet.cxx:2256
Definition: EdbPVRec.h:148
void SetCouplesPeriodic(int istart, int iperiod)
Definition: EdbPVRec.cxx:982
TObjArray * eTracks
Definition: EdbPVRec.h:161
int PropagateTrack(EdbTrackP &tr, bool followZ, float probMin=0.05, int ngapMax=3, int design=0)
Definition: EdbPVRec.cxx:2569
void SetScanCond(EdbScanCond *scan)
Definition: EdbPVRec.h:171
void FitTracks(float p=10., float mass=0.139, TObjArray *gener=0, int design=0)
Definition: EdbPVRec.cxx:1893
void SetID(int id)
Definition: EdbPattern.h:309
int PID() const
Definition: EdbPattern.h:320
void AddPattern(EdbPattern *pat)
Definition: EdbPattern.cxx:1707
Definition: EdbTrackFitter.h:17
EdbPVRec * gAli
Definition: check_vertex.C:14
float ProbMinP
Definition: check_vertex.C:22
float mass
Definition: check_vertex.C:21
float momentum
Definition: check_vertex.C:20
int ngapmax
Definition: check_vertex.C:24

◆ DoVertexing()

int EdbDSRec::DoVertexing ( )

631 {
633 // performing vertexing
634 if(gEDBDEBUGLEVEL>=1) printf("%d tracks for vertexing\n", eVRec.eEdbTracks->GetEntries() );
635 int nvtx = eVRec.FindVertex();
636 if(gEDBDEBUGLEVEL>=1) printf("%d 2-track vertexes was found\n",nvtx);
637 if(nvtx == 0) return 0;
639 for(int nv=0;nv<eVRec.Nvtx();nv++){
641 if(v->Flag()<0)delete v;
642 else eVTX.Add(v);
643 }
644 return eVRec.Nvtx();
645}
TObjArray eVTX
Definition: EdbDataStore.h:85
void SetOwnVertices(bool own=true)
Definition: EdbDataStore.h:42
Int_t ProbVertexN()
Definition: EdbVertex.cxx:1426
Int_t Nvtx() const
Definition: EdbVertex.h:287
Int_t FindVertex()
Definition: EdbVertex.cxx:1065
EdbVertex * GetVertex(Int_t &i)
Definition: EdbVertex.h:288
Definition: EdbVertex.h:69
Int_t Flag() const
Definition: EdbVertex.h:124
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ FillECovPV()

void EdbDSRec::FillECovPV ( EdbPatternsVolume pv,
EdbScanCond cnd = 0 
)

prepare segments' cov matrix

740 {
741 assert(pv);
742 EdbPattern* pat=0;
743 EdbSegP* seg;
744 for(int i=0;i<pv->Npatterns();i++){
745 pat=pv->GetPattern(i);
746 pat->SetSegmentsPID();
747 for(int j=0;j<pat->GetN();j++){
748 seg=pat->GetSegment(j);
749 seg->SetID(j);
750 FillECovSeg(seg);
751 }
752 }
753}
void FillECovSeg(EdbSegP *seg, EdbScanCond *cnd=0)
Definition: EdbDataStore.cxx:731
void SetID(int id)
Definition: EdbSegP.h:128
Int_t GetN() const
Definition: EdbPattern.h:65

◆ FillECovSeg()

void EdbDSRec::FillECovSeg ( EdbSegP seg,
EdbScanCond cnd = 0 
)

if cnd=0 - decide which to use

731 {
732 assert(seg);
733 seg->SetErrors();
734 seg->SetErrorP(1.);
735 if(cnd==0)
736 cnd=(seg->DZ()>100)?&eCond_b:&eCond_m;
737 cnd->FillErrorsCov(seg->TX(),seg->TY(),seg->COV());
738}
void SetErrors()
Definition: EdbSegP.h:90
Float_t DZ() const
Definition: EdbSegP.h:154
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
TMatrixD & COV() const
Definition: EdbSegP.h:123
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
void SetErrorP(float sp2)
Definition: EdbSegP.h:94

◆ FillECovTrks()

void EdbDSRec::FillECovTrks ( )

721 {
722 EdbTrackP* trk=0;
723 TObjArray trx1;
724 for(int nt=0;nt<Nt();nt++){
725 trk=GetTrack(nt);
726 for(int ns=0;ns<trk->N();ns++) FillECovSeg(trk->GetSegment(ns));
727 trk->FitTrackKFS(0);
728 }
729}
int Nt()
count methods:
Definition: EdbDataStore.h:37
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:195
int FitTrackKFS(bool zmax=false, float X0=5810., int design=0)
Definition: EdbPattern.cxx:722

◆ FillErrorsCOV()

void EdbDSRec::FillErrorsCOV ( )

755 {
758 FillECovTrks();
759}
void FillECovPV(EdbPatternsVolume *, EdbScanCond *cnd=0)
prepare segments' cov matrix
Definition: EdbDataStore.cxx:740
void FillECovTrks()
Definition: EdbDataStore.cxx:721

Member Data Documentation

◆ eCond_b

EdbScanCond EdbDSRec::eCond_b

◆ eCond_m

EdbScanCond EdbDSRec::eCond_m

◆ eMomEst

EdbMomentumEstimator EdbDSRec::eMomEst

◆ eVRec

EdbVertexRec EdbDSRec::eVRec

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