BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//_____________________________________________________________________
//
//  2 dimensional track finder module
//
//_____________________________________________________________________
//
//  $Id: BrDC2DTrackFinder.cxx,v 1.9 2002/02/20 13:36:46 ufstasze Exp $
//  $Author: ufstasze $
//  $Date: 2002/02/20 13:36:46 $
//
//

#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef ROOT_TString
#include "TString.h"
#endif

#include "BrDC2DTrackFinder.h"

#include "BrGeometryDbManager.h"
#include "BrParameterDbManager.h"

#ifndef BRAT_BrDCHitPosition
#include "BrDCHitPosition.h"
#endif

#include "BrDCHit.h"
#include "BrDetectorHit.h"
#include "BrTrackHit2D.h"
#include "BrTrackHit.h"
#include "BrMath.h"

ClassImp(BrDC2DTrackFinder);

//////////////////////////////////////////////////////////////////////////////////
//
// The BrDC2DTrackFinder class finds 2D tracks from DC's
//
//////////////////////////////////////////////////////////////////////////////////

 BrDC2DTrackFinder::BrDC2DTrackFinder()
{
  // DO NOT USE
  SetState(kSetup);

  fTrackHits   = 0;
  fTrackHits2D = 0;
}

//_____________________________________________________________
 BrDC2DTrackFinder::BrDC2DTrackFinder(const Char_t *Name,const Char_t *Title,
				     BrDCClusterFinder *clusterFinder,
				     const Float_t cut_dx, 
				     const Float_t cut_dy, 
				     const Float_t asigma, 
				     const Int_t max_miss_main, 
				     const Int_t max_miss_sub, 
				     const Int_t max_miss_other, 
				     Bool_t EventByEvent)
  : BrModule(Name,Title)
{

  fClusterFinder = clusterFinder;

  fDCCutMatchDx = cut_dx;
  fDCCutMatchDy = cut_dy;
  fAsigma       = asigma;
  
  fDCCutMaxMissMain  = max_miss_main;
  fDCCutMaxMissSub   = max_miss_sub;
  fDCCutMaxMissOther = max_miss_other;
  fEventByEvent      = EventByEvent;
  
}

//_____________________________________________________________
 BrDC2DTrackFinder::~BrDC2DTrackFinder(){
}

//_____________________________________________________________
 void BrDC2DTrackFinder::Init() 
{
  // create tables of selected hits
  fSelectedDCHits = new BrDataTable("hits", "hits");
  if(!strcmp("T3",GetName()))fNumpv=3;
  else if (!strcmp("T4",GetName()) || !strcmp("T5",GetName()))fNumpv=3;
  else cout<<"BrDC2DTrackFinder::Init(): check your BrDCTrackingModule constructor"<<endl;

  // specific options for matching hits to track.
  // If fCheckSingleHits is TRUE we look for single
  // (one plane) hit but only in the case when we coudn't
  // find the combined hit.
  // It was found that for T3, this option significantly
  // slows down the program execution.

  if(!strcmp("T4",GetName()))fCheckSingleHits = kTRUE;
  if(!strcmp("T5",GetName()))fCheckSingleHits = kTRUE;

  if(fBfsAngle>5 || fRunNo>6166){
    if(!strcmp("T3",GetName()))fCheckSingleHits = kTRUE;
    if(!strcmp("T3",GetName()))fStrictParam=1.0;
    if(!strcmp("T4",GetName()))fStrictParam=0.7;
    if(!strcmp("T5",GetName()))fStrictParam=0.7;
  }else{
    if(!strcmp("T3",GetName()))fCheckSingleHits = kTRUE;
    if(!strcmp("T3",GetName()))fStrictParam=2.0;
    if(!strcmp("T4",GetName()))fStrictParam=0.7;
    if(!strcmp("T5",GetName()))fStrictParam=0.7;
  }

  fAsig=4.0;

  // limitation on the track angle
  // The defaults seem to be 14.7 mrad in bend plane
  // i.e 1.0 deg +- .5deg*4 
  // One should look at the actual angles for real found tracks
  // as a diagnostic (at end of finding tracks?)
  //
  if(!strcmp("T3",GetName())){
    fLimitAngle = kTRUE;
    fOffset[0] = -0.0147;
    fOffset[1] = -0.0015;
    fSigma[0]  =  0.0076;
    fSigma[1]  =  0.0057;
  }
  if(!strcmp("T4",GetName())){
    fLimitAngle = kTRUE;
    fOffset[0] = -0.016;
    fOffset[1] =  0.00064;
    fSigma[0]  =  0.012;
    fSigma[1]  =  0.0043;
  }
  if(!strcmp("T5",GetName())){
    fLimitAngle = kTRUE;
    fOffset[0] = -0.0143;
    fOffset[1] =  0.0039;
    fSigma[0]  =  0.0191;
    fSigma[1]  =  0.0046;
  }

}

