BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//-------------------------------------------------------------------------
//
// BrBfsTrackingModule
//
// Tracking module for the BFS. The module will setup the module that 
// combines the local tracks found front and behind the D4 or D3 magnet. 
// The module does not instigate
// the local tracking of T2/T3, T4 or T5. This is the responsibility of the 
// framework this will also allow a reconstruction job to substitute 
// the combine track part by another module without any penalty in
// code organisation.
//
// The BrGeometryDbManager much be setup such the 'hidden' 
// calls to the manager  to get the geometry information will work.
//
// Note: When compilers develops, the creation of submodules should be 
//       enclosed in try and catch blocks since the initialization is 
//       done in the constructors and thus has no other mean of 
//       indicating failures.
//
// Requires: BrDetectorTracks Table in inputnode
//

//-------------------------------------------------------------------------
//
// $Id: BrBfsTrackingModule.cxx,v 1.24 2002/08/15 14:21:46 ufstasze Exp $
// $Author: ufstasze $
// $Date: 2002/08/15 14:21:46 $
// $Copyright: 2001 Brahms Collaboration 
//
//-------------------------------------------------------------------------


#ifndef ROOT_TString
#include "TString.h"
#endif
#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef BRAT_BrBfsTrackingModule
#include "BrBfsTrackingModule.h"
#endif
#ifndef BRAT_BrClonesArray
#include "BrClonesArray.h"
#endif
#ifndef BRAT_BrVector3D
#include "BrVector3D.h"
#endif
#ifndef BRAT_BrMath
#include "BrMath.h"
#endif
#ifndef BRAT_BrModuleMatchTrack
#include "BrModuleMatchTrack.h"
#endif
#ifndef BRAT_BrDetectorVolume
#include "BrDetectorVolume.h"
#endif
#ifndef BRAT_BrMagnetVolume
#include "BrMagnetVolume.h"
#endif
#ifndef BRAT_BrEventNode
#include "BrEventNode.h"
#endif
#ifndef BRAT_BrDataTable
#include "BrDataTable.h"
#endif
#ifndef BRAT_BrDetectorTrack
#include "BrDetectorTrack.h"
#endif
#ifndef BRAT_BrSpectrometerTracks
#include "BrSpectrometerTracks.h"
#endif
#ifndef BRAT_BrGeometryDbManager
#include "BrGeometryDbManager.h"
#endif
#ifndef BRAT_BrParameterDbManager
#include "BrParameterDbManager.h"
#endif
#ifndef BRAT_BrUnits
#include "BrUnits.h"
#endif
#ifndef BRAT_BrRunInfoManager
#include "BrRunInfoManager.h"
#endif
#ifndef BRAT_BrBbVertex
#include "BrBbVertex.h"
#endif
#ifndef BRAT_BrSpectrometerTracks
#include "BrSpectrometerTracks.h"
#endif

ClassImp(BrBfsTrackingModule);


//__________________________________________________________________
 BrBfsTrackingModule::BrBfsTrackingModule()
{
  //// default constructor DO NOT USE
  fCombineT2T4 = 0;
  fCombineT3T4 = 0;
  fCombineT4T5 = 0;

  fT2Volume = 0;
  fT3Volume = 0;
  fT4Volume = 0;
  fT5Volume = 0;
  fD3Volume = 0;
  fD4Volume = 0;
  fH1Volume = 0;
  fH2Volume = 0;

  fT2T4Tracks  = 0;
  fT3T4Tracks  = 0;
  fT4T5Tracks  = 0;

  fH1Params    = 0;
  fH2Params    = 0;

  fFrontMode    = 0;
  fCleanUp      = 0;
  fVirtualBack  = kTRUE;
  fVirtualFront = kTRUE;
  fUseDb        = kTRUE;
  fArmsAligned  = kFALSE;
  fKeepMatched  = kFALSE;
  
  fD3FieldScaleFactor = 1.0;
  fD4FieldScaleFactor = 1.0;
  SetFrontMode(kT3,kTRUE);
  SetFrontMode(kT2,kTRUE);
}


//__________________________________________________________________
 BrBfsTrackingModule::BrBfsTrackingModule(const Char_t* name, const Char_t* title)
  : BrModule(name, title)
{
  // ------------------------------------
  // Normal constructor.
  // This Name & title is informational.
  // ------------------------------------
  // ----- track matching modules
  fCombineT4T5 = new BrModuleMatchTrack("T4_T5", "T4 - T5 matching module",
					"T4", "T5", "D4");
  fCombineT2T4 = new BrModuleMatchTrack("T2_T4", "T2 - T4 matching module",
					"T2", "T4", "D3");
  fCombineT3T4 = new BrModuleMatchTrack("T3_T4", "T3 - T4 matching module",
					"T3", "T4", "D3");
  fT2Volume = 0;
  fT3Volume = 0;
  fT4Volume = 0;
  fT5Volume = 0;
  fD3Volume = 0;
  fD4Volume = 0;
  fH1Volume = 0;
  fH2Volume = 0;

  fT2T4Tracks  = 0;
  fT3T4Tracks  = 0;
  fT4T5Tracks  = 0;

  fH1Params    = 0;
  fH2Params    = 0;

  fFrontMode    = 0;
  fCleanUp      = 0;
  fVirtualBack  = kTRUE;
  fVirtualFront = kTRUE;
  fUseDb        = kTRUE;
  fArmsAligned  = kFALSE;
  fKeepMatched  = kFALSE;
  
  fD3FieldScaleFactor = 1.0;
  fD4FieldScaleFactor = 1.0;

  SetFrontMode(kT3,kTRUE);
  SetFrontMode(kT2,kTRUE);
}

//__________________________________________________________________
 BrBfsTrackingModule::~BrBfsTrackingModule()
{
  //// default destructor
  if(fCombineT2T4)
    delete fCombineT2T4;
  if(fCombineT3T4)
    delete fCombineT3T4;
  if(fCombineT4T5)
    delete fCombineT4T5;
}

