BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
////////////////////////////////////////////////////////////////
//
// The TBrDigitizeDC is the base module for digitization of drift
//  Chambers. The algorithms are based on those in Sonata and Sonata++
//
//
////////////////////////////////////////////////////////////////

//  $Id: BrDigitizeDC.cxx,v 1.8 2002/06/24 16:45:44 ufstasze Exp $
//
#include "BrDigitizeDC.h"
#include "BrDCPlane.h"
#include "BrDetectorParamsDC.h"
#include "BrDetectorVolume.h"
#include "BrGeometryDbManager.h"
#include "BrParameterDbManager.h"
#include "BrDcDig.h"
#include "BrEventNode.h"
#include "BrMath.h"
#include "TRandom.h"
#include "BrDCHitPosition.h"
#include "BrTableNames.h"
#include "BrDetectorDC.h"
#include "BrDriverDC.h"
#ifndef BRAT_BrVector3D
#include "BrVector3D.h"
#endif
#ifndef BRAT_BrPlane3D
#include "BrPlane3D.h"
#endif
#ifndef BRAT_BrLine3D
#include "BrLine3D.h"
#endif
//#include "BrDCTrack.h"
 
#include "TROOT.h"
#include "TCanvas.h"
#include "TView.h"

#include <iostream.h>

ClassImp(BrDigitizeDC);
 BrDigitizeDC::BrDigitizeDC()
    :BrModule()
{ 
  fDetectorNode = NULL;
  fDCDriver = 0;
}

 BrDigitizeDC::BrDigitizeDC(const Char_t *Name, const Char_t *Title)
  :BrModule(Name,Title)
{
  fParams_p       = NULL;
  for(Int_t i=0;i<4;i++) fVolumeParams_p[i] = NULL;
  fNumVol = 0;
  fDetectorNode = NULL;
  fVolume;
  fDCDriver = 0;
  //fPosres = (Float_t)0.02;
  fPosres = (Float_t)0.012;

  fBragGeoMode=kFALSE;

  if(!strcmp("T3",GetName()))fSortDighits=470; // this parameters are used to decide whether
  if(!strcmp("T4",GetName()))fSortDighits=190; // dighit should be sort or not
  if(!strcmp("T5",GetName()))fSortDighits=190;

}

//____________________________________________________________________________________
 void BrDigitizeDC::DefineHistograms(){
  //
  // Called from Module::Book
  //
  // Create histograms

  TDirectory* saveDir = gDirectory;
  TDirectory* histDir = gDirectory->mkdir(Form("Digitized%s", GetName())); 
  histDir->cd();

  Char_t HistName[32],HistTitle[64];

  sprintf(HistName,"Ghits%s",GetName());
  sprintf(HistTitle,"Geant Hits for %s",GetName()); 
  d_Ghits = new TH1F(HistName,HistTitle,40, (float)0., (float)40.);

  sprintf(HistName,"Nhits%s",GetName());
  sprintf(HistTitle,"Dig Hits for %s",GetName()); 
  d_Hhits = new TH1F(HistName,HistTitle,500, (float)0., (float)1000.);
  for(Int_t idet=0;idet<3;idet++) {
    sprintf(HistName,"Nhits%s_%d",GetName(),idet);
    sprintf(HistTitle,"Hits Number in %d-%s",idet+1,GetName()); 
    fHits[idet] = new TH1F(HistName,HistTitle,501, -0.5, 500.5);
  }
  
  sprintf(HistName,"DZ%s",GetName());
  d_Dz    = new TH1F(HistName,"Distribution of Z",310,-15.0,15.0);

  for(Int_t im=0;im<3;im++){
    for(Int_t ip=0;ip<fDCDriver->GetMaxPlaneNum();ip++){
      sprintf(HistName,"Diffx_m%d_p%d",im+1,ip+1);  
      sprintf(HistTitle,"Spread of drift distans for %s in %d %d (module,plane)",GetName(),im+1,ip+1);
      fhDiffx[im][ip] = new TH1F(HistName,HistTitle,1000,-1.5,1.5);
    }
  }
  
  sprintf(HistName,"Spread%s",GetName());  
  sprintf(HistTitle,"Spread of drift distans for %s",GetName());
  fhSpread = new TH1F(HistName,HistTitle,1000,-1.5,1.5);
  
  sprintf(HistName,"Drift%s",GetName());  sprintf(HistTitle,"Drift distans in %s",GetName());  
  d_Drift = new TH1F(HistName,HistTitle,1000,-.5,1.5);
  
  gDirectory = saveDir;

}