//_____________________________________________________________
 void BrDC2DTrackFinder::Find2DTracks(BrDCViewList *viewList,
				     TObjArray *trackHits, 
				     TObjArray *trackHits2D) {
  fTrackHits    = trackHits;
  fTrackHits2D  = trackHits2D;

  if(fEventByEvent)HitsOut(viewList); // hit are written down only for event by event mode

  Int_t frontPlane = 1;
  Int_t backPlane  = 3;
  Int_t checkPlane = 2;
  GenerateInitial2DTracks(viewList,frontPlane,backPlane,checkPlane);
  CheckIntermediatePositions(viewList,checkPlane);
  frontPlane = 1;
  backPlane  = 2;
  checkPlane = 3;
  GenerateInitial2DTracks(viewList,frontPlane,backPlane,checkPlane);
  CheckIntermediatePositions(viewList,checkPlane);
  frontPlane = 2;
  backPlane  = 3;
  checkPlane = 1;
  GenerateInitial2DTracks(viewList,frontPlane,backPlane,checkPlane);
  CheckIntermediatePositions(viewList,checkPlane);
  CleanTracks();
}

//_____________________________________________________________
 void BrDC2DTrackFinder::HitsOut(BrDCViewList *viewList) {
  // this function is only for debugging purpose to write
  // down detector combined hits 
  ofstream OUTHit;
  Char_t file[64];
  sprintf(file,"./visu/%sclusters.hit",GetName());
  OUTHit.open(file,ios::out);
  OUTHit.setf(ios::showpoint);
  BrDCHit *hit_p;
  for(Int_t iv=0;iv<3;iv++) {
    Int_t NPlaneViews = viewList->GetNumPlaneViews(iv);
    //cout<<"iv , NPV: "<<iv<<" "<<NPlaneViews<<endl;
    for(Int_t j=0;j<NPlaneViews;j++) {      
      TObjArray *detectorHitsPlane = viewList->GetDetectorHits(iv,j+1); 
      Int_t Nhits = detectorHitsPlane->GetEntries();
      for(Int_t i=0;i<Nhits;i++) {
        hit_p = (BrDCHit*)detectorHitsPlane->At(i);
        Float_t u=hit_p->GetPos()[0];
        Float_t z=hit_p->GetPos()[2];
        Int_t iview = hit_p->GetImod() - (hit_p->GetImod()/10)*10; // decode iview value
        Int_t idet = hit_p->GetImod()/100 - 1;
	Int_t PrimHits=hit_p->GetNumCombined();
	OUTHit<<u<<" "<<z<<" "<<iview<<" "<<idet<<" "<<PrimHits<<endl;
        //cout<<u<<" "<<z<<" "<<iview<<" "<<idet<<" "<<endl;
      }
    }
  }
  OUTHit.close();
}