//__________________________________________________________________
 void BrBfsTrackingModule::Init()
{
  // -----------------------------------------
  // Initialze the matching module and setup 
  // local information for volumes
  // -----------------------------------------

  // -----------------------
  
  BrGeometryDbManager*  geoMan = BrGeometryDbManager::Instance();
  BrParameterDbManager* parMan = BrParameterDbManager::Instance();
  
  fCombineT4T5->Init();
  fCombineT2T4->Init();
  fCombineT3T4->Init();

  // check that both T2 and T3 tracks are not "cleaned up"
  if (TESTBIT(fCleanUp, kT2) && TESTBIT(fCleanUp, kT3)) {
    Abort("Init", "You should not clean up both T2 and T3 tracks! Aborting...");
    return;
  }

  // --- check DB information
  if (fUseDb) {
    BrRunInfoManager* runMan =  BrRunInfoManager::Instance();
    const BrRunInfo* run = runMan->GetCurrentRun();
    if (run->GetRunNo() == -1) {
      Abort("Init", "Run number is -1, check it out: "
	    "register a valid run number!");
      return;
    }
    
    if (run->GetFFSAngle() == run->GetBFSAngle())
      fArmsAligned = kTRUE; 
  }
  
  
  if (!fArmsAligned && fCombineT2T4) {
    Warning("Init", "It makes no sense to run this module "
	    "for T2-T4-T5 when FFS and BFS are not aligned."
	    "Disable T2-T4 matching module and use T3-T4..."); 
    SetFrontMode(kT2, kFALSE);
    delete fCombineT2T4;
  }
  
  // TOF2 stuff
  fH2Volume = (BrDetectorVolume*)geoMan->GetDetectorVolume("BrDetectorVolume","TOF2");
  fH2Params = (BrDetectorParamsTof*)parMan->
    GetDetectorParameters("BrDetectorParamsTof","TOF2");

  // Tracking chambers, magnets, and modules
  fT2Volume = fCombineT2T4->GetFrontVolume();
  fT3Volume = fCombineT3T4->GetFrontVolume();
  fT4Volume = fCombineT4T5->GetFrontVolume();  
  fT5Volume = fCombineT4T5->GetBackVolume();   

  fD3Volume = fCombineT2T4->GetMagnetVolume();
  fD4Volume = fCombineT4T5->GetMagnetVolume();

  fCombineT2T4->SetDebugLevel(DebugLevel());
  fCombineT2T4->SetVerbose(Verbose());
  fCombineT3T4->SetDebugLevel(DebugLevel());
  fCombineT3T4->SetVerbose(Verbose());

  if (fArmsAligned) {
    fH1Volume = (BrDetectorVolume*)geoMan->
      GetDetectorVolume("BrDetectorVolume","TOF1");
    fH1Params = (BrDetectorParamsTof*)parMan->
      GetDetectorParameters("BrDetectorParamsTof","TOF1");
  }

  fCombineT4T5->SetDebugLevel(DebugLevel());
  fCombineT4T5->SetVerbose(Verbose());
}

//___________________________________________________________________
 void BrBfsTrackingModule::SetMagnetFields()
{
  // private method

  // ----------------------------------------------------
  // get field information from database
  // update info each time a field is needed (event loop)
  // ----------------------------------------------------
  if (!fUseDb)
    return;

  BrRunInfoManager* runMan =
    BrRunInfoManager::Instance();

  const BrRunInfo* run = runMan->GetCurrentRun();
  if (run->GetRunNo() == -1) {
    Error("SetMagnetFields", "Run number is -1, check it out");
    return;
  }
  
  // get magnet name
  fD3Volume->SetCurrent(run->GetD3Set(), run->GetD3Pol());
  fD4Volume->SetCurrent(run->GetD4Set(), run->GetD4Pol());
  
  fD3Volume->SetField(fD3Volume->GetField()*fD3FieldScaleFactor);
  fD4Volume->SetField(fD4Volume->GetField()*fD4FieldScaleFactor);
  
  if (DebugLevel() > 10)
    cout << "  D3 field set from database: " << fD3Volume->GetField() << endl
	 << "  D4 field set from database: " << fD4Volume->GetField() << endl;
}

//__________________________________________________________________________
 void BrBfsTrackingModule::Clear()
{
  // --------------------
  // Clear datastructures
  // --------------------

  fT2T4Tracks = 0;
  fT3T4Tracks = 0;
  fT4T5Tracks = 0;

  if (TESTBIT(fFrontMode, kT2)) 
    fCombineT2T4->Clear();
  if (TESTBIT(fFrontMode, kT3)) 
    fCombineT3T4->Clear();

  fCombineT4T5->Clear();
}

//__________________________________________________________________________
 void BrBfsTrackingModule::DefineHistograms()
{
  // -----------------------------------------------------------
  // define histograms for modules owned by this 'super module'
  // -----------------------------------------------------------

  TDirectory* saveDir = gDirectory;
  TDirectory* histDir = gDirectory->mkdir("Tracking_BFS");
  histDir->cd();

  // --- book track matching histograms
  if (TESTBIT(fFrontMode, kT2)) 
    fCombineT2T4->Book();
  if (TESTBIT(fFrontMode, kT3)) 
    fCombineT3T4->Book();
  fCombineT4T5->Book();

  // --- H1 histograms  
  if (fArmsAligned) {
    if (TESTBIT(fFrontMode, kT2)) {
      fT2H1Proj = 
	new TH2F("hT2toH1Proj", "T2 track projection onto H1 plane",
		 200, -fH1Params->GetTofWidth()/2.-5, fH1Params->GetTofWidth()/2.+5,
		 200, -fH1Params->GetSlatHeight()/2.-5, fH1Params->GetSlatHeight()/2.+5);
      fT2H1Slat = new TH1F("hT2toH1Slat", "H1 slats pointed by T2 tracks" , 50, -4.5, 45.5);
    }
    
    if (TESTBIT(fFrontMode, kT3)) {
      fT3H1Proj = 
	new TH2F("hT3toH1Proj", "T3 track projection onto H1 plane",
		 200, -fH1Params->GetTofWidth()/2.-5, fH1Params->GetTofWidth()/2.+5,
		 200, -fH1Params->GetSlatHeight()/2.-5, fH1Params->GetSlatHeight()/2.+5);
      fT3H1Slat = new TH1F("hT3toH1Slat", "H1 slats pointed by T3 tracks" , 50, -4.5, 45.5);
    }
    

    fT4H1Proj = 
      new TH2F("hT4toH1Proj", "T4 track projection onto H1 plane",
	       200, -fH1Params->GetTofWidth()/2.-5, fH1Params->GetTofWidth()/2.+5,
	       200, -fH1Params->GetSlatHeight()/2.-5, fH1Params->GetSlatHeight()/2.+5);
    fT4H1Slat = new TH1F("hT4toH1Slat", "H1 slats pointed by T4 tracks", 50, -4.5, 45.5);
  }
  
  // --- if arms not aligned: vertex stuff (from T3 to BB)

  if (!fArmsAligned) {
    fTrkVX  = new TH1F("hTrkVX", "Track vertex in the X direction", 400, -20, 20);
    fTrkVY  = new TH1F("hTrkVY", "Track vertex in the Y direction", 400, -20, 20);
    fTrkVXY = new TH2F("hTrkVXY","Track vertex in the X, Y plane", 400, -20, 20, 400, -20, 20);
    fTheta  = new TH1F("hTheta", "Polar angle calculated from vertex plane", 200, 0, 20);
    fPhi    = new TH1F("hPhi",   "Azimuthal angle calculated from vertex plane", 300, -30, 30);
  }
  
  // --- H2 stuff

  fT4H2Proj  = new TH2F("hT4toH2Proj", "T4 track projection onto H2 plane",
			200, -fH2Params->GetTofWidth()/2.-5, fH2Params->GetTofWidth()/2.+5,
			200, -fH2Params->GetSlatHeight()/2.-5, fH2Params->GetSlatHeight()/2.+5);
  fT4H2Slat  = new TH1F("hT4toH2Slat", "H2 slats pointed by T4 tracks",
			42, -4.5, 37.5);
  fT5H2Proj  = new TH2F("hT5toH2Proj", "T5 track projection onto H2 plane",
			200, -fH2Params->GetTofWidth()/2.-5, fH2Params->GetTofWidth()/2.+5,
			200, -fH2Params->GetSlatHeight()/2.-5, fH2Params->GetSlatHeight()/2.+5);
  fT5H2Slat  = new TH1F("hT5toH2Slat", "H2 slats pointed by T5 tracks",
			42, -4.5, 37.5);

  // --- momentum stuff  
  fD3Mom   = new TH1F("hD3Momentum",   "Momentum from D3", 400, -20, 20);
  fD4Mom   = new TH1F("hD4Momentum",   "Momentum from D4", 400, -20, 20);
  fD3D4Mom = new TH1F("hD3D4Momentum", "Average Momentum from D3 and D4", 400, -20, 20);
  fMomDiff = new TH1F("hMomentumDifference",   "#DeltaMomentum D4 - D3", 200, -10, 10);
  fMomCorr = new TH2F("hMomentumCorrelation","D4 Momentum vs D3 momentum",
		      300, -15, 15, 300, -15, 15);
  
  fMomDevVsMom = new TH2F("hDpVsP",   "#DeltaMomentum D4 - D3 versus average momentum", 
			  500, -25, 25, 200, -10, 10);
  
  
  // --- path length stuff
  // if FFS and BFS aligned: path from H2 to H1, of not: H2 to primary vtx
  fPathLength = new TH1F("hPathLength", "Path Length from H2 to track vertex or H1", 
			 1500, 500, 2000);
  // From H2 to front chamber
  fPartialLength = new TH1F("hPartialLength", "Path Length from H2 to front chamber", 1000, 100, 1500);
  
  // if we extend D3 tracks to D4 in case no D4 tracks are found
  if (fVirtualBack) {
    fD4Virtual = new TH1F("hD4Virtual", "Number of virtual tracks in D4 (extension of D3 tracks)",
			  10, 0.5, 10.5);
    fD4Tracks  = new TH1F("hD4Tracks", "1: real tracks, 2: virtual track (D3 extension)",
			  20, 0.5, 2.5);
  }
  
    if (fVirtualFront) {
    fD3Virtual = new TH1F("hD3Virtual", "Number of virtual tracks in D3 (extension of D4 tracks down to T3 only)",
			  10, 0.5, 10.5);
    fD3Tracks  = new TH1F("hD3Tracks", "1: real tracks, 2: virtual track (D4 extension down to T3 only)",
			  20, 0.5, 2.5);
  }
  
  if (fVirtualFront || fVirtualBack)
    fBfsTracksL  = new TH2F("hLengthVsVirtual", "Length vs 1: real tracks, 2: virtual track (D3 and D4 extension)",
			    20, 0.5, 2.5, 4000, 800, 2400);
  
  gDirectory = saveDir;
}