//__________________________________________________________________________
 void BrDigitizeDC::SetDetectorParamsDC(BrDetectorParamsDC* par)
{ 
  // Set the DC parameters with the values of the external parameters
  // referenced by pointer
  // Note that no owner ship is transferred. A copy is made.
  // 
  // 
	
  if( fParams_p != 0)
    delete fParams_p;
  fParams_p = new BrDetectorParamsDC("Dummy","DCINT");
  *fParams_p = *par;
  // The conversion cast is needed because fParams_p is of
  // type BrDetectorParams.
}

 void BrDigitizeDC::SetDetectorParamsDC(const BrDetectorParamsDC& par)
{
  // Set the TPC parameters with the values of the external parameters.
  // Note that no owner ship is transferred. A copy is made.
  // There is no real need for the copy code here. Both Name and Title
  // will be copied (it seems)     
  if( fParams_p != 0)
    delete  fParams_p;
  fParams_p = new BrDetectorParamsDC("Dummy","ParamsTPC");
  *fParams_p = par;
  // The conversion cast is needed because fParams_p is of
  // type BrDetectorParams.
}

 void BrDigitizeDC::ListDetectorParameters() const{
  // List the current value of the digitization parameters on
  // standard out. 
  if( fParams_p != NULL)
    fParams_p->ListParameters();
  else
    cout << "No parameters set for this Digitisation module" << endl;
}


 void BrDigitizeDC::SetDetectorVolume(BrDetectorVolume *vol) {
  Int_t i;
  for(i=0;i<fNumVol;i++) {
    if(!strcmp(fVolumeParams_p[i]->GetName(),vol->GetName())) {
      fVolumeParams_p[i] = vol;
      return;
    }
  }
  if(fNumVol<4) {
    fVolumeParams_p[fNumVol] = vol;
    fNumVol++;
  }
  else 
    cout<<"Too many volume parameters, this one not added"<<endl;
}

//________________________________________________________________
 BrDetectorVolume *BrDigitizeDC::GetDetectorVolume(const Char_t *name) {
  Int_t i;
  for(i=0;i<fNumVol;i++) {
    if(!strcmp(fVolumeParams_p[i]->GetName(),name)) {
      return fVolumeParams_p[i];
    }
  }
  cout<<"Detector Volume "<<name<<" does not exist"<<endl;
  return 0;
}

//________________________________________________________________
 void BrDigitizeDC::ListEventStatistics() {
  BrModule::ListEventStatistics();
  for(Int_t i=0;i<3;i++)
    cout << "NumGeantHits for " << GetName() << " A" << i << " is "
	 << fNumGeantHits[i] << endl;
  cout<<"NumDigHits = "<<fNumDigHits<<endl;
}

//________________________________________________________________
 BrDigitizeDC::~BrDigitizeDC(){
  //
  // Is there a problem with the fact the the fParams_p can point to different kinds
  // of datastructures i.e. potential memory leak. probably yes since fParams_p is
  // defined as BrDetectorParams nor BrDetectorParamsTPC.
  // By second thought it is probably ok. As long as all parameters are defined/set
  // Using the SetDetectorParamsTPC / DC /BB etc..
  //
  if( fParams_p != 0) delete fParams_p;
  for(Int_t i=0;i<4;i++) {
    if( fVolumeParams_p[i] != 0) delete fVolumeParams_p[i];
  }
  if(fDetectorNode) delete fDetectorNode;
}