//_____________________________________________________________
 void BrDC2DTrackFinder::GenerateInitial2DTracks(BrDCViewList
						*viewList,
						Int_t frontPlane,
						Int_t backPlane,
						Int_t checkPlane) 
{
  //Main and submain view 2d tracking
  //=================================
  //
  //Define local 2d tracks by taking hit_p for the first and last
  //planes in main view and sub view

  Int_t i,j;
  Int_t iv;
  Float_t u1,u2,z1,z2,du1,du2,au,dau;
  Float_t tra2d_pos[2];
  BrDCHit *hit1_p,*hit2_p;
  BrTrackHit2D  *tra2d_p;
  for(iv=0;iv<2;iv++) {
     Int_t numViews = viewList->GetNumPlaneViews(iv);
     TObjArray *detectorHitsFirstPlane = viewList->GetDetectorHits(iv,frontPlane);
     //TObjArray *detectorHitsLastPlane  = viewList->GetDetectorHits(iv,numViews);
     TObjArray *detectorHitsLastPlane  = viewList->GetDetectorHits(iv,backPlane);

     Int_t numHitsFirstPlane = detectorHitsFirstPlane->GetEntries();
     Int_t numHitsLastPlane  = detectorHitsLastPlane->GetEntries();
     
     if(fEventByEvent){
       if(iv==0)cout<<"X: First-Last Plane: "<<numHitsFirstPlane<<" "<<numHitsLastPlane<<endl;
       if(iv==1)cout<<"Y: First-Last Plane: "<<numHitsFirstPlane<<" "<<numHitsLastPlane<<endl;
     }

     for(i=0;i<numHitsFirstPlane;i++) {
        hit1_p = (BrDCHit*)detectorHitsFirstPlane->At(i);
        u1  = hit1_p->GetPos()[0];
        z1  = hit1_p->GetPos()[2];
        du1 = hit1_p->GetDpos()[0];
        for(j=0;j<numHitsLastPlane;j++) {
           hit2_p = (BrDCHit*)detectorHitsLastPlane->At(j);
           u2  = hit2_p->GetPos()[0];
           z2  = hit2_p->GetPos()[2];
           du2 = hit2_p->GetDpos()[0];
           au = (u2-u1)/(z2-z1);
           dau = (Float_t)TMath::Sqrt(du1*du1 + du2*du2)/TMath::Abs(z2-z1);

/*
	BrClonesArray &tra2d = *fTrackHits2D;
	Int_t ntra2d = fTrackHits2D->GetLast()+1;
	tra2d_p = new(tra2d[ntra2d]) BrTrackHit2D();
*/

	   if(fLimitAngle)
	     if(TMath::Abs(au-fOffset[iv])>fAsig*fSigma[iv])continue;

           tra2d_p = new BrTrackHit2D();
           fTrackHits2D->Add(tra2d_p);

           tra2d_pos[0] = hit1_p->GetPos()[2];
           tra2d_pos[1] = hit1_p->GetPos()[0];
           tra2d_p->SetPos(tra2d_pos);
           tra2d_p->SetVec(au);
           tra2d_p->SetDpos(du1);
           tra2d_p->SetDvec(dau);
           tra2d_p->SetNhit(2);
           tra2d_p->SetFlg(iv);        //iv = 0 if main view, 1 if sub view.
           tra2d_p->SetStat(checkPlane);
           tra2d_p->AddDetectorHit(hit1_p);
           tra2d_p->AddDetectorHit(hit2_p);
           InsertTrahit(hit1_p,tra2d_p);
           InsertTrahit(hit2_p,tra2d_p);
           }		//loop over j
        }			//loop over i
     }			//loop over iv  
}