//_________________________________________________________________________
 void BrBfsTrackingModule::Begin() 
{
  // Set magnetic fields before each sequence
  SetMagnetFields();
}

//_________________________________________________________________________
 void BrBfsTrackingModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  // ------------------------------------------------------------- 
  //  BFS tracking Event:
  //  Combines T3(T2) and T4, T4 and T5
  //  Match combined tracks via T4
  //  Fill table of BfsTracks and put them to output node  
  // -------------------------------------------------------------
  
  //---------
  // get BB vertex
  
  BrBbVertex* vtx = 0;
  if (!fArmsAligned) {
    vtx = (BrBbVertex*)inNode->GetObject("BB Vertex");
    if (!vtx) vtx = (BrBbVertex*)outNode->GetObject("BB Vertex");
    if (!vtx) Info(50, "Event", "Could not find primary BB Vertex");
  }
  
  // ---- clear clones arrays
  Clear();

  // ---- combine front and back via magnet
  if (TESTBIT(fFrontMode, kT2)) {
    fCombineT2T4->CombineFrontBack(inNode);
    fT2T4Tracks  = fCombineT2T4->GetMatchedTracks();
  }
  if (TESTBIT(fFrontMode, kT3)) {
    fCombineT3T4->CombineFrontBack(inNode);
    fT3T4Tracks  = fCombineT3T4->GetMatchedTracks();
  }
  
  // ---- get matched tracks found
  fCombineT4T5->CombineFrontBack(inNode);
  fT4T5Tracks   = fCombineT4T5->GetMatchedTracks();
  
  if (DebugLevel() > 5) {
    cout << " --- Number of matched tracks: " << endl;
    if (fT2T4Tracks)
      cout << "     T2 - T4 : " << fT2T4Tracks->GetEntries() << endl;
    if (fT3T4Tracks)
      cout << "     T3 - T4 : " << fT3T4Tracks->GetEntries() << endl;
    cout << "     T4 - T5 : " << fT4T5Tracks->GetEntries() << endl;
  }

  // ---- check correlation with H2 geometry and H1 if arms aligned
  if (HistOn()) {
    if (fT2T4Tracks) CheckD3Tracks(fT2T4Tracks, kT2);
    if (fT3T4Tracks) CheckD3Tracks(fT3T4Tracks, kT3);
    CheckD4Tracks();
  }
  
  // ------ clean up , i.e remove T2_T4 track or T3_T4 track from the front 
  // track array if they share the same T4 track
  Bool_t* keepT2T4 = 0;
  Bool_t* keepT3T4 = 0;

  if (fT2T4Tracks) {
    keepT2T4 = new Bool_t [fT2T4Tracks->GetEntries()];
    for (Int_t i=0; i < fT2T4Tracks->GetEntries(); i++)
      keepT2T4[i] = kTRUE;
  }
  
  if (fT3T4Tracks) {
    keepT3T4 = new Bool_t [fT3T4Tracks->GetEntries()];
    for (Int_t i=0; i < fT3T4Tracks->GetEntries(); i++)
      keepT3T4[i] = kTRUE;
  }
  
  if (fCleanUp && fT2T4Tracks && fT3T4Tracks) {
    if (TESTBIT(fCleanUp,kCUT2))
      RemoveDuplicates(fT2T4Tracks, fT3T4Tracks, keepT2T4);
    if (TESTBIT(fCleanUp,kCUT3))
      RemoveDuplicates(fT3T4Tracks, fT2T4Tracks, keepT3T4);
  }

  // -------- match combined tracks
  // check matched tracks T3/T2-T4 && T4-T5

  BrDataTable* bfsTable = new BrDataTable("BfsTracks");

  // DO: here, combining both matched tracks: the criterium is that they should
  // share the same T4 track (might be a bit too simple...)

  if (fT2T4Tracks) CombineMatched(bfsTable, fT2T4Tracks, kT2, keepT2T4);
  if (fT3T4Tracks) CombineMatched(bfsTable, fT3T4Tracks, kT3, keepT3T4);

  if (fVirtualBack) {
    if (fT2T4Tracks) ExtendD3MatchedTracks(bfsTable, fT2T4Tracks, kT2, keepT2T4);
    if (fT3T4Tracks) ExtendD3MatchedTracks(bfsTable, fT3T4Tracks, kT3, keepT3T4);
  }

  if (fVirtualFront)
    ExtendD4MatchedTracks(bfsTable);

  if (fT2T4Tracks) delete [] keepT2T4;
  if (fT3T4Tracks) delete [] keepT3T4;
  
  // 
  // Add selected tracks to the outNode
  //
  
  // fill tables with matched tracks for efficiency module
  if(fKeepMatched){
    BrDataTable* T2T4Tracks = new BrDataTable("MatchedT2T4");
    BrDataTable* T3T4Tracks = new BrDataTable("MatchedT3T4");
    BrDataTable* T4T5Tracks = new BrDataTable("MatchedT4T5");
    Int_t NumT2T4=0;
    if(fT2T4Tracks)NumT2T4=fT2T4Tracks->GetEntries();
    for(Int_t i = 0; i < NumT2T4; i++){
      BrMatchedTrack* mtr = (BrMatchedTrack*)fT2T4Tracks->At(i);
      BrMatchedTrack* newmtr = (BrMatchedTrack*)mtr->Clone();
      T2T4Tracks->Add(newmtr);
    }
    Int_t NumT3T4=0;
    if(fT3T4Tracks)NumT3T4=fT3T4Tracks->GetEntries();
    for(Int_t i = 0; i < NumT3T4; i++){
      BrMatchedTrack* mtr = (BrMatchedTrack*)fT3T4Tracks->At(i);
      BrMatchedTrack* newmtr = (BrMatchedTrack*)mtr->Clone();
      T3T4Tracks->Add(newmtr);
    }
    for(Int_t i = 0; i < fT4T5Tracks->GetEntries(); i++){
      BrMatchedTrack* mtr = (BrMatchedTrack*)fT4T5Tracks->At(i);
      BrMatchedTrack* newmtr = (BrMatchedTrack*)mtr->Clone();
      T4T5Tracks->Add(newmtr);
    }
    
    if(T2T4Tracks)outNode->AddObject(T2T4Tracks);
    if(T3T4Tracks)outNode->AddObject(T3T4Tracks);
    outNode->AddObject(T4T5Tracks);
    
  }

  if (bfsTable->GetEntries())
    outNode->AddObject(bfsTable);
  else {
    delete bfsTable;
    return;
  }
  
  Info(50, "Event", " ----- Loop over %d BFS tracks ", bfsTable->GetEntries()); 

  for (Int_t i = 0; i < bfsTable->GetEntries(); i++) {
    BrBfsTrack* bfstrack = (BrBfsTrack*)bfsTable->At(i);
    bfstrack->SetTrackId(i);
    Double_t d3P = bfstrack->GetD3Momentum();
    Double_t d4P = bfstrack->GetD4Momentum();


    // a lot of things are done in the method TrackToVertex
    // set some important info for PID (H2 to BB, H2 to H1)
    // if the pointer vtx = 0, the calculation won't go so far
    // and will stop at the D1 front plane if FFS and BFS are on the same angle
    // to the front chamber otherwise

    TrackToVertex(bfstrack, vtx);

    if(HistOn()) {
      fD3Mom->Fill(d3P);
      fD4Mom->Fill(d4P);
      fMomDiff->Fill(d4P - d3P);
      fMomCorr->Fill(d3P, d4P);

      if (fVirtualBack || fVirtualFront) {
	Int_t d4virt = ((bfstrack->GetD4Status() >> 4) & 1); // gives 0 or 1
	Int_t d3virt = ((bfstrack->GetD3Status() >> 4) & 1); // gives 0 or 1
	Int_t bfsvirt = (d4virt | d3virt) + 1; // give 1 or 2
	fBfsTracksL->Fill(bfsvirt, bfstrack->GetPathLength());
      }
    }
    
    if (Verbose() >= 50)
      bfstrack->Print("A");
  }
}