//________________________________________________________________
 void BrDigitizeDC::Init(){
  // The init member function/method is used to define histograms
  // The histograms is added to the list of hist objects in this
  // module.
  if(DebugLevel())
    cout << "Entering Init " << GetName() << endl;
  Int_t i;
  Char_t volbuf[64];
  BrGeometryDbManager  *gGeomDb;
  BrParameterDbManager *gParamDb;
  BrDetectorVolume     *vol;
  
  // First Set Volumes
  // Create instance of data base managers
  gGeomDb  = BrGeometryDbManager::Instance();
  gParamDb = BrParameterDbManager::Instance();
  
  fParams_p = (BrDetectorParamsDC*) 
    gParamDb->GetDetectorParameters("BrDetectorParamsDC",GetName());
  
  //Get Master volume for DC's
  strcpy(volbuf,GetName());
  vol = (BrDetectorVolume*) 
    gGeomDb->GetDetectorVolume("BrDetectorVolume",volbuf);
  fVolume=(BrDetectorVolume*) 
    gGeomDb->GetDetectorVolume("BrDetectorVolume",GetName());
  if(vol) {
    SetDetectorVolume(vol);
    // Get DC Sub Modules
    for(i=0;i<3;i++) {
      sprintf(volbuf,"%sA%d",GetName(),i+1);
      vol = (BrDetectorVolume*) 
	gGeomDb->GetDetectorVolume("BrDetectorVolume",volbuf);
      if (DebugLevel() > 3) 
	cout << volbuf << " " << vol << endl;
      SetDetectorVolume(vol);
    }
  }

  // define pointer to BrDriverDC
  if(!fDCDriver) {
    if(DebugLevel())cout << "Creating fDCDriver" << endl;
    fDCDriver = new  BrDriverDC(GetName(),GetName());   //Pawel
  }
  fDCDriver->SetSetupStatus(kFALSE);
  //   if(!fDCTrack) {
  //  fDCTrack = new  BrDCTrack((Char_t*)GetName(),(Char_t*)GetName());   //Pawel
  //}
  
}

 void BrDigitizeDC::Event(BrEventNode* InputTable,BrEventNode* OutputTable)
  // Event method for DC digitization
  //
{
  Digitize(InputTable,OutputTable);
}

 Bool_t BrDigitizeDC::Digitize(BrEventNode* InputTable,BrEventNode* OutputTable)
{
  //Digitize DC events
  //return kTRUE if successful
  //return kFALSE if error

  static const Double_t rad = .01745329;
  Int_t i;
  Int_t ihit,iplane;
  Float_t ax,ay,z,cosT,sinT;
  Float_t r,xc,yc,uc,au;
  Float_t xwire,u2,z2;
  
  BrDCPlane* dcpl;
  Int_t numplanes;

  Int_t NumHits;
  Int_t idet;

  Int_t wireArr[10];     //Pawel
  Float_t driftArr[10];  //Pawel
  Float_t detZ;
  // Float_t SubModulePosition[3] = {(Float_t)-30.0,(Float_t)0.0,(Float_t)30.0};
  

  BrDetectorParamsDC *par;
  BrDetectorVolume *vol[4];

  BrGeantHit* ghit_p;
  BrDcDig *dighit_p;

  /*
  Int_t NSameWire[3][10][48];
  for(idet=0;idet<3;idet++) {
    for(Int_t ip=0;ip<10;ip++){
      for(Int_t iw=0;iw<48;iw++)NSameWire[idet][ip][iw]=0;
    }
  }
  */
  

//Get stuff from the input table
  Char_t TableName[32];
  for(idet=0;idet<3;idet++ ) {
    // sprintf(TableName,"GeantHits %sA%dx0",GetName(),idet+1);
    sprintf(TableName,"%s %sA%d",BRTABLENAMES kGeantHits,GetName(),idet+1);
    fGeantHits[idet] = InputTable->GetDataTable(TableName);
    if(fGeantHits[idet] == 0 ) {
      if(DebugLevel() > 0){
	cout << "No hits for " << GetName() << "Module number "<<idet+1<<endl;
      }
      return kFALSE;
    }
  }
  //NumHits = GeantHits->GetEntries();
  
// GetObjectList returns TObjectArray

//Prepare the output table
//sprintf(TableName,"%s DigitizedDCx0",GetName());
  sprintf(TableName,"%s %s",BRTABLENAMES kDigitizedDC,GetName());
  BrDataTable *DigitizedHits = new BrDataTable(TableName);
  OutputTable->AddDataTable(DigitizedHits);

  par = GetDetectorParamsDC();
  if(par == NULL){
    Warning("Digitize"," %s Parameters must be set ",GetName());
    return kFALSE;
  }

  vol[3] = GetDetectorVolume(GetName());
  if(vol[3] == NULL){
    Warning("Digitize"," %s Volume must be set ",GetName());
    return kFALSE;
  }

//Write some data to file for debugging
  /*
Char_t outStr[80];
ofstream outFile;
outFile.open(GetName(),ios::out);
outFile<<"idet iplane        uc          x          y         z"<<endl;
  */
  Int_t DigitHits[3]; // to see howmany hits in each module
  
  //fDCTrack->Event(InputTable);
  if(DebugLevel()>1)
    cout<< "Detector: "<<GetName()<<endl;
  
  for(idet=0;idet<3;idet++) {
    DigitHits[idet]=0;
    
    fDCDriver->AtSubModule(idet);         // Pawel
    
    NumHits = fGeantHits[idet]->GetEntries();
    if(HistOn())
      d_Ghits->Fill(NumHits);
    
    if(DebugLevel()>0) {
      cout<<"Starting Digitization for "<<GetName()<< "NumHits = "<<NumHits<<endl;
      for(ihit=0;ihit<NumHits;ihit++) {
	cout<<(BrGeantHit*)fGeantHits[idet]->At(ihit);
      }
    }
    
    //Initialize things to make g++ happy (or at least not sad)
    au    = (Float_t)0.0;
    xwire = (Float_t)0.0;
    u2    = (Float_t)0.0;
    z2    = (Float_t)0.0;
    
    //detZ = SubModulePosition[idet];	//Need to change!!!
    detZ = fDCDriver->GetSubModulePosition();
    
    numplanes = par->GetNplane();
    
    if(numplanes == 0) {
      cout<<"numplanes = 0!!!, there is some problem with parameters; check it out!!!"<<endl;
      return kFALSE;
    }
    for( ihit=0;ihit<NumHits;ihit++) {
      ghit_p = (BrGeantHit*)fGeantHits[idet]->At(ihit);

      if(fBragGeoMode){
	SetSubModuleTrackParams(ghit_p);
	fDCDriver -> SetSubModuleTrack(fTrPar);
      }else{
	SetDetectorTrackParams(ghit_p);
	fDCDriver -> SetDetectorTrack(fTrPar);
      }
      
      Float_t ax=fTrPar[0];
      Float_t ay=fTrPar[1];
      Float_t x0=fTrPar[2];
      Float_t y0=fTrPar[3];
      //if(fDCTrack -> MatchThis(fDCDriver->GetTrack()) && idet){
      // Float_t *wsk = fDCTrack -> GetTrack(idet);
      //fTrPar[0]=*wsk++;fTrPar[1]=*wsk++;fTrPar[2]=*wsk++;fTrPar[3]=*wsk;
      //fDCDriver -> SetSubModuleTrack(fTrPar);
      //cout<<"Changed : "<<idet<<"  ";fDCDriver->ShowTrack();
      //}else{
      //cout<<"Original: "<<idet<<"  ";fDCDriver->ShowTrack();
      //}
      
      //cout<<"Modul= "<<idet<<" "; fDCDriver -> ShowTrack();          
      fDCDriver -> GoThroughModule(wireArr, driftArr);                     //Pawel  
      
      // Not that follow this point 'wireArr' and 'driftArr' contain wire numbers and
      // the drift distances, respectively. Each wire number refers to the closest wire 
      // in the plane to a given track.
      if(HistOn()){
	d_Dz->Fill(ghit_p->LocalPosIn()[2]);
	d_Dz->Fill(ghit_p->LocalPosOut()[2]);
      }
      
      for(iplane=0;iplane<fDCDriver->GetMaxPlaneNum();iplane++) {
        fDCDriver->AtPlane(iplane);                               //Pawel 
	Int_t iwire=wireArr[iplane];                             //Pawel                   
	if(fDCDriver->WireStatus(iwire)) {  // Wire Status returns null when there is
	  // no such a wire in the chamber
	  //NSameWire[idet][iplane][iwire] += 1;
	  dcpl = par->GetPlaneAt(iplane);
	  if(fBragGeoMode)z = fDCDriver->SubModulePlanePosition();                       //Pawel
	  else z = fDCDriver->GetPlanePosition(idet,iplane); 
	  cosT = (Float_t)TMath::Cos(fDCDriver->GetAngle() * rad);
	  sinT = (Float_t)TMath::Sin(fDCDriver->GetAngle() * rad);
	  if(DebugLevel() > 4) {
	    cout<<"Digitize z position is "<<z<<endl;
	  }
          xc = ax*z + x0;
          yc = ay*z + y0;
	  uc = xc * cosT - yc * sinT; // This 3 lines are for debugging porpose
	  
	  r = gRandom->Rndm();
	  //r = 0.0; // for debugg
	  //            cout<<" Effic= "<<par->GetEff()<<" "<<r<<endl;
	  if(r<= par->GetEff()) {
	    
	    DigitHits[idet]++; // one more hit
	    Float_t spread = fPosres*BrMath::FGauss(); 
	    fhSpread->Fill(spread);
	    //Float_t spread = 0.0; // for debugg
	    Float_t xdrift = TMath::Abs(driftArr[iplane] + spread); //Pawel
	    Float_t time =  xdrift /par->GetDriftv();
	    Int_t adc_time = (Int_t)(time / dcpl->GetTdcConv(iwire) + dcpl->GetTdcOffset(iwire));
	    
	    
	    dighit_p = new BrDcDig();
	    DigitizedHits->Add(dighit_p);
	    dighit_p->SetID(DigitizedHits->GetEntries());
	    dighit_p->SetSubModule(idet+1);	//+1??to be like FORTRAN
	    dighit_p->SetPlane(iplane+1); //should it be +1 to be like FORTRAN
	    dighit_p->SetWire(iwire+1);
	    dighit_p->SetTime(adc_time);
	    dighit_p->SetWidth((Int_t)(par->GetTwoPar() / par->GetDriftv()));
	    
	    //Writing file for debugging
	    //sprintf(outStr,"  %d       %d       %6.2f       %6.2f     %6.2f      %6.2fn",idet,iplane,uc,xc,yc,z+detZ);
	    //outFile<<outStr;
	    
	    //             Here is a small analysis to see how well original position is
	    //             reproduced.
	    xwire = fDCDriver->GetWirePosition(iwire);    //Pawel
	    //if(!strcmp("T5",GetName()) && iplane>1 && iplane<4)
	    //cout<<"Digitize2: "<<idet<<" "<<xwire+xdrift<<" "<<xwire-xdrift<<" "<<iplane<<
	    //                       " "<<fDCDriver->PlanePosition()<<endl;
	    Float_t reco_time;
	    reco_time = (dighit_p->GetTime() - dcpl->GetTdcOffset(iwire)) * dcpl->GetTdcConv(iwire);
	    xdrift = par->GetDriftv() * reco_time;
	    
	    if(HistOn()){
	      d_Drift->Fill(xdrift);              
	      
	      Float_t xspread; 
	      if(uc > xwire)xspread = xwire + xdrift - uc;
	      else xspread = xwire - xdrift - uc;
	      fhDiffx[idet][iplane]->Fill(xspread);              
	    }	      	      
	  } //if r<Eff
	} //if(WireStatus
      } //iplane
    } //ihit
  }  //idet
  //outFile.close();
  
  for(i=0;i<3;i++)fNumGeantHits[i] = fGeantHits[i]->GetEntries();
  
  //  Int_t SameMax=0;  
  //for(idet=0;idet<3;idet++) {
  //  for(Int_t ip=0;ip<10;ip++){
  //    for(Int_t iw=0;iw<48;iw++) {//
  //	if(NSameWire[idet][ip][iw]>1)SameMax += 1;
  //  }
  // }
  //}
  fNumDigHits = DigitizedHits->GetEntries();

  if(HistOn()){
    d_Hhits->Fill(fNumDigHits);
    for(i=0;i<3;i++)fHits[i]->Fill(DigitHits[i]); // number of hits for individual modules
  } 
  //cout<<GetName()<<" Digitize: "<<fNumDigHits<<endl;
  //if(fNumDigHits>fSortDighits) {
    SortDighit(DigitizedHits);
    fNumDigHits = DigitizedHits->GetEntries();
    //}
 
  if(DebugLevel()>4) {
    cout<<"Digitize3: "<<fNumDigHits<<endl;
    for(i=0;i<fNumDigHits;i++) {
      dighit_p = (BrDcDig*)DigitizedHits->At(i);
      idet = dighit_p->GetSubModule()-1;
           fDCDriver->AtSubModule(idet);
           iplane = dighit_p->GetPlane()-1;
           fDCDriver->AtPlane(iplane);                   
           dcpl   = par->GetPlaneAt(iplane);
           Int_t iwire = dighit_p->GetWire()-1;
           xwire = fDCDriver->GetWirePosition(iwire);
           Float_t time = dcpl->GetRealTime(dighit_p->GetTime(),iwire);
           Float_t xdrift = par->GetDriftv() * time;
           if(!strcmp("T5",GetName()) && iplane>1 && iplane<4)
	     cout<<"Digitize3: "<<idet<<" "<<xwire+xdrift<<" "<<xwire-xdrift<<" "<<iplane<<
	                             " "<<fDCDriver->PlanePosition()<<endl;
           //cout<<dighit_p;
    }    
   }
  return kTRUE;
}