//_____________________________________________________________
 void BrDC2DTrackFinder::CheckIntermediatePositions(BrDCViewList *viewList,
						   Int_t checkPlane) 
{
  //Loop over intermediate positions in the local tracks and compare
  //to the initially defined local tracks.
  //This might create additional ghost 2d tracks

#define max_saved 1024
  
  Int_t i,j,k;
  Int_t iv;
  Int_t NumTra2DD;
  Int_t nmmiss,nok;
  Int_t NoSaved,NoSavedPl,it;
  Float_t u1,u3,z1,z3,au,dau,du1,du3,uc,duc;
  
  Float_t dddd[2];
  Int_t maxmiss[3];
  
  BrDataTable SavedId;
  BrTrackHit2D *tra2d_p1;
  BrTrackHit2D *tra2d_p2;
  BrDCHit *hit_p;
  
  maxmiss[0] = fDCCutMaxMissMain;
  maxmiss[1] = fDCCutMaxMissSub;
  maxmiss[2] = fDCCutMaxMissOther;
  
  dddd[0] = fDCCutMatchDx;
  dddd[1] = fDCCutMatchDy;
  
  NumTra2DD = fTrackHits2D->GetLast()+1;
  
  if(NumTra2DD>800){
    maxmiss[0]=0; // be very strict in this case 
    maxmiss[1]=0;
  }
  
  //Loop over 2D track hits found so far
  if(fEventByEvent)cout<<"Number of Initial 2DTracks: "<<NumTra2DD<<endl; 
  for(i=0;i<NumTra2DD;i++) {
    NoSaved = 1;
    tra2d_p1 = (BrTrackHit2D*)fTrackHits2D->At(i);
    SavedId.Clear();
    SavedId.Add(tra2d_p1);
    if(tra2d_p1->GetStat() == checkPlane) {
      tra2d_p1->SetStat(0);
      iv = tra2d_p1->GetFlg();			//x view (iv=0)or y view (iv=1)
      nmmiss = 0;
      u1  = tra2d_p1->GetPos()[1];
      z1  = tra2d_p1->GetPos()[0];
      au  = tra2d_p1->GetVec();
      dau = tra2d_p1->GetDvec();
      du1 = tra2d_p1->GetDpos();
      //cout<<" Initial 2DTrack: "<<GetName()<<" "<<iv<<" "<<au<<" "<<u1-au*z1<<endl;
      Int_t numViews = viewList->GetNumPlaneViews(iv);
      
      //      Loop over intermediate planes of this view (iv)
      //for(mview=2;mview<numViews;mview++) {
      Int_t mview=checkPlane;
      nok = 0;
      NoSavedPl = NoSaved;
      TObjArray *detectorHitsThisView = viewList->GetDetectorHits(iv,mview);
      Int_t numHitsThisView = detectorHitsThisView->GetEntries();
      for(k=0;k<numHitsThisView;k++) {
	hit_p = (BrDCHit*)detectorHitsThisView->At(k);
	u3  = hit_p->GetPos()[0];
	z3  = hit_p->GetPos()[2];
	du3 = hit_p->GetDpos()[0];
	uc = u1 + au*(z3-z1);
	//if(!strcmp("T5",GetName()))
	//cout<<" In (mview,k) Loop: "<<mview<<" "<<k<<" "<<u3<<" "<<uc<<endl;
	//duc = dddd[iv]; // original version
	duc = du3 + du1;
	Float_t lowerLimit = uc-fAsigma*duc;
	Float_t upperLimit = uc+fAsigma*duc;
	//if(!strcmp("T5",GetName()))
	//cout<<" Low, Upp: "<<lowerLimit<<" "<<upperLimit<<endl;
	if(BrMath::BTWN(lowerLimit,u3,upperLimit)) {
	  //Coming here means that we found a hit that is consistent
	  //with this 2D track
	  nok++;
	  if(nok == 1) {
	    if(NoSaved>max_saved) printf("NoSaved = %dn",NoSaved);
	    //It is the first one, so increment nhit 
	    //and add detector hit to the list.
	    for(it=0;it<NoSaved;it++) {
	      BrTrackHit2D *tra2d_pp = (BrTrackHit2D*)SavedId.At(it);
	      tra2d_pp->IncNhit();
	      tra2d_pp->AddDetectorHit(hit_p);
	      InsertTrahit(hit_p,tra2d_pp);
	    }
	  }
	  //else if(nok > 1) {
	  else if(nok == 1) { 
	    //                  More than one hit -- create a ghost local track
	    if(2*NoSaved > max_saved) {
	      cout<<"**** Error : cannot save GHOSTS **"<<endl;
	    }
	    else {
	      for(j=0;j<NoSavedPl;j++) {
		tra2d_p2 = new BrTrackHit2D(tra2d_p1);
		fTrackHits2D->Add(tra2d_p2);
		
		if(NoSaved+j>max_saved) 
		  printf("NoSaved1 = %dn",NoSaved+j);
		SavedId.AddAt(tra2d_p2,NoSaved+j);
		CopyTrahit1((BrTrackHit2D*)SavedId.At(j),tra2d_p2);
	      }
	      
	      //            Restore the current hitcmb picked at this plane
	      //            before inserting into the trahit structure for
	      //            this new ghost track(s)
	      for(j=0;j<NoSavedPl;j++) {
		BrTrackHit2D* tra2d_pp = (BrTrackHit2D*)SavedId.At(NoSaved+j);
		InsertTrahit(hit_p,tra2d_pp);
	      }
	      
	      //            Update total number of tracks which belong to
	      //            this group of ghosts.
	      NoSaved += NoSavedPl;
	    }
	  }
	}
      }   //loop over k
      
      ////////////////////////////////////////////////////////////
      for(Int_t ipv=0;ipv<fNumpv;ipv++){
	if(nok==0 && fCheckSingleHits){
	  //try to find something in single planes
	  Int_t idet=checkPlane-1;
	  Int_t iview = iv;
	  Int_t ipl = iv*3 + ipv; 
	  z3 = (fClusterFinder->GetDCDriver())->GetPlanePosition(idet,ipl);
	  uc = u1 + au*(z3-z1);
	  fClusterFinder->SelectDCHits(idet,ipl,fSelectedDCHits);
	  Int_t numSelectedHits = fSelectedDCHits->GetEntries();
	  for(Int_t k=0;k<numSelectedHits;k++) {
	    BrDCHitPosition *hitpos_p = (BrDCHitPosition*)fSelectedDCHits->At(k);
	    BrDcDigInfo *diginfo_p = hitpos_p->GetDigDCInfo();     
	    if(diginfo_p->GetUsed() == 0) {
	      u3  = hitpos_p->GetUhit();
	      du3 = (fClusterFinder->GetDCParams())->GetPosres();
	      duc = du3/fStrictParam; // to be rather strict (was 2.5)
	      Float_t lowerLimit = uc-fAsigma*duc;
	      Float_t upperLimit = uc+fAsigma*duc;
	      if(BrMath::BTWN(lowerLimit,u3,upperLimit)) {
		//Coming here means that we found a hit that is consistent
		//with this 2D track, so first create ne BrDCHit:
		BrClonesArray &hit = *fClusterFinder->GetDetectorHits();
		Int_t nhit = (fClusterFinder->GetDetectorHits())->GetLast()+1;
		hit_p = new(hit[nhit]) BrDCHit();
		
		Float_t hit_pos[3];
		Float_t hit_dpos[3];
		hit_pos[0]  = u3;
		hit_pos[1]  = (Float_t)0.0;
		hit_pos[2]  = z3;
		
		hit_p->SetPos(hit_pos);
		hit_dpos[0] = du3;
		hit_dpos[1] = (Float_t)0.0;
		hit_p->SetDpos(hit_dpos);

		//cout<<"2D: "<<idet<<" "<<ipl<<" "<<iview<<endl;

		hit_p->SetImod((idet + 1)*100 + iview+1);
		hit_p->SetStat(1);
		
		// save primary hit
		hit_p -> SetNumCombined(1);
		hit_p -> SetU(u3,0,0);
		hit_p -> SetZ(z3,0,0);
		hit_p -> SetPlane(ipl,0,0);
		hit_p -> SetWire(hitpos_p->GetNwire(),0,0);
		hit_p -> SetTdc(hitpos_p->GetTdc(),0,0);
		
		nok++;
		if(nok == 1) {
		  if(NoSaved>max_saved) printf("NoSaved = %dn",NoSaved);
		  //It is the first one, so increment nhit 
		  //and add detector hit to the list.
		  for(it=0;it<NoSaved;it++) {
		    BrTrackHit2D *tra2d_pp = (BrTrackHit2D*)SavedId.At(it);
		    tra2d_pp->IncNhit();
		    tra2d_pp->AddDetectorHit(hit_p);
		    InsertTrahit(hit_p,tra2d_pp);
		  }
		  diginfo_p->IncUsed(); 
		}
	      }        //if(BrMath::BTWN(lowerLimit,u3,upperLimit))
	    }
	  }            // loop over k
	  fSelectedDCHits->Clear();
	}              // if(nok==0 && fSingleHitsCheckingOption)
      }                // loop over ipv 
      
      ///////////////////////////////////////////////////////////
      
      
      // No hits found for this view.
      // Update the missed hit pointer
      if(nok == 0) {
	nmmiss++;
	if(nmmiss>maxmiss[iv]) {
	  //for(j=0;j<NoSaved;j++) RemoveTrack((BrTrackHit2D*)SavedId.At(j));
	  for(j=0;j<NoSaved;j++)((BrTrackHit2D*)SavedId.At(j))->MarkAsBad();
	  break;
	}
      }           
      //if(DebugLevel() > 2) cout<<"nmiss for tra2d = "<<nmmiss<<endl;
    }
  }
  
  //Clean up tra2d table
  //if(DebugLevel()>5) cout<<"Calling DeleteBadTrack()"<<endl;
  DeleteBadTrack();
  
  SavedId.Clear();
}