//___________________________________________________________________
 void BrBfsTrackingModule::RemoveDuplicates(TClonesArray* ar1, 
					   TClonesArray* ar2, 
					   Bool_t* keep)
{
  // private

  // ----------------------------------------------------------------
  // loop over both clones arrays of matched tracks and tag duplicates
  // to remove
  // ----------------------------------------------------------------

  // chech first only matched tracks that are ok
  Bool_t checked[ar2->GetEntries()];
  for (Int_t j = 0; j < ar2->GetEntries(); j++)
    checked[j] = kFALSE;

  for (Int_t i = 0; i < ar1->GetEntries(); i++) {
    Int_t t4Id1 = ((BrMatchedTrack*)ar1->At(i))->GetBackTrack()->GetID();
    if (((BrMatchedTrack*)ar1->At(i))->GetStatus()!=BrMatchedTrack::kOk) {
      keep[i] = kFALSE;
      continue;
    }

    for (Int_t j = 0; j < ar2->GetEntries(); j++) {
      if (checked[j])
	continue;

      Int_t t4Id2 = ((BrMatchedTrack*)ar2->At(j))->GetBackTrack()->GetID();
      if (((BrMatchedTrack*)ar2->At(j))->GetStatus()!=BrMatchedTrack::kOk) 
	continue;

      if (t4Id1 == t4Id2) {
	keep[i] = kFALSE;
	checked[j] = kTRUE;
      }
    }
  }
}

//___________________________________________________________________
 void BrBfsTrackingModule::CheckD3Tracks(TClonesArray* ar, 
					EFrontMode tid)
{

  // private method: only for checking matched tracks
  // if histos are booked.
  // this method has no other action whatsoever to the
  // BFS track. 

  // ---------------------
  // look at Front-Middle tracks
  // ---------------------
  BrDetectorVolume* frontVolume = fT2Volume;
  if (tid == kT3)
    frontVolume = fT3Volume;

  Int_t n = ar->GetEntries();
  Debug(1, "CheckD3Tracks", " --- Middle matched tracks %d", n);

  BrPlane3D tofPlane(0,0,0, 0,1,0, 1,0,0);
  
  //Loop over tracks
  
  for (Int_t i = 0; i < n; i++) {
    BrMatchedTrack *mtc = (BrMatchedTrack*)ar->At(i);
    
    // ----- get detector tracks (T3/T2 and T4)
    BrDetectorTrack* front = mtc->GetFrontTrack();
    BrDetectorTrack* back  = mtc->GetBackTrack();
    Double_t momentum      = mtc->GetMomentum();

    // ----- check T4 track projection on H2:
    BrLine3D t4Line   = fT4Volume->LocalToGlobal(back->GetTrackLine());
    BrLine3D d4EnLine = fD4Volume->GlobalToLocal(t4Line);
    BrLine3D d4ExLine = fD4Volume->SwimForward(d4EnLine, momentum);
    t4Line = fD4Volume->LocalToGlobal(d4ExLine);
    BrLine3D h2Line(fH2Volume->GlobalToLocal(t4Line));

    BrVector3D proj(tofPlane.GetIntersectionWithLine(h2Line));
    Int_t pSlat = fH2Params->GetSlatNo(proj(0));

    fT4H2Proj->Fill(proj.GetX(), proj.GetY());
    fT4H2Slat->Fill(pSlat);
    
    // ----- check T3/T2 track projection on H1:
    if (!fArmsAligned)
      continue;
    
    BrLine3D frLine   = frontVolume->LocalToGlobal(front->GetTrackLine());
    BrLine3D h1Line(fH1Volume->GlobalToLocal(frLine));
    proj  = tofPlane.GetIntersectionWithLine(h1Line);
    pSlat = fH1Params->GetSlatNo(proj(0));

    if (tid == kT2) {
      fT2H1Proj->Fill(proj.GetX(), proj.GetY());
      fT2H1Slat->Fill(pSlat);
    }
    else {
      fT3H1Proj->Fill(proj.GetX(), proj.GetY());
      fT3H1Slat->Fill(pSlat);
    }
  }
}

