BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
// 
// Performe matching between FfsTracks and Td1 counter.
// 
// 

//
// 

//____________________________________________________________________
//
// $Id: BrTd1MatchingModule.cxx,v 1.3 2002/04/02 19:20:53 videbaek Exp $
// $Author: videbaek $
// $Date: 2002/04/02 19:20:53 $
// $Copyright: (C) 2002 BRAHMS Collaboration <brahmlib@rhic.bnl.gov>
//
#ifndef BRAT_BrTd1MatchingModule
#include "BrTd1MatchingModule.h"
#endif

#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif
#ifndef ROOT_TNtuple
#include "TNtuple.h"
#endif
#ifndef ROOT_TH1
#include "TH1.h"
#endif
#ifndef ROOT_TH2
#include "TH2.h"
#endif

#include <BrIostream.h>

#ifndef BRAT_BrTd1Rdo
#include "BrTd1Rdo.h"
#endif

#ifndef BRAT_BrEventNode
#include "BrEventNode.h"
#endif
#ifndef BRAT_BrEvent
#include "BrEvent.h"
#endif
#ifndef BRAT_BrEventNode
#include "BrEventNode.h"
#endif

#ifndef BRAT_BrDataTable
#include "BrDataTable.h"
#endif

#ifndef BRAT_BrDetectorList
#include "BrDetectorList.h"
#endif

#ifndef BRAT_BrPlane3D
#include "BrPlane3D.h"
#endif

#ifndef BRAT_BrTrack
#include "BrTrack.h"
#endif

#ifndef BRAT_BrRunInfoManager
#include "BrRunInfoManager.h"
#endif
#ifndef BRAT_BrSpectrometerTracks
#include "BrSpectrometerTracks.h"
#endif

#ifndef BRAT_BrGeometryDbManager
#include "BrGeometryDbManager.h"
#endif

#ifndef BRAT_BrCoordinateSystem
#include "BrCoordinateSystem.h"
#endif

#ifndef BRAT_BrMagnetVolume
#include "BrMagnetVolume.h"
#endif

#ifndef BRAT_BrTd1TrackMatch
#include "BrTd1TrackMatch.h"
#endif

#include <vector>

//____________________________________________________________________
ClassImp(BrTd1MatchingModule);

//____________________________________________________________________
 BrTd1MatchingModule::BrTd1MatchingModule()
{
  // Default constructor. DO NOT USE
  SetState(kSetup);
}
//____________________________________________________________________
 BrTd1MatchingModule::BrTd1MatchingModule(const Char_t* name, const Char_t* title)
  : BrModule(name, title)
{
  // Named Constructor
  SetState(kSetup);
  fDx=2.4;
  fDy=2.4;
  fXoffset=0.0;
  fYoffset=0.5;
}

//____________________________________________________________________
 void BrTd1MatchingModule::DefineHistograms()
{
  // Define histograms. They are:
  // <fill in here>

  if (GetState() != kInit) {
    Stop("DefineHistograms", "Must be called after Init"); 
    return;  
  }

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

  // Make histograms here 

  fhProjectedPosition = new TH2F("projPosition","Tracks projected to TD1",40, -5, 15, 40, -10, 10);
  fhXYDiff = new TH2F("xyDiff","Difference between tracks and slats",40, -20, 20, 40, -20, 20);
  fhYYDiff1 = new TH2F("yyDiff1","Difference between track vs slat y",80, -10, 10, 80, -10, 10);
  fhYYDiff2 = new TH2F("yyDiff2","Difference between track vs slat y",80, -10, 10, 80, -10, 10);
  fhYYDiff3 = new TH2F("yyDiff3","Difference between track vs slat y",80, -10, 10, 80, -10, 10);
  fhYDiff1 = new TH1F("yDiff1","Difference between track vs slat y",80, -10, 10);
  fhYDiff2 = new TH1F("yDiff2","Difference between track vs slat y",80, -10, 10);
  fhYDiff3 = new TH1F("yDiff3","Difference between track vs slat y",80, -10, 10);
  fhXSlatDiff1 = new TH1F("xDiff1","Difference between Track and Slatpos 1",100,-10,10);
  fhXSlatDiff2 = new TH1F("xDiff2","Difference between Track and Slatpos 2",100,-10,10);
  fhXSlatDiff3 = new TH1F("xDiff3","Difference between Track and Slatpos 3",100,-10,10);
  fhXSlat1 = new TH1F("Slat1","Track pos for Slat1 hit",100,-10,10);
  fhXSlat2 = new TH1F("Slat2","Track pos for Slat2 hit",100,-10,10);
  fhXSlat3 = new TH1F("Slat3","Track pos for Slat3 hit",100,-10,10);
  fhAcpXSlat1 = new TH1F("AcpSlat1","Track pos for Slat1 hit",100,-10,10);
  fhAcpXSlat2 = new TH1F("AcpSlat2","Track pos for Slat2 hit",100,-10,10);
  fhAcpXSlat3 = new TH1F("AcpSlat3","Track pos for Slat3 hit",100,-10,10);

  gDirectory = saveDir;
}