//_________________________________________________________________________
 void BrDC2DTrackFinder::InsertTrahit(BrDCHit *hit_p,BrTrackHit2D* tra2dd)
{
//Create a new TrackHit.  Indicate which combined hit it comes from as well
//as which 2D track.  Also add this track to the 2D track from which it came
  BrTrackHit* trahit_p;

/*
  BrClonesArray &trahit = *fTrackHits;
  Int_t ntrahit = fTrackHits->GetLast()+1;
  trahit_p = new(trahit[ntrahit]) BrTrackHit();
*/
  trahit_p = new BrTrackHit();
  fTrackHits->Add(trahit_p);

  trahit_p->SetImod(hit_p->GetImod());
  trahit_p->SetLoctra(0);
  //do we really need these two or does it just arrive because
  //of differences from FORTRAN version?  Do it for the moment.
  trahit_p->SetTra2d(tra2dd);
  trahit_p->SetDetectorHit(hit_p);
}

//_______________________________________________________________________________
 void BrDC2DTrackFinder::CopyTrahit1(BrTrackHit2D* oldTrack,
				 BrTrackHit2D* newTrack)
{
  //Copy trahit entries excluding the last detector hit
  
  Int_t itr;
  
  Int_t NumTrackHitsFound;
  TObjArray trackHitsFound;  //List of track hits having 2DTrack of oldTrack
  BrTrackHit* trahit_p;
  BrDCHit* hit_p;
  Int_t NumTrackHits;
  
//First, figure out where the last entry is
  NumTrackHits = fTrackHits->GetLast()+1;
  for(itr=0;itr<NumTrackHits;itr++) {
     trahit_p = (BrTrackHit*)fTrackHits->At(itr);
     if(trahit_p->GetTra2d() == oldTrack) {
        trackHitsFound.Add(trahit_p);
        }
     }
 
  NumTrackHitsFound = trackHitsFound.GetEntries();
  for(itr=0;itr<NumTrackHitsFound-1;itr++) {
    //hit_p = ((BrTrackHit*)trackHitsFound.At(itr))->GetDetectorHit();
    hit_p = (BrDCHit*)((BrTrackHit*)trackHitsFound.At(itr))->GetDetectorHit();
   
     if(hit_p != 0) {
/*
        BrClonesArray &trahit = *fTrackHits;
        Int_t ntrahit = fTrackHits->GetLast()+1;
        trahit_p = new(trahit[ntrahit]) BrTrackHit();
*/
        trahit_p = new BrTrackHit();
        fTrackHits->Add(trahit_p);
             
        trahit_p->SetDetectorHit(hit_p);
        trahit_p->SetTra2d(newTrack);     
        }
     }
  trackHitsFound.Clear();	//clear so destructor does not delete trahits
}