//___________________________________________________________________
 void BrBfsTrackingModule::CheckD4Tracks()
{
  // private method

  // ---------------------
  // look at T4-T5 tracks
  // ---------------------
  
  Int_t n = fT4T5Tracks->GetEntries();
  Debug(1, "CheckD4Tracks", " ---  T4-T5 matched tracks %d", n); 
  
  BrPlane3D tofPlane(0,0,0, 0,1,0, 1,0,0);

  //Loop over tracks
  
  for (Int_t i = 0; i < n; i++) {
    BrMatchedTrack *mtc = (BrMatchedTrack*)fT4T5Tracks->At(i);
    
    // ----- get detector tracks (T4-T5)
    BrDetectorTrack* front = mtc->GetFrontTrack();
    BrDetectorTrack* back  = mtc->GetBackTrack();
    Double_t momentum      = mtc->GetMomentum();

    // ----- check T5 track projection on H2:
    BrLine3D t5Line   = fT5Volume->LocalToGlobal(back->GetTrackLine());
    BrLine3D h2Line(fH2Volume->GlobalToLocal(t5Line));
    BrVector3D proj(tofPlane.GetIntersectionWithLine(h2Line));
    Int_t pSlat = fH2Params->GetSlatNo(proj(0));

    fT5H2Proj->Fill(proj.GetX(), proj.GetY());
    fT5H2Slat->Fill(pSlat);
    
    if (!fArmsAligned)
      continue;
    
    // ----- check T4 track projection on H1:
    BrLine3D t4Line   = fT4Volume->LocalToGlobal(front->GetTrackLine());
    BrLine3D d3ExLine = fD3Volume->GlobalToLocal(t4Line);
    BrLine3D d3EnLine = fD3Volume->SwimBack(d3ExLine, momentum);
    t4Line = fD3Volume->LocalToGlobal(d3EnLine);
    BrLine3D h1Line(fH1Volume->GlobalToLocal(t4Line));

    proj  = tofPlane.GetIntersectionWithLine(h1Line);
    pSlat = fH1Params->GetSlatNo(proj(0));

    fT4H1Proj->Fill(proj.GetX(), proj.GetY());
    fT4H1Slat->Fill(pSlat);
  }
}

//___________________________________________________________________
 void BrBfsTrackingModule::TrackToVertex(BrBfsTrack* trk, BrBbVertex* vtx) 
{
  // evaluate path length from H2 to primary vertex and other useful things 
  // like projections on TOF planes, pointed slats, etc.

  Double_t length = 0;
  BrVector3D wrong(-999, -999, -999);

  // length from T4 to T5
  BrMatchedTrack*  t4t5  = trk->GetBfsBackTrack();
  length += t4t5->GetPathLength();

  // other stuff
  BrMatchedTrack*  frt4  = trk->GetBfsFrontTrack();
  BrDetectorTrack* t5Trk = t4t5->GetBackTrack();
  BrDetectorTrack* frTrk = frt4->GetFrontTrack();
  BrDetectorTrack* t4Trk = t4t5->GetFrontTrack();

  // ----- track projection on H2:
  BrPlane3D tofPlane(0,0,0, 0,1,0, 1,0,0);

  BrLine3D t5Line   = fT5Volume->LocalToGlobal(t5Trk->GetTrackLine());
  BrLine3D h2Line(fH2Volume->GlobalToLocal(t5Line));
  BrVector3D h2Proj(tofPlane.GetIntersectionWithLine(h2Line));
  trk->SetProjOnTof(h2Proj);

  // ----- pointed slat ------
  Int_t pslat = fH2Params->GetSlatNo(h2Proj[0]);
  trk->SetPointedSlat(pslat);

  BrVector3D gH2Proj;
  fH2Volume->LocalToGlobal(h2Proj, gH2Proj, 0);
  
  // ----- line from T5 track position to H2 ------
  BrVector3D gT5Pos;
  fT5Volume->LocalToGlobal(t5Trk->GetPos(), gT5Pos, 0);
  
  length += (gH2Proj - gT5Pos).Norm();
  length += fH2Params->GetSlatPos(pslat)[2]; // approximately ok

  // ----- length from T4 to front chamber
  length += frt4->GetPathLength();
  trk->SetPartialPath(length);

  if (HistOn())
    fPartialLength->Fill(length);

  // -------- length from front chamber to track vertex
  BrDetectorVolume* frontVolume = fT3Volume;
  if (trk->T2WasUsed())
    frontVolume = fT2Volume;
  
  BrVector3D frPos;
  frontVolume->LocalToGlobal(frTrk->GetPos(), frPos, 0);

  BrLine3D frLine = frontVolume->LocalToGlobal(frTrk->GetTrackLine());

  if (fArmsAligned) { // H1 calculations
    // -------- projection on TOF1
    BrLine3D h1Line(fH1Volume->GlobalToLocal(frLine));
    BrVector3D h1Proj(tofPlane.GetIntersectionWithLine(h1Line));
    trk->SetProjOnTof1(h1Proj);

    // -------- pointed slat on TOF1
    pslat = fH1Params->GetSlatNo(h1Proj[0]);
    trk->SetTof1PointedSlat(pslat);

    // --------- length from front chamber to H1
    BrVector3D gH1Proj;
    fH1Volume->LocalToGlobal(h1Proj, gH1Proj, 0);
    if (trk->T2WasUsed()) // substract to get H2 to H1 length
      length -= (frPos - gH1Proj).Norm();
    else                  // add
      length += (frPos - gH1Proj).Norm();

    trk->SetPathLength(length);
    trk->SetTrackVertex(gH1Proj);

    if (HistOn())
      fPathLength->Fill(length);

    return; // don't go further, the rest depend on the FFS tracking
  }

  // ----------- swim back to vertex (if vtx exits) if T3 faces beam line directly
  BrVector3D trkVtx(0,0,0);
  BrVector3D bbVtx(0, 0, 0);
  if (vtx) 
    bbVtx.SetZ(vtx->GetZ0());

  BrVector3D dir(0, 0, 1);
  BrPlane3D bbPlane(bbVtx, dir);

  // ----------- track length
  length += (bbVtx - frPos).Norm();
  trk->SetPathLength(length);

  // ----------- intersection with BB vtx plane
  trkVtx = bbPlane.GetIntersectionWithLine(frLine);
  trk->SetTrackVertex(trkVtx);

  // ----------- theta, phi
  BrVector3D trkDir(frPos - bbVtx);
  Float_t theta = TMath::ACos(trkDir(2)/trkDir.Norm());
  Float_t phi   = TMath::ATan(trkDir(1)/trkDir(0));
  trk->SetTheta(theta);
  trk->SetPhi(phi);
  
  if (Verbose() > 15) 
    trk->Print("A");

  if (HistOn()) {
    fTheta->Fill(theta/BrUnits::degree);
    fPhi->Fill(phi/BrUnits::degree);
    fTrkVX->Fill(trkVtx(0));
    fTrkVY->Fill(trkVtx(1));
    fTrkVXY->Fill(trkVtx(0), trkVtx(1));
    fPathLength->Fill(length);
  }
}