//____________________________________________________________________
 void BrTd1MatchingModule::Init()
{
  // Job-level initialisation
  SetState(kInit);
  BrGeometryDbManager* geomDb = BrGeometryDbManager::Instance();
  fD1Volume  = (BrMagnetVolume*)geomDb->
    GetDetectorVolume("BrMagnetVolume", "D1"); 
  if(!fD1Volume)
    Abort("Init","No Magnet Volume : D1");

  
  fT1Volume  = (BrDetectorVolume*)geomDb->
    GetDetectorVolume("BrDetectorVolume", "T1");
 
  fTd1Volume  = (BrDetectorVolume*)geomDb->
    GetDetectorVolume("BrDetectorVolume", "TD1");
   
  fD1System = const_cast<BrCoordinateSystem*> (fD1Volume->GetCoordinateSystem());
  fT1System = const_cast<BrCoordinateSystem*> (fT1Volume->GetCoordinateSystem());

}

//____________________________________________________________________
 void BrTd1MatchingModule::Begin()
{
  // Run-level initialisation
  SetState(kBegin);
  SetMagnetFields();
}


//___________________________________________________________________
 void BrTd1MatchingModule::SetMagnetFields()
{
  // ----------------------------------------------------
  // get field information from database
  // update info each time a field is needed (event loop)
  // ----------------------------------------------------

  //  if(!fUseDb)
  //  return;
  BrRunInfoManager* runMan =
    BrRunInfoManager::Instance();
  
  if (!runMan) {
    Abort("SetMagnetFields",
	  "Could not get an instance of BrRunInfoManager");
    return;
  }
  
  const BrRunInfo* run = runMan->GetCurrentRun();
  if(!run){
    Abort("SetMagnetField", "Could not get current run");
    return;
  }

  if (run->GetRunNo() == -1) {
    Abort("SetMagnetFields", "Run number is -1, check it out");
    return;
  }
  
  // get magnet settings

  fD1Volume->SetCurrent(run->GetD1Set(), run->GetD1Pol());
  cout << "  D1 field:   " << fD1Volume->GetField() 
       << "  D1 Current: " << run->GetD1Set() 
       << endl;

}