//_______________________________________________________________________________
 void BrDC2DTrackFinder::RemoveTrack(BrTrackHit2D* tra2dd)
{

  //This routine should be more or less debugged although
  //the results don't match exactly SONATA after the first time
  //through.  I think for the moment that is because of the
  //different way it is done here and in using ADAMO.  This should
  //be monitored
  Int_t NumTrackHits;
  Int_t i;
  BrTrackHit* trahit_p;
  BrDCHit* hit_p;

  NumTrackHits = fTrackHits->GetLast()+1;
  for( i=0;i<NumTrackHits;i++) {
     trahit_p = (BrTrackHit*)fTrackHits->At(i);
     if(trahit_p->GetTra2d() == tra2dd) {

//      This is a shortcut over the original version copied from Sonata.
//      It looped over combined hits, compared each detector hit to the hit
//      in trahit.  I see no reason to do that.  To see how that was done,
//      look at this method in BrLocalTrackDC
        hit_p = (BrDCHit*)trahit_p->GetDetectorHit();
        hit_p->SetStat(BrMath::MOD(hit_p->GetStat(),16) 
                                   + (hit_p->GetStat()/16-1)*16);
//      Now Delete this track.
        fTrackHits->RemoveAt(i);
     }
  }
  fTrackHits->Compress();

  tra2dd->MarkAsBad();
}