//___________________________________________________________
 void BrDigitizeDC::SetDetectorTrackParams(BrGeantHit* hit)
{
  // define track associated with a hit in the local detector frame

  BrPlane3D Plane(0,0,0, 0,1,0, 1,0,0);

  BrLine3D globalLine = hit->GetGlobalLine();
  BrLine3D localLine = fVolume->GlobalToLocal(globalLine);
  
  BrVector3D proj(Plane.GetIntersectionWithLine(localLine));
  if(localLine.GetDirection()[2]==0){
    cout<<"BrDigitizeDC::SetDetectorTrackParams::can not calculate track direction"<<endl;
    return;
  }
  Float_t ax = localLine.GetDirection()[0]/localLine.GetDirection()[2];
  Float_t ay = localLine.GetDirection()[1]/localLine.GetDirection()[2];
  Float_t x0 = proj.GetX();
  Float_t y0 = proj.GetY();

  fTrPar[0]=ax; fTrPar[1]=ay; fTrPar[2]=x0; fTrPar[3]=y0;
}

//___________________________________________________________
 void BrDigitizeDC::SetSubModuleTrackParams(BrGeantHit* hit)
{
  // define track associated with a hit in the local sub-module frame

  Float_t ax = ( hit->LocalPosOut()[0] - hit->LocalPosIn()[0] ) / 
    ( hit->LocalPosOut()[2] - hit->LocalPosIn()[2] );
  Float_t ay = ( hit->LocalPosOut()[1] - hit->LocalPosIn()[1] ) / 
    ( hit->LocalPosOut()[2] - hit->LocalPosIn()[2] );
  
  Float_t x0 = ( hit->LocalPosOut()[0] + hit->LocalPosIn()[0] ) * 0.5;  //Pawel
  Float_t y0 = ( hit->LocalPosOut()[1] + hit->LocalPosIn()[1] ) * 0.5;  //Pawel
  fTrPar[0]=ax; fTrPar[1]=ay; fTrPar[2]=x0; fTrPar[3]=y0;
}