//___________________________________________________________________
 void BrBfsTrackingModule::Finish() 
{
  // call finish from combine modules
  
  fCombineT4T5->Finish();
  if (TESTBIT(fFrontMode, kT3))
    fCombineT3T4->Finish(); 
  if (TESTBIT(fFrontMode, kT2))
    fCombineT2T4->Finish(); 
}

//___________________________________________________________________
 void BrBfsTrackingModule::CombineMatched(BrDataTable* bfsTab, 
					 TClonesArray* ar, 
					 EFrontMode tid,
					 Bool_t* keep)
{
  // private method

  // --------------------------------------------
  // This method combines T2/T3-T4 and T4-T5 tracks
  // with (for now) a very simple condition:
  // when the T4 track is in common. If there are
  // more than one candidate, the momentum diff  
  // will determine the choice (choose the minimum)
  // One could also check the geometrical properties
  // (maybe later on, once this class works...)
  // --------------------------------------------

  const Int_t d3trk = ar->GetEntries();
  const Int_t d4trk = fT4T5Tracks->GetEntries();
  
  if (HistOn()) {
    if (fVirtualBack)
      fD4Tracks->Fill(1., fT4T5Tracks->GetEntries());
    if (ar == fT3T4Tracks && fVirtualFront)
      fD3Tracks->Fill(1., ar->GetEntries());
  }
  
  for (Int_t i = 0; i < d4trk; i++) {
    BrMatchedTrack*  d4t = (BrMatchedTrack*)fT4T5Tracks->At(i);
    // check if it's a good matched track
    if (d4t->GetStatus()!=BrMatchedTrack::kOk)
      continue;

    BrDetectorTrack* t4f = d4t->GetFrontTrack();
    Double_t         d4P = d4t->GetMomentum();
    
    for (Int_t j = 0; j < d3trk; j++) {
      BrMatchedTrack*  d3t = (BrMatchedTrack*)ar->At(j);
      // check if it's a good matched track
      if (d3t->GetStatus()!=BrMatchedTrack::kOk)
	continue;

      BrDetectorTrack* t4b = d3t->GetBackTrack();
      Double_t         d3P = d3t->GetMomentum();

      if (!keep[j]) 
	continue;

      // hopefully, the track Id has been set correctly!!
      if (t4b->GetID() != t4f->GetID()) 
      	continue;
      
      BrBfsTrack* trk = new BrBfsTrack();
      trk->SetBfsFrontTrack(d3t);
      trk->SetBfsBackTrack(d4t);
      trk->SetT2WasUsed(kFALSE);
      if (tid == kT2)
	trk->SetT2WasUsed(kTRUE);

      bfsTab->Add(trk);

      Debug(15, "CombineMatched", "D4 momentum : %fnD3 momentum : %f", d4P, d3P);
      if (DebugLevel() >= 15)
	trk->Print();
      
      if (HistOn()) {
	fD3D4Mom->Fill((d4P + d3P)/2.);
	fMomDevVsMom->Fill((d4P + d3P)/2., d4P-d3P);
      }
    }
  }
}

//___________________________________________________________________
 void BrBfsTrackingModule::ExtendD3MatchedTracks(BrDataTable* bfsTab, 
						TClonesArray* ar,
						EFrontMode tid,
						Bool_t* keep)
{
  // private method

  // ------------------------------------------------------------------
  // The purpose of this method is to check that D3 matched tracks that
  // don't belong to a BFS track (cf. CombinedMatched) can still be
  // used by swimming them forward through D4. In this case, a virtual
  // D4 matched track is created and a BFS track is added to the table
  // ------------------------------------------------------------------

  BrPlane3D Plane(0,0,0, 0,1,0, 1,0,0);
  BrPlane3D frontEdge  = fD4Volume->GetEffectiveEdgeEntrancePlane();
  

  const Int_t d3trk = ar->GetEntries();
  const Int_t bfst  = bfsTab->GetEntries();
  Int_t      nD4Trk = fT4T5Tracks->GetEntries();   
  
  // ----- loop over D3 matched tracks
  for (Int_t i = 0; i < d3trk; i++) {
    BrMatchedTrack*  d3t = (BrMatchedTrack*)ar->At(i);
    // check if track is good
    if (d3t->GetStatus()!=BrMatchedTrack::kOk)
      continue;

    Bool_t used = kFALSE;
    
    if (!keep[i])
      continue;

    // ---- check if BFS Track already used this matched track
    for (Int_t bfs = 0; bfs < bfst; bfs++) {
      BrBfsTrack* btrk = (BrBfsTrack*)bfsTab->At(bfs);
      BrMatchedTrack* m = btrk->GetBfsFrontTrack();
      if (m->GetTrackId() == d3t->GetTrackId()) {
	used = kTRUE;
	break;
      }
    }
    
    if (used)  
      continue; // check next D3 matched track
    
    // ---- if NOT used, swim forward through D4 and set a virtual D4 track      

    // ---- calculate d3t to T5 volume
    Double_t momentum = d3t->GetMomentum();
    BrDetectorTrack* T4Track = d3t->GetBackTrack();

    BrLine3D Line   = fT4Volume->LocalToGlobal(T4Track->GetTrackLine());
    BrLine3D EnLine = fD4Volume->GlobalToLocal(Line);
    BrLine3D ExLine = fD4Volume->SwimForward(EnLine, momentum);

    if (fD4Volume->GetSwimStatus(EnLine,ExLine,fCombineT4T5->GetFiducialCutDx(),
				 fCombineT4T5->GetFiducialCutDy()) == kFALSE) 
      continue; // Pawel's choice 
    // DO: do we still want them (like FFS through D1) ?

    BrMatchedTrack* d4t = new BrMatchedTrack();

    Line = fD4Volume->LocalToGlobal(ExLine);
    BrLine3D t5Line = fT5Volume->GlobalToLocal(Line);
    BrDetectorTrack* t5track = new BrDetectorTrack();
    t5track->SetID((Int_t)t5track);

    BrVector3D projT5(Plane.GetIntersectionWithLine(t5Line));

    Float_t slope[3];
    slope[0] = t5Line.GetDirection()[0]/t5Line.GetDirection()[2];
    slope[1] = t5Line.GetDirection()[1]/t5Line.GetDirection()[2];
    slope[2] = 1.;

    t5track->SetPos(projT5);
    t5track->SetAlpha(slope);

    d4t->SetFrontTrack(T4Track);
    d4t->SetBackTrack(t5track);
    d4t->SetMomentum(momentum);
    d4t->SetStatus(BrMatchedTrack::kVirtual);
    d4t->SetDang(0);
    d4t->SetDy(0);
    d4t->SetDaly(0);
    d4t->SetTrackId(nD4Trk);
    d4t->SetExit(ExLine.GetOrigin());
    d4t->SetEntrance(frontEdge.GetIntersectionWithLine(EnLine));
    d4t->SetPathLength(VirtualLength(d4t, fT4Volume, fD4Volume, fT5Volume));

    BrBfsTrack* trk = new BrBfsTrack();
    trk->SetBfsFrontTrack(d3t);
    trk->SetBfsBackTrack(d4t);
    trk->SetT2WasUsed(kFALSE);
    if (tid == kT2)
      trk->SetT2WasUsed(kTRUE);

    if (HistOn())
      fD3D4Mom->Fill(momentum);

    bfsTab->Add(trk);
    nD4Trk++;
  }

  // ---- number of virtual tracks created
  if (HistOn()) {
    fD4Virtual->Fill(nD4Trk - fT4T5Tracks->GetEntries());
    fD4Tracks->Fill(2., nD4Trk - fT4T5Tracks->GetEntries());
  }
}