//______________________________________________________
 void BrDC2DTrackFinder::DeleteBadTrack()
{
  int i;
  Int_t NumTra2D;
  BrTrackHit2D* tra2d_p;
  /*
  if(DebugLevel()>5) 
    cout<<"DeleteBadTrack: num entries = "
	<<fTrackHits2D->GetEntries()<<endl;
   */
  NumTra2D = fTrackHits2D->GetLast()+1;
  for(i=0;i<NumTra2D;i++) {
     tra2d_p = (BrTrackHit2D*)fTrackHits2D->At(i);
     /*
     if(DebugLevel()>5) 
        cout<<"i = "<< i<<"GetStat()= "<<tra2d_p->GetStat()<<endl;
     */
     if( tra2d_p->IsBad()) {
       tra2d_p = (BrTrackHit2D*) fTrackHits2D->RemoveAt(i);
       delete tra2d_p;
     }
  }
  fTrackHits2D->Compress();
}

//_____________________________________________________________
 void BrDC2DTrackFinder::CleanTracks() 
{

  BrDataTable SavedId;
  Int_t idelArr[50];

  Int_t NumTra2DD = fTrackHits2D->GetLast()+1;
  Int_t Ndel=0;
  for(Int_t i=0;i<NumTra2DD;i++) {
    BrTrackHit2D* track2d_p1 = (BrTrackHit2D*)fTrackHits2D->At(i);
    for(Int_t j=i;j<NumTra2DD;j++) {
      if(i==j)continue;
      BrTrackHit2D* track2d_p2 = (BrTrackHit2D*)fTrackHits2D->At(j);
      // check if tracks have the same hit list
      Int_t NumHits_p1=track2d_p1->GetDetectorHitEntries();
      Int_t NumHits_p2=track2d_p2->GetDetectorHitEntries();
      Int_t nsame=0;
      if(NumHits_p1==NumHits_p2){
	for(Int_t ih1=0;ih1<NumHits_p1;ih1++){
	  for(Int_t ih2=0;ih2<NumHits_p2;ih2++){
	    if(track2d_p1->GetDetectorHitAt(ih1)==track2d_p2->GetDetectorHitAt(ih2)){
	      nsame++;
	    }
	  }
	}
      }
      if(nsame==NumHits_p1 && nsame){
	//find identical track
	//cout<<"identical tracks at "<<i<<" and "<<j<<endl;
	Bool_t AddToDelList=kTRUE;
	for(Int_t idl=0;idl<Ndel;idl++)if(j==idelArr[idl])AddToDelList=kFALSE;
	if(AddToDelList){
	  idelArr[Ndel]=j;
	  Ndel++;
	} 
      }
    }
  }
  //if(Ndel)cout<<"Number of tracks to deleted is "<<Ndel<<endl;
  for(Int_t i=0;i<Ndel;i++){
    //cout<<"idelArr: "<<i<<" "<<idelArr[i]<<endl;
    BrTrackHit2D* track2d= (BrTrackHit2D*)fTrackHits2D->At(idelArr[i]);
    //RemoveTrack(track2d);
    track2d->MarkAsBad();
  }
  DeleteBadTrack();
  
 SavedId.Clear();
}