//____________________________________________________________________
 void BrTd1MatchingModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  // 
  // Input:
  //  Td1Rdo datatables.
  //  FFsTracks
  // Output:
  //
  //
  // NOTE: this module assumes that you have the data in the innode.
  // If this is not the case add the module to an AppendContainer rather than
  // the MainModule and you are in business.
  //
  // Method:
  // project tracks (swim back through D1) and looks for slat hit. There is a complication
  // not taken care of namely that the is a physical overlap between the slat 2 and 1 or 3
  // by several mm so there are valid combinations of one track to two slats.
  // The cuts are deliberately loose at this point.
  //
  SetState(kEvent);
  
  BrTd1Rdo* td1 = (BrTd1Rdo*)  inNode->GetObject("Td1Rdo");
  BrDataTable* tracks = inNode->GetDataTable("FfsTracks");
  if(DebugLevel()>10)
    inNode->ListObjects();

  if(!tracks)
    return;
  if(!td1)
    return;

  
  // get midplane for projections.
  // This is in Global system.
  //
  BrPlane3D Td1Plane = fTd1Volume->GetMidPlane();

  int nTracks = tracks->GetEntries();
  
  vector<Int_t> PointedSlat(nTracks);

  for(int i=0;i<nTracks;i++) PointedSlat[i]=0;

  for(int ntr=0;ntr<nTracks;ntr++){
    BrGlbTrack* d2Track = (BrGlbTrack*) tracks->At(ntr);
    BrMatchedTrack*  matchTrack= (BrMatchedTrack*) d2Track->GetMatchedTrack();

    // Fix me - have geometry in DB or files.
    // These coordinates are in the D1 local system (and approximate)
    //
    
    BrDetectorTrack* t1Trk  = matchTrack->GetFrontTrack();
    
    Double_t momentum = matchTrack->GetMomentum();
    
    BrLine3D gTrkLine  = fT1Volume->LocalToGlobal(t1Trk->GetTrackLine());
    BrLine3D d1ExTrkLine = fD1Volume->GlobalToLocal(gTrkLine);
    
    //project to exit of D1 (local to D1)
    BrPlane3D  d1ExPlane = fD1Volume->GetEffectiveEdgeExitPlane();
    BrVector3D d1ExPoint =
      d1ExPlane.GetIntersectionWithLine(d1ExTrkLine);

    d1ExTrkLine.SetOrigin(d1ExPoint);
    BrLine3D gD1ExTrkLine = fD1Volume->LocalToGlobal(d1ExTrkLine);
    BrVector3D gD1ExPoint = gD1ExTrkLine.GetOrigin();
    
    // swim back track through D1
    BrLine3D d1EnTrkLine =
      fD1Volume->SwimBack(d1ExTrkLine, momentum);

    BrLine3D gD1EnTrkLine = fD1Volume->LocalToGlobal(d1EnTrkLine);

    BrVector3D gtrackPosition = Td1Plane.GetIntersectionWithLine(gD1EnTrkLine);
    BrVector3D trackPosition;
    fTd1Volume->GlobalToLocal(gtrackPosition, trackPosition,0);

    // Found projected track to TD1 plane
    // Now compare and see if a match is found in 'Slat pos' and 'Time i.e. y
    //
    if(HistOn())
      fhProjectedPosition->Fill(trackPosition(0),trackPosition(1));
    for(int j=0;j<3;j++){
      if(td1->GetStatus(j+1)){
	float x = trackPosition(0)- td1->GetXpos(j+1);
	float y = trackPosition(1)- td1->GetYpos(j+1);

	if(HistOn()){
	  fhXYDiff->Fill(x,y);
	  if(j==0)
	    {
	      fhYYDiff1->Fill(trackPosition(1),td1->GetYpos(j+1));
	      fhYDiff1->Fill(td1->GetYpos(j+1));
	      fhXSlatDiff1->Fill(trackPosition(0)-td1->GetXpos(j+1));
	      fhXSlat1->Fill(trackPosition(0));
	      if(TMath::Abs(x-fXoffset)<fDx &&
		 TMath::Abs(y-fYoffset)<fDy)
		fhAcpXSlat1->Fill(trackPosition(0));
	    }
	  if(j==1){
	    fhYYDiff2->Fill(trackPosition(1),td1->GetYpos(j+1));
	    fhYDiff2->Fill(td1->GetYpos(j+1));
	    fhXSlatDiff2->Fill(trackPosition(0)-td1->GetXpos(j+1));
	    fhXSlat2->Fill(trackPosition(0));
	    if(TMath::Abs(x-fXoffset)<fDx &&
	       TMath::Abs(y-fYoffset)<fDy)
	      fhAcpXSlat2->Fill(trackPosition(0));
	  }
	  
	  
	  if(j==2){
	    fhYYDiff3->Fill(trackPosition(1),td1->GetYpos(j+1));
	    fhYDiff3->Fill(td1->GetYpos(j+1));
	    fhXSlatDiff3->Fill(trackPosition(0)-td1->GetXpos(j+1));
	    fhXSlat3->Fill(trackPosition(0));
	    if(TMath::Abs(x-fXoffset)<fDx &&
	       TMath::Abs(y-fYoffset)<fDy)
	      fhAcpXSlat3->Fill(trackPosition(0));
	  }
	}
	//
	// performe the matching check.
	//
	if(TMath::Abs(x-fXoffset)<fDx &&
	   TMath::Abs(y-fYoffset)<fDy){
	  PointedSlat[ntr]=j+1;
	  break;
	}
	
      }
    } // end slat loop
  } // end track loop
  
  //
  // Insert Check for Double Hits either tracks and or from Possible loops for double Hits...
  //

  BrDataTable* tab = new BrDataTable("FfsMatchTd1","rdo");
  Int_t nmatch=0;
  for(int ntr=0;ntr<nTracks;ntr++){
    BrGlbTrack* d2Track = (BrGlbTrack*) tracks->At(ntr);
    if(PointedSlat[ntr]>0){
      BrTd1TrackMatch * obj = new BrTd1TrackMatch();
      Int_t slat = PointedSlat[ntr];
      Float_t length = 23.5;
      if(slat==2) length=31.0;
      obj->SetMatching(d2Track->GetTrackId(), slat, td1->GetTime(slat), length);
      tab->Add(obj);
      nmatch++;
    }
    
  }
  if(tab->GetEntries())
    outNode->AddDataTable(tab);
  else
    delete tab;
  
}

//____________________________________________________________________
 void BrTd1MatchingModule::End()
{
  // Run-level finalisation
  SetState(kEnd);
}

//____________________________________________________________________
 void BrTd1MatchingModule::Finish()
{
  // Job-level finalisation
  SetState(kFinish);
}

//____________________________________________________________________
 void BrTd1MatchingModule::Print(Option_t* option) const
{
  // Print module information
  // See BrModule::Print for options.
  // In addition this module defines the Option:
  // <fill in here>

  TString opt(option);
  opt.ToLower(); 
  
  BrModule::Print(option); 
  if (opt.Contains("d")) 
   cout << endl 
         << "  Original author: Flemming Videbaek" << endl
         << "  Last Modifications: " << endl 
         << "    $Author: videbaek $" << endl  
         << "    $Date: 2002/04/02 19:20:53 $"   << endl 
         << "    $Revision: 1.3 $ " << endl  
         << endl 
         << "-------------------------------------------------" << endl;
}

//____________________________________________________________________
//
// $Log: BrTd1MatchingModule.cxx,v $
// Revision 1.3  2002/04/02 19:20:53  videbaek
// Additional diagnostic histograms
//
// Revision 1.2  2002/03/21 14:12:56  videbaek
// Code had an unguarded HistOn() for an histogram increment.
// Fixed.
//
// Revision 1.1  2002/03/20 19:17:56  videbaek
// Added module for matching between ffsTracks and Td1 start counter used
// in the pp. Running. The module is similar in concept to the TofMatching modules
// but utilizes a different intermediat object.
//

//
//
//

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