//___________________________________________________________________
 void BrBfsTrackingModule::ExtendD4MatchedTracks(BrDataTable* bfsTab)
{
  // private method

  // ------------------------------------------------------------------
  // The purpose of this method is to check that D4 matched tracks that
  // don't belong to a BFS track (cf. CombinedMatched) can still be
  // used by swimming them back through D3. In this case, a virtual
  // D3 matched track is created and a BFS track is added to the table
  // Note: I only consider a swim back to T3 (not T2)
  // ------------------------------------------------------------------

  BrPlane3D Plane(0,0,0, 0,1,0, 1,0,0);
  BrPlane3D backEdge  = fD3Volume->GetEffectiveEdgeExitPlane();
  

  const Int_t d4trk = fT4T5Tracks->GetEntries();
  const Int_t bfst  = bfsTab->GetEntries();
  Int_t nD3Trk = 0;
  if (fT3T4Tracks) nD3Trk += fT3T4Tracks->GetEntries();   
  
  // ----- loop over D4 matched tracks
  for (Int_t i = 0; i < d4trk; i++) {
    BrMatchedTrack*  d4t = (BrMatchedTrack*)fT4T5Tracks->At(i);
    if (d4t->GetStatus()!=BrMatchedTrack::kOk)
      continue;
    Bool_t useIt = kTRUE;
    
    // ---- check if BFS Track already used this matched track
    for (Int_t bfs = 0; bfs < bfst; bfs++) {
      BrBfsTrack* btrk = (BrBfsTrack*)bfsTab->At(bfs);
      BrMatchedTrack* m = btrk->GetBfsBackTrack();

      if (m->GetTrackId() == d4t->GetTrackId()) {
	useIt = kFALSE;
	break;
      }
    }
    
    if (!useIt)  
      continue; // check next D4 matched track
    
    // ---- if NOT used, swim back through D3 and set a virtual D3 track      

    // ---- calculate d4t to T3 volume
    Double_t momentum = d4t->GetMomentum();
    BrDetectorTrack* T4Track = d4t->GetFrontTrack();

    BrLine3D Line   = fT4Volume->LocalToGlobal(T4Track->GetTrackLine());
    BrLine3D ExLine = fD3Volume->GlobalToLocal(Line);
    BrLine3D EnLine = fD3Volume->SwimBack(ExLine, momentum);

    if (fD3Volume->GetSwimStatus(EnLine,ExLine,fCombineT3T4->GetFiducialCutDx(),
				 fCombineT3T4->GetFiducialCutDy()) == kFALSE) 
      continue; // Pawel's choice 
    // DO: do we still want them (like FFS through D1) ?

    BrMatchedTrack* d3t = new BrMatchedTrack();

    Line = fD3Volume->LocalToGlobal(EnLine);
    BrLine3D t3Line = fT3Volume->GlobalToLocal(Line);
    BrDetectorTrack* t3track = new BrDetectorTrack();
    t3track->SetID((Int_t)t3track);

    BrVector3D projT3(Plane.GetIntersectionWithLine(t3Line));

    Float_t slope[3];
    slope[0] = t3Line.GetDirection()[0]/t3Line.GetDirection()[2];
    slope[1] = t3Line.GetDirection()[1]/t3Line.GetDirection()[2];
    slope[2] = 1.;
    
    t3track->SetPos(projT3);
    t3track->SetAlpha(slope);

    d3t->SetFrontTrack(t3track);
    d3t->SetBackTrack(T4Track);
    d3t->SetMomentum(momentum);
    d3t->SetStatus(BrMatchedTrack::kVirtual);
    d3t->SetDang(0);
    d3t->SetDy(0);
    d3t->SetDaly(0);
    d3t->SetTrackId(nD3Trk);
    d3t->SetExit(backEdge.GetIntersectionWithLine(ExLine));
    d3t->SetEntrance(EnLine.GetOrigin());
    d3t->SetPathLength(VirtualLength(d3t, fT3Volume, fD3Volume, fT4Volume));
      
    BrBfsTrack* trk = new BrBfsTrack();
    trk->SetBfsFrontTrack(d3t);
    trk->SetBfsBackTrack(d4t);
    trk->SetT2WasUsed(kFALSE);

    if (HistOn())
      fD3D4Mom->Fill(momentum);

    bfsTab->Add(trk);
    nD3Trk++;
  }

  // ---- number of virtual tracks created
  if (HistOn()) {
    Int_t nReal = 0;
    if (fT3T4Tracks) nReal += fT3T4Tracks->GetEntries();   
    fD3Virtual->Fill(nD3Trk - nReal);
    fD3Tracks->Fill(2., nD3Trk - nReal);
  }
}

//_____________________________________________________________________________________
 Float_t BrBfsTrackingModule::VirtualLength(BrMatchedTrack* mtrk, 
					   BrDetectorVolume* frontVol, 
					   BrMagnetVolume* magVol,
					   BrDetectorVolume* backVol)
{
  // private method

  // ------------------------------------------------------------------
  // Evaluate length from T3/T4 middle plane to T4/T5 middle plane for
  // virtual matched tracks created in ExtendBfsFront
  // (stolen from BrModuleMatchTrack)
  // ------------------------------------------------------------------

  Float_t length = (mtrk->GetExit() - mtrk->GetEntrance()).Norm();
  
  BrVector3D gFrPos; frontVol->LocalToGlobal(mtrk->GetFrontTrack()->GetPos(), gFrPos, 0);
  BrVector3D gBkPos; backVol->LocalToGlobal(mtrk->GetBackTrack()->GetPos(), gBkPos, 0);
  
  BrVector3D frPos; magVol->GlobalToLocal(gFrPos, frPos, 0);
  BrVector3D bkPos; magVol->GlobalToLocal(gBkPos, bkPos, 0);
  
  length += (frPos - mtrk->GetEntrance()).Norm();
  length += (bkPos - mtrk->GetExit()).Norm();

  Debug(5, "VirtualLength", "%s matched track virtual length: %f", magVol->GetName(),length);

  return length;
}