//______________________________________________________________
 void BrDigitizeDC::SortDighit(BrDataTable* DigitizedHits)
{

  BrDetectorParamsDC *par;

  Int_t jhit,khit;
  Int_t NumDigHits;

  BrDcDig* dighit_p1;
  BrDcDig* dighit_p2;

  Int_t iwire1,iwire2;
  Float_t time1,time2;
  BrDCPlane *dcpl1,*dcpl2;

  par = GetDetectorParamsDC();

  NumDigHits = DigitizedHits->GetEntries();

  DigitizedHits->Sort(NumDigHits);

  for( jhit = 0;jhit < NumDigHits-1;jhit++) {
     dighit_p1 = (BrDcDig*)DigitizedHits->At(jhit);
     iwire1 = dighit_p1->GetWire()-1;
     dcpl1 = par->GetPlaneAt(dighit_p1->GetPlane()-1);
     time1 = (dighit_p1->GetTime() - dcpl1->GetTdcOffset(iwire1)) * dcpl1->GetTdcConv(iwire1);
     for( khit = jhit+1;khit < NumDigHits;khit++ ) {
        dighit_p2 = (BrDcDig*)DigitizedHits->At(khit);
        iwire2 = dighit_p1->GetWire()-1;
        dcpl2 = par->GetPlaneAt(dighit_p2->GetPlane()-1);
        time2 = (dighit_p2->GetTime() - dcpl2->GetTdcOffset(iwire2)) * dcpl2->GetTdcConv(iwire2);

        if( dighit_p1->GetSubModule() == dighit_p2->GetSubModule() &&
            dighit_p1->GetPlane() == dighit_p2->GetPlane() &&
            dighit_p1->GetWire() == dighit_p2->GetWire() ) {
	  if( time1 + dighit_p1->GetWidth() > time2 ) {
              dighit_p1->SetWidth((Int_t)time2 + dighit_p2->GetWidth() - (Int_t)time1);
              //dighit_p2->SetIfuse(dighit_p1->GetID());
              DigitizedHits->Remove(dighit_p2);
              DigitizedHits->Compress();
              NumDigHits--;
              }
           }
        }		//loop over khit
     }			//loop over jhit
}

 void BrDigitizeDC::DrawGeantHits() {
  Int_t i;
  BrDetectorVolume *vol[4];
  Char_t volbuf[64];
  Char_t canvasname[64];
  Char_t nodename[64];
  TCanvas *canvas;

  sprintf(canvasname,"%sCanvas",GetName());

  if(!fDetectorNode) {
    canvas = (TCanvas*)gROOT->FindObject(canvasname);
    if(!canvas) canvas = new TCanvas(canvasname,canvasname,50,50,400,400);
    sprintf(volbuf,"DC %s",GetName());
    vol[3] = GetDetectorVolume(GetName());
    for(i=0;i<3;i++) {
      sprintf(volbuf,"%sA%d",GetName(),i+1);
      vol[i] = GetDetectorVolume(volbuf);
    }
    sprintf(nodename,"%sDetectorNode",GetName());

    fDetectorNode = new BrDetectorDC(nodename,nodename,GetName());//,vol[3],vol[0],vol[1],vol[2]);

    fDetectorNode->Draw();
    TView *view = canvas->GetView();

    // Set the view to dimensions we know for DC's; this could be made more general
    // but following looks pretty good.
    view->SetRange((Float_t)-35.0,(Float_t)-35.0,(Float_t)-35.0,(Float_t)35.0,(Float_t)35.0,(Float_t)35.0);
    Int_t irep;
    view->SetView(-(Float_t)90.0,(Float_t)90.0,(Float_t)0.0,irep);


  }

  fDetectorNode->Clear();

  canvas = (TCanvas*)gROOT->FindObject(canvasname);
  //Check if canvas still exists; if not, don't draw the hits.
  if(canvas) {
    canvas->cd();
    fDetectorNode->DrawGeantHits(fGeantHits);
    canvas->Modified();
  }
}