//  $Log: BrDC2DTrackFinder.cxx,v $
//  Revision 1.9  2002/02/20 13:36:46  ufstasze
//  fCheckSingleHits and fStrictParam can be set according to spectr. angle and run number.
//
//  Revision 1.8  2002/02/15 16:05:15  ufstasze
//  fCheckSingleTrack set to kTRUE
//
//  Revision 1.7  2002/01/26 16:50:49  videbaek
//  Found a nasty memory leak in the DeleteTrack.
//  The bad trackcandidate was removed from the list that was an TObjArray.
//  This does not call the destructor of the object. This in contract to the
//  behaviour of TClonesArray where the destructor is called.
//  Nore the original code was written using TClonesArray.
//
//  Revision 1.6  2001/11/19 17:02:09  ufstasze
//  Added Slope limit on the 2D initial tracks.
//
//  Revision 1.5  2001/08/23 15:50:53  staszel
//  Now we don't use RemoveTracks method. We just
//  mark track as bad. Is significantly saves cpu and
//  almost removed memory leak.
//
//  Revision 1.4  2001/08/23 10:18:03  staszel
//  CheckSingleHits option has been introduced.
//
//  Revision 1.3  2001/08/10 16:23:58  staszel
//  Now we create more initial tracks and tracking efficiency in better.
//
//  Revision 1.2  2001/08/07 13:48:25  ouerdane
//  Cleaned up a lot the DC tracking modules and utilities:
//   Made all utilities deriving from BrModule
//   Set the class version number to 0 (NOT 1 like before)
//   Removed + in the LinkDef file for all modules
//   Added ClassDef in the View Combiner class (didn't exist)
//   Structurized the histogramming section in TDirectories
//
//  Revision 1.1.1.1  2001/06/21 14:55:11  hagel
//  Initial revision of brat2
//
//  Revision 1.8  2001/06/20 11:55:01  staszel
//  BrDetectorHit replaced by BrDCCHit.
//
//  Revision 1.7  2001/06/15 15:37:58  staszel
//  Changes having to do with development of dc tracking package
//
//  Revision 1.6  2001/01/19 17:53:27  staszel
//  Changes having to do with development of BrDCTrackingModule
//
//  Revision 1.5  2000/12/14 15:13:01  staszel
//  none
//
//  Revision 1.4  2000/10/25 17:58:05  hagel
//  Further development of BrDCTrackingModule
//
//  Revision 1.3  2000/10/05 02:21:31  hagel
//  Further development of BrDCTrackingModule
//
//  Revision 1.2  2000/09/29 02:15:59  hagel
//  Changes having to do with development of BrDCTrackingModule
//
//  Revision 1.1  2000/09/19 21:10:55  hagel
//  Initial placeholder revision


This page automatically generated by script docBrat by Christian Holm

Copyright ; 2002 BRAHMS Collaboration <brahmlib@rcf.rhic.bnl.gov>
Last Update on by

Validate HTML
Validate CSS