//__________________________________________________________________
 void BrBfsTrackingModule::Print(Option_t* option) const 
{
  // 
  TString opt(option);
  opt.ToLower(); 
  
  BrModule::Print(option); 
  if (opt.Contains("d")) 
    cout << endl 
         << "  Original author: Djamel Ouerdane" << endl
         << "  Last Modifications: " << endl 
         << "    $Author: ufstasze $" << endl  
         << "    $Date: 2002/08/15 14:21:46 $"   << endl 
         << "    $Revision: 1.24 $ " << endl  
         << endl 
         << "-------------------------------------------------" << endl
	 << "  D3 Field factor:     " << fD3FieldScaleFactor << endl 
	 << "  D4 Field factor:     " << fD4FieldScaleFactor << endl
	 << "  FFS and BFS aligned? " << (fArmsAligned ? "yes" : "no")
	 << endl  
	 << "  Using DB?            " << (fUseDb ? "yes" : "no") 
	 << endl
	 << "  Using virtual back?  " << (fVirtualBack ? "yes" : "no") 
	 << endl
	 << "  Using T2?            " 
	 << (TESTBIT(fFrontMode, kT2) ?  "yes" : "no") << endl
	 << "  Using T3?            " 
	 << (TESTBIT(fFrontMode, kT3) ?  "yes" : "no") << endl
	 << endl;
}

//__________________________________________________________________
//
//  $Log: BrBfsTrackingModule.cxx,v $
//  Revision 1.24  2002/08/15 14:21:46  ufstasze
//  In ExtendD4MatchedTracks and ExtendD3MatchedTracks we use fiducial cuts taken from BrModuleMatchTracks
//
//  Revision 1.23  2002/08/12 16:20:16  ufstasze
//  Modified section for filling matched track tables
//
//  Revision 1.22  2002/06/25 08:58:33  ufstasze
//  Now, if fKeepMatched is true the matched tracks are writen to the outNode
//
//  Revision 1.21  2002/04/15 08:15:33  bjornhs
//  Added declaration of all volumes indep. of status of T3,T4,T5 combinations.
//  The conditional declaration caused a seg. violation when setting e.g.
//      bfs->SetFrontMode(BrBfsTrackingModule::kT2, kTRUE);
//      bfs->SetFrontMode(BrBfsTrackingModule::kT3, kFALSE);
//
//  Revision 1.20  2002/04/10 11:09:38  ekman
//  Only combine matched tracks with status=BrMatchedTrack::kOk. And only clean up tracks
//  with status kOk.
//
//  Revision 1.19  2002/04/09 14:01:13  ekman
//  Fixed small bug (fVirtualBack was not checked before filling hist -> segment. fault).
//
//  Revision 1.18  2002/04/09 02:08:35  ouerdane
//  extend D4 matched tracks to D3 matched if no D3 matched tracks available
//
//  Revision 1.17  2002/04/03 13:05:36  ouerdane
//  removed correction of local track IDs
//
//  Revision 1.16  2002/04/02 14:42:17  ouerdane
//  a lot of updates added into this somewhat messy module, cf brahms-dev-l email for details
//
//  Revision 1.15  2002/03/06 20:52:38  ekman
//  Added debug message.
//
//  Revision 1.14  2002/02/23 15:41:09  videbaek
//  Found yet another memory leak.
//  The BfsTable was neither deleted nor put onto node if no track exsisted.
//  Please guys, when you write code check with the gObjectTable->Print()
//  and .rootrc settings that there are no Tobject leaks in the modules, and structures.
//
//  Revision 1.13  2001/12/20 15:55:56  ouerdane
//  Added methods SetCombineMode() and ExtendBfsFrontToBfs() originally written by Pawel
//  and debugged for production.
//
//  The idea is that we shouldn't restrict ourselves with combining D3 matched tracks
//  to D4 matched tracks since T5 had a module off for a long time, thus decreasing
//  its efficiency. When a D3 track can't find a D4 track to create a BFS track,
//  it will be swimmed forward through D4. If the swim status is good, a virtual D4
//  track (matched track)is created (reason why I introduced a new bit in the status of the
//  matched track; this bit can be retrieved this way:
//
//    Int_t virt = (d4Track->GetStatus() >> 4) & 1;
//
//  if virt is 0, then it's a real track, if 1, it's a virtual track.
//
//  This procedure gives around twice more BFS tracks. The correlation with H2 hits will remove
//  most of the ghosts (hopefully).
//
//  Revision 1.12  2001/12/13 14:12:55  ouerdane
//  fixed a small bug
//
//  Revision 1.11  2001/12/13 14:10:41  ouerdane
//  call the default value of UseDb in the constructor
//
//  Revision 1.10  2001/12/13 13:36:03  ouerdane
//  completed tracking down to vertex by swimming tracks through D1 and D2
//
//  Revision 1.9  2001/12/07 21:36:59  videbaek
//  Add UseDb method to bypass RunDb. For particular Geant use.
//  Bypass (in FFS) swimstatus if momentum undetermined, but D1 is on. This happens for
//  some zero filed runs. SwinStatus would give an error message per track if not bypassed.
//  complete analysis in ffs even if no Vtx from BB. Set vz to 0 in Ntuples, etc.
//  Needed and useful for Geant + pp running when this comes up.
//
//  Revision 1.8  2001/11/20 13:14:31  ufstasze
//  Added SetD3FieldScaleFactor method and fD3FieldScaleFactor new private member.
//
//  Revision 1.7  2001/11/13 23:02:20  ouerdane
//  track Id set properly
//
//  Revision 1.6  2001/11/05 06:42:21  ouerdane
//  added new classes replacing older spectrometer track modules
//
//  Revision 1.5  2001/10/02 20:22:26  ouerdane
//  Added method SetUseNewBb to use the new BrBbVertex or not.
//  Fixed minor things consequently in the method PathLength()
//
//  Revision 1.4  2001/10/02 02:01:36  ouerdane
//  added track length evaluation from final tof plane to primary vertex for FFS,
//  from H2 to H1 for BFS if arms are not aligned.
//  moved SetMagnetField to event method for simple reason:
//  when dealing with mutiple files (sequences) of more than one run, one would like to
//  update the field if changed. In principle it should in Begin but one might have a file
//  mixing more than one run (which I strongly warn not to do...)
//
//  Revision 1.3  2001/09/10 15:04:10  staszel
//  Corrected histogram name (fMomDevVsMom).
//
//  Revision 1.2  2001/09/07 14:05:13  staszel
//  Front and Back Detector tracks are written to disk and
//  can be access from bfsTrack.
//
//  Revision 1.1  2001/09/06 09:26:13  staszel
//  Initial release.
//

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