// $Log: BrDigitizeDC.cxx,v $
// Revision 1.8  2002/06/24 16:45:44  ufstasze
// Added fBragGeoMode option to be independent from the brag geometry
//
// Revision 1.7  2002/06/24 16:30:50  staszel
// Added fBragGeoOption
//
// Revision 1.6  2002/06/10 16:38:12  ufstasze
// Splitted decoding to DecodeDat for data and DecodeDig for digitize geant simul.
//
// Revision 1.5  2002/06/09 15:48:32  ouerdane
// added histogram directory for more clarity
//
// Revision 1.4  2002/01/03 19:51:26  cholm
// Prepared to use BrTableNames class (or perhaps BrDetectorList) for table names
//
// Revision 1.3  2001/08/12 13:32:55  cholm
// Changed CTOR prototype to comply with BrModule CTOR (const Char_t*'s instead
// of just Char_t*'s).
//
// Revision 1.2  2001/06/22 17:40:04  cholm
// Changes t odo with data classes renaming.
//
// Revision 1.1.1.1  2001/06/21 14:55:06  hagel
// Initial revision of brat2
//
// Revision 1.44  2001/04/23 18:17:02  videbaek
// Inititalize fDCDriver to 0 in constructor.
// Only increment histograms when HistOn().
//
// Revision 1.43  2001/04/18 20:27:20  videbaek
// Introduced Defined Histograms - and removed defs of those from Init()
//
// Revision 1.42  2001/01/19 17:49:07  staszel
// Changes having to do with development of BrDCTrackingModule
//
// Revision 1.41  2000/11/23 01:32:57  brahmlib
// Cleaned up code in general. BrDriverDC now usess ClassImp/ClassDef.
//
// Revision 1.40  2000/11/22 16:37:40  videbaek
// Moved Init() out of constructor.
// Changed cout << to Warning(..)
//
// Revision 1.39  2000/10/25 17:58:06  hagel
// Further development of BrDCTrackingModule
//
// Revision 1.38  2000/10/25 15:14:53  staszel
// Now pointer to BrBriverDC is defined in the Init method +
// small changes to fit to new BrDriverDC.
//
// Revision 1.37  2000/10/11 14:40:27  staszel
// Changes in Digitize method has been done to implement new class BrDriverDC.
//
// Revision 1.36  2000/09/29 02:15:59  hagel
// Changes having to do with development of BrDCTrackingModule
//
// Revision 1.35  2000/04/06 20:19:15  cholm
// Corrected some `errors'. Many scope of loop variable i, and cast to int,
// both in ReadASCIILine
//
// Revision 1.34  2000/03/17 03:05:15  hagel
// Implemented constant table names
//
// Revision 1.33  1999/12/30 19:43:41  videbaek
// leanup of code comments - a few places
//
// Revision 1.32  1999/08/25 18:33:46  videbaek
// change BTWN to reference BrMath::BTWN
//
// Revision 1.31  1999/08/25 18:26:42  videbaek
// Remoived ref to btwn.h include file
//
// Revision 1.30  1999/04/15 15:33:03  hagel
// Eliminate possibility of negative time in TDC
//
// Revision 1.29  1999/04/07 21:27:03  hagel
// Added diagnostic histogram h_Diffx
//
// Revision 1.28  1999/03/22 17:37:59  hagel
// Add fPosres as data member.  Use it instead of par->GetPosres() to
// simulate DC position resolution.  That decouples this from the
// DC tracking routine.  In principle, it should be the same however.
// Having it decoupled makes some things easier to simulate.
//
// Revision 1.27  1999/03/16 21:04:12  hagel
// 1. Make BrDCPlane inherit from BrPlane3D.  BrDCPlane how has information about
//    where it is located.
// 2. Change BrDigitizeDC and BrLocalTrackDC to make use of new BrDCPlane information.
// 3. Move Decode(...) from BrDigitizeDC to BrLocalTrackDC.
// 4. Change cut parameters to be member data with defaults rather than hardwired.
//
// Revision 1.26  1999/03/07 00:00:43  hagel
//
// Revision 1.25  1999/02/28 22:14:37  hagel
// Changes to make event display more robust and seamless
//
// Revision 1.24  1999/02/24 16:36:03  hagel
// Implemented BrParameterDbManager
//
// Revision 1.23  1999/02/11 22:27:25  hagel
// Changes to make run on Solaris
//
// Revision 1.22  1999/02/11 03:11:43  hagel
// Debug Event Display for Linux
//
// Revision 1.21  1999/01/15 15:36:02  videbaek
// Changes from merging fv and kh files. should be cosmetic only
// Corrected makeNT (non cygnus makefile for win95)
//
// Revision 1.20  1998/12/17 23:27:14  hagel
// Begin to use BrGeometryDbManager in DC stuff
//
// Revision 1.19  1998/12/14 23:37:53  hagel
// cleanup of Initial version of Event Display stuff for DC's
//
// Revision 1.18  1998/12/09 21:00:57  hagel
// Add support for drawing GeantHits
//
// Revision 1.17  1998/12/08 21:32:28  hagel
// Add return kFALSE in Digitize
//
// Revision 1.16  1998/12/04 21:38:40  videbaek
// Changed local track to return properly in no digitize banks present.
// Routine Digitize retruns a bool_t
// Fix type SHLIB->$(SHLIB) ..
//
// Revision 1.15  1998/09/27 16:57:37  alv
// removed EventStatisticsStart and EventStatisticsEnd
// (inherit them from BrModule)
// In ListEventStatistics
//   use BrModule::ListEventStatistics to list CPU time
//
// Revision 1.14  1998/09/08 15:24:36  alv
// Added fNumVol=0 in constructor.
// Corrected calculation of au in Digitize (sinT term was +au, not -ay).
//
// Revision 1.13  1998/08/24 22:52:10  hagel
// Make fDcpl indicate index instead of pointer
//
// Revision 1.12  1998/07/03 15:47:29  hagel
// Clean up g++ warnings
//
// Revision 1.11  1998/07/02 17:54:15  hagel
// Add const to GetDetectorVolume
//
// Revision 1.10  1998/07/01 21:26:03  videbaek
//  fixed some comments
//
//

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