BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
//____________________________________________________________________
//
//	BrTileDigModule
//
// Digitisation Module class for scintilaltor tiles. You have the
// option to use calibrations rather than the ASCII file parameters.
// Set UseCalibrations to kTRUE. 
//

//____________________________________________________________________
//
// $Id: BrTileDigModule.cxx,v 1.14 2002/04/22 16:58:36 cholm Exp $
// $Author: cholm $
// $Date: 2002/04/22 16:58:36 $
// $Copyright: (c) 2001 BRAHMS Collaboration <brahmlib@rcf.rhic.bnl.gov>
//
#include <iostream.h>
#include <fstream.h>
#ifndef ROOT_TSystem
#include "TSystem.h"
#endif
#ifndef ROOT_TMath
#include "TMath.h"
#endif
#ifndef ROOT_TRandom
#include "TRandom.h"
#endif
#ifndef BRAT_BrGeometryDbManager
#include "BrGeometryDbManager.h"
#endif
#ifndef BRAT_BrParameterDbManager
#include "BrParameterDbManager.h"
#endif
#ifndef BRAT_BrCalibrationManager
#include "BrCalibrationManager.h"
#endif
#ifndef BRAT_BrDetectorList
#include "BrDetectorList.h"
#endif
#ifndef BRAT_BrTileDigModule
#include "BrTileDigModule.h"
#endif
#ifndef BRAT_BrTileParameters
#include "BrTileParameters.h"
#endif
#ifndef BRAT_BrDetectorVolume
#include "BrDetectorVolume.h"
#endif
#ifndef BRAT_BrEventNode
#include "BrEventNode.h"
#endif
#ifndef BRAT_BrDataTable
#include "BrDataTable.h"
#endif
#ifndef BRAT_BrGeantHit
#include "BrGeantHit.h"
#endif
#ifndef BRAT_BrTileDig
#include "BrTileDig.h"
#endif
#ifndef   BRAT_BrVector3D
#include "BrVector3D.h"
#endif
#ifndef BRAT_BrTableNames
#include "BrTableNames.h"
#endif
#ifndef BRAT_BrPathManager
#include "BrPathManager.h"
#endif
#ifndef BRAT_BrTileCalibration
#include "BrTileCalibration.h"
#endif
#ifndef ROOT_TDirectory
#include "TDirectory.h"
#endif

//__________________________________________________________________
ClassImp(BrTileDigModule);


//__________________________________________________________________
 BrTileDigModule::BrTileDigModule()
{
  // default
  // Constructor. 
  SetState(kSetup);
  fParamsMultTile_p = 0;
  fCalibration      = 0;
  fGapLimit         = 4096;
  UseCalibrations();
}
//____________________________________________________________________
 BrTileDigModule::BrTileDigModule(const Char_t *Name, const Char_t *Title)
  :BrModule(Name,Title)
{
  // Normal Named constructor
  // 
  SetState(kSetup);
  fParamsMultTile_p = 0;
  fCalibration      = 0;
  fGapLimit         = 4096;
  UseCalibrations();
}

//___________________________________________________________________________
 BrTileDigModule::~BrTileDigModule()
{
};

//___________________________________________________________________________
 void BrTileDigModule::Init()
{
  // Initialize entry point. This is not called at creation time but
  // This would be a good place to setup histograms, statistics variables etc.
  // This could also be done at the time of the constructor. This
  // might likely be more Usefull since one has the opportunity to
  // setup variables before this call (e.g. detector parameters or
  // alike). 
  SetState(kInit);
#if 0
  BrGeometryDbManager *gGeomDb = BrGeometryDbManager::Instance();
  fVolumeParamsMultTile_p = (BrDetectorVolume*) 
    gGeomDb->GetDetectorVolume("BrDetectorVolume", GetName());
#endif 

  BrParameterDbManager *gParamDb = BrParameterDbManager::Instance();
  fParamsMultTile_p = (BrTileParameters*) 
    gParamDb->GetDetectorParameters("BrTileParameters",GetName());


  Int_t ActiveDetTable[6][8];
  Int_t side, loc;
  Int_t i, j;
  for (i=0; i<6; i++) {
    for (j = 0; j<8; j++) {
      // initiallize all channel to # 99 (which does not exist)
      ActiveDetTable[i][j]=99; 
    }
  }

  // 
  const Char_t *bratsys = BrPathManager::Instance()->GetDataDir();
  // const Char_t *bratsys = gSystem->Getenv("BRATSYS");
  char multdetector_file[128];
  sprintf(multdetector_file,"%s/params/ActiveMultDet.dat",bratsys);
  ifstream activeDet(multdetector_file);
  if (!activeDet) {
      cout << "ActiveMultDet.dat file could not be found." << endl;
      cout << "Use the default setting for converting Geant channel"
	   << endl << " to Brat data table channel number" << endl;;
      cout << "User need to check the output file carefully." << endl;
      Int_t k = 0;
      for (i = 0; i<6; i++) {
	for (j = 0; j<8; j++) {
	  ActiveDetTable[i][j]=k++;
	}
      }
  }
  else  {
    Int_t  k;
    Char_t HistName[32], mname[5];
    while ( activeDet >> HistName >> k) {
      if(HistName[0]=='t') {
	sscanf(HistName,"%4s%d_%d",mname,&side,&loc);
	ActiveDetTable[side][(loc-1)]=k;
      }
    }
  }
  
  // Now, with ActiveDetTable, convert the geant sub volue number to 
  // actual FAST Bus ADC channle (location) number
  for (i = 0; i < 48; i++)  {
    side = (Int_t)floor((i/8));
    loc = (Int_t)fmod(i,8);
    GeantFBConvTable[i]=ActiveDetTable[side][loc];
  }
  //cout << "end of Init() " << endl;

  if (DebugLevel() > 9)
    Print("T");  

  
  if (fUseCalibrations) {
    // Calibrations
#ifdef BR_MULT_CAL_TMP
    fCalibration = BrTileTmpCalibration::Instance();
#else
    BrCalibrationManager *manager =
      BrCalibrationManager::Instance();
    
    fCalibration = static_cast<BrTileCalibration*>
      (manager->Register("BrTileCalibration",
			 BrDetectorList::GetDetectorName(kBrTile)));
#endif
    if (!fCalibration)
      Failure("Init", "didn't get calibration from manager");
    
    fCalibration->Use("pedestal");
    fCalibration->Use("pedestalWidth");
    fCalibration->Use("adcGap");
    fCalibration->Use("adcGain");
    // fCalibration->Use("conversionFunc");
    // These are geometry things proper, but put here for now, until a
    // concensus  of how to do the geometry database.
    fCalibration->Use("ringMap");
    fCalibration->Use("rowMap");
    fCalibration->Use("ringPostion");
    fCalibration->Use("tilt");
  }
}

//___________________________________________________________________________
 void BrTileDigModule::DefineHistograms() 
{                     
  // Define  histograms 
  if (GetState() != kInit) {
    Stop("DefineHistograms", "Must be called after Init"); 
    return;  
  }

  if (!fParamsMultTile_p) {
    Failure("DefineHistograms", "no detector parameters");
    return;
  }

  // Get some parameters 
  Int_t    nSingles   = fParamsMultTile_p->GetNoModules();

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


  // 
  fSingleVsHits         = new TH2D("singleVsHits", "Single Hits", 
				   nSingles, 1, nSingles+1, 
				   200, 0, 200);
  fSingleVsEnergy       = new TH2D("singleVsEnergy", 
				   "Single energies",
				   nSingles, 1, nSingles+1, 
				   100, 0, 0.2);
  fSingleVsEnergyPerHit = new TH2D("singleVsEnergyPerHit", 
				   "Single energies per hit",
				   nSingles, 1, nSingles+1, 
				   100, 0, 0.2);
  fSingleVsAdc          = new TH2D("singleVsAdc", "Single ADCs", 
				   nSingles, 1, nSingles+1, 
				   100, 0, 33500);
  fSingleVsAdcPerHit    = new TH2D("singleVsAdcPerHit", 
				   "Single ADCs per hit", 
				   nSingles, 1, nSingles+1, 
				   100, 0, 33500);
  
  gDirectory = saveDir;
}

//___________________________________________________________________________
 void BrTileDigModule::Event(BrEventNode* inNode, BrEventNode* outNode)
{
  // Event method to be called once per event. 
  // This performes the actual slow simulation
  // of the hits on the MultTile.
  //
  // Description of method:
  // Each tile can in fact be hit by multiple hits
  // 1) First a loop over all hits is made to see how many
  //    hits there are per tube. If the is only a single (charged) hit
  //	  The adc and tdc value is calculated from the adc is calculated from
  //    de/dx and flight time
  // 2) For multiple hits each end of the tube is dealt with seperately.
  //	  The fastest signals is found. If the following hits are within
  //    the resolving time / gate width (100n sec) the adc values are
  //    added to those from the first hit.
  // 
  // Some debug stuff 
  SetState(kEvent);
  Debug(5, "Event", "Inside Digitize of %s", GetName());
  // if (DebugLevel() > 4)
  //   inNode->ListObjects();


  // Get simulation parameters
  BrTileParameters* par = GetTileParameters();
  if(!par) {
    Stop("Event","parameters are not set");
    return ;
  }
  if(DebugLevel() > 4)
    par->ListParameters();

  // No of Tubes for a given kind
  Int_t nTubes = par->GetNoModules();

  // Prepare the output table
  BrTileDig* digit = new BrTileDig(BRTABLENAMES kDigTiles,"sim");
  if (!digit) {
    Stop("Event", "could not make output digit");
    return;
  }
  outNode->AddObject(digit);
  
  	
  // The name for input table object for tile 
  TString tableName("GeantHits TILE");
  Debug(5, "Event", "* Inspecting %s", tableName.Data());
  	
  BrDataTable* hits = (BrDataTable*)inNode->GetObject(tableName.Data());    
  if(!hits) {
    Stop("Event", "no hits found (%s)", tableName.Data());
    return;
  }

  // number of hits in the "inNode" for tile
  Int_t numHits = hits->Entries();
  Info(10, "Event", "found %d hits", numHits);

  // Loop variables; 
  Int_t i = 0;

  // Make temporary arrays for the TDC and ADC signals, and initialise
  // them 
  Float_t  badTdc = 10000000;
  UShort_t nHits[nTubes];           // Just for debugging
  Float_t  adc[nTubes];             // The accumelated energy
  Float_t  tdc[nTubes];             // First hit time 
  for (i = 0; i < nTubes; i++) {
    nHits[i] = 0;
    adc[i]   = 0; 
    tdc[i]   = badTdc;
  }
  
 
  // Loop over hits 
  BrGeantHit* hit = 0; 
  for (Int_t i = 0; i < numHits; i++) {
    // Get the ith hit 
    hit = (BrGeantHit*)hits->At(i); 
    if (!hit) {
      Abort("Event", "hit in table is 0, should not happen");
      return;
    }

    // Find the corresponding Fastbus channel 
    Int_t channel = GeantFBConvTable[hit->Isub()-1];
    if (channel == 99) {
      Debug(20, "Event", 
	    "* Active map says tile %d (Fastbus # %d) isn't on",
	    hit->Isub(), channel); 
      // This channel is not on, so we skip this hit
      continue;
    }
    if (fUseCalibrations) {
      if (fCalibration->GetRingMap(channel) < 0 ||
	  fCalibration->GetRingMap(channel) > 8) {
	Debug(10, "Event", "Ring map says tile %d isn't on", 
	      channel);
	continue;
      }
    }

    Debug(10, "Event", "* Found a hit at %d (%d)" , hit->Isub(),
	  channel);
    

    // Count the number of hits
    nHits[channel]++; 

    // Check if this hit was before any of the previous hits in this
    // channel 
    tdc[channel] = TMath::Min(hit->Tof(), tdc[channel]);
    
    // Accumelate the ADC signal from
    adc[channel] += TMath::Max(hit->Dedx(), Float_t(0));
  }
  
  // Now having looped over all the hits, we ready to make the final
  // output digit. 
  for (i = 0; i < nTubes; i++) {
    // Get two random numbers in [0, 1) 
    Float_t a, b;
    gRandom->Rannor(a,b);

    if (fUseCalibrations) {
      if (fCalibration->GetRingMap(i) < 0 ||
	  fCalibration->GetRingMap(i) > 8) {
	Debug(10, "Event", "Ring map says tile %d isn't on", 
	      i);
	continue;
      }
    }

    if (HistOn()) {
      fSingleVsHits->Fill(i,nHits[i]); 
      fSingleVsEnergy->Fill(i,adc[i]);
      if (nHits[i] != 0) 
	fSingleVsEnergyPerHit->Fill(i, adc[i], nHits[i]);
    }

    if (nHits[i] != 0)
      Debug(10, "Event", "Ch: %3d   H: %4d  E: %8f  T:  %8f", 
	    i, nHits[i],  adc[i], tdc[i]);

    // If no time signal was found in this channel, set it to zero. 
    if (tdc[i] == badTdc) 
      tdc[i] = 0;
    // smear the time signal 
    tdc[i] += a * fParamsMultTile_p->GetSigmaTime();
    tdc[i] /= fParamsMultTile_p->GetTdcConv();
    tdc[i] += fParamsMultTile_p->GetTdcOffset(); 

    if (!fUseCalibrations) {
      // Smear the energy signal
      adc[i]  *= fParamsMultTile_p->GetAdcGain();
      adc[i]  += fParamsMultTile_p->GetAdcOffset();
    }
    else {
      adc[i]  *= 1000 / fCalibration->GetAdcGain(i); 
      adc[i]  += gRandom->Gaus(fCalibration->GetPedestal(i), 
			       fCalibration->GetPedestalWidth(i));
      if (adc[i] > fGapLimit) {
	Debug(15, "Event", "Ch: %3d above gap %8f", i, adc[i]);
	adc[i] += fCalibration->GetAdcGap(i);
	Debug(15, "Event", "Ch: %3d corrected gap %8f", i, adc[i]);
      }
    }
    
    // Check if the ADC was saturated, and if so, set it to the
    // saturation value 
    adc[i] = TMath::Min(adc[i], Float_t(32767));

    if (HistOn()) {      
      fSingleVsAdc->Fill(i, adc[i]);
      if (nHits[i] != 0) 
	fSingleVsAdcPerHit->Fill(i, adc[i], nHits[i]);
    }
    
    Debug(10, "Event", "Ch: %3d A: %8f  T:  %8f", 
	  i, adc[i], tdc[i]);

    // Finally set the output
    digit->SetId(i+1, i+1);
    digit->SetAdc(i+1, Int_t(adc[i]));
  }
  if(DebugLevel() > 9) {
    cout << " * Output from BrTileDigModule:" << endl;
    digit->Print("R");
  }
}  

//______________________________________________________________________
 void BrTileDigModule::ListDetectorParameters() const {
  // List the current value of the digitization parameters on
  // standard out. 
  if( fParamsMultTile_p != NULL){
    cout <<    "The name is " << GetName() <<"n";
    fParamsMultTile_p->ListParameters();
  }
  else
    cout << "No parameters set for this Digitisation module" << endl;
}

//____________________________________________________________________________
 void BrTileDigModule::Print(Option_t* option) const
{
  // Standard information printout. 
  // 
  // Options: See BrModule::Print
  //   T       Print conversion tables 
  TString opt(option);
  opt.ToLower();
  
  BrModule::Print(option);
  if (opt.Contains("d"))
    cout << endl
         << "  Original author: Hironori Ito, Version 0.2" << endl
         << "    $Author: cholm $" << endl
         << "    $Date: 2002/04/22 16:58:36 $"   << endl
         << "    $Revision: 1.14 $ " << endl
         << endl
         << "-------------------------------------------------" << endl
	 << "  " << (!fUseCalibrations ? "Not using" : "Using") 
	 << " real calibrations" << endl;
  if (opt.Contains("t")) {
    cout << "  BRAG subid to Fastbus ADC channel conversion: " << endl;
    for (Int_t i = 0; i < 48 ; i++) {
      cout << "    " << setw (3) << i + 1 << ": " 
	   << setw(3) << GeantFBConvTable[i] << flush;
      if ((i + 1) % 8 == 0)
	cout << endl;
    }
    cout << endl;
  }
}

//      $Log: BrTileDigModule.cxx,v $
//      Revision 1.14  2002/04/22 16:58:36  cholm
//      Some fixes for the modes using the calibrations rahter than
//      parameters.  Also changed a lot of
//
//        if (DebuLevel() >  ...)
//           cout << ....
//
//      to
//
//        Debug(..., ..., ...)
//
//      and similar for Info.
//
//      Revision 1.13  2002/04/19 21:33:15  hito
//      Minor change for compiling with gcc3.
//
//      Revision 1.12  2002/01/03 19:51:41  cholm
//      Prepared to use BrTableNames class (or perhaps BrDetectorList) for table names
//
//      Revision 1.11  2001/10/08 11:28:56  cholm
//      Changed to use new DB access classes
//
//      Revision 1.10  2001/09/20 13:42:38  cholm
//      Removed the provisions to use temp ASCII files, since it's no longer
//      needed. Also added some doc.
//
//      Revision 1.9  2001/08/30 12:40:10  cholm
//      Added some more documentation.
//
//      Revision 1.8  2001/08/30 12:37:38  cholm
//      Added option to use calibration gains, pedestals, and so on, rather
//      than the ones in the ASCII parameter files.  This is to make the
//      digitised data more closely resemple actual data.  It gives much
//      better results, especially since the dual range ADC gab is simulated
//      too.
//
//      Revision 1.7  2001/08/24 18:51:43  cholm
//      Better histogram limits. Higher debug level for messages.
//
//      Revision 1.6  2001/08/21 12:43:00  cholm
//      Increased some debugging levels
//
//      Revision 1.5  2001/08/12 13:15:29  cholm
//      Fixed a few things.
//
//      Revision 1.4  2001/08/12 09:43:42  cholm
//      Filled in the DefineHistograms method.  My mistake I overlooked this
//      file when commiting changes - sorry.
//
//      Revision 1.2  2001/08/10 14:27:12  cholm
//      Cleaned up this module a lot.  The event method was very inefficient,
//      and there seemed to be a problem with setting the ID in the digit, so
//      that too many channels had ID of 0 (zero).
//
//      Revision 1.1.1.1  2001/06/21 14:55:07  hagel
//      Initial revision of brat2
//
//      Revision 1.5  2001/06/04 13:37:27  cholm
//      Changes to use BrPathManager, BrVersion, and perpare to use BrFileTag.
//
//      Revision 1.4  2001/05/03 18:45:50  cholm
//      Fixed a few superflous warnings, and uninitialised array in BrTileParameters
//
//      Revision 1.3  2001/03/07 12:20:22  cholm
//      * Made the method BrModule::Info() const obsolete in favour of
//        BrModule::Print(Option_t* option="B") const.
//
//      Revision 1.2  2001/02/02 19:50:51  videbaek
//      Made Warning contingent on Verbose
//
//      Revision 1.1  2001/01/29 20:27:22  cholm
//      Renamed class BrDigitizeMultTile to BrTileDigModule to fit with naming
//      scheme.
//
//      Revision 1.10  2001/01/17 02:01:48  hagel
//      Fixed float to int warning that made into error
//
//      Revision 1.9  2000/11/22 21:16:44  videbaek
//      Moved geometry and parameter calls from constructor to Initm method
//
//      Revision 1.8  2000/09/21 20:49:44  cholm
//      Several changes, since this was really messy code. Most important, moved
//      the arrays of constant size to pointers.
//
//      Revision 1.7  2000/08/19 21:01:44  hito
//      Minor update.
//
//      Revision 1.6  2000/08/17 14:50:43  hito
//      *** empty log message ***
//
//      Revision 1.5  2000/07/26 22:08:37  hito
//      Mapping of digitization for si and tile were changed to reflect the raw data
//      and the change in gbrahms.
//
//      Revision 1.4  2000/05/24 14:11:01  cholm
//      Use name for object from "BrTableNames.h"
//
//      Revision 1.3  2000/03/21 21:22:25  cholm
//      Several changes: A few hacks where needed to compile on Digital Unix, noticably in my - Christian Holm - classes for the DB and Exceptions. Proberly still need to do something to some of Konstantins stuff. Also, I removed a lot of warnings fromt the compiler, by teddying up the code. Several places, old C-style indecies in for loops were assumed. Corrected. Unused variables have been commented out.
//
//      Revision 1.2  2000/03/17 16:13:21  cholm
//      Took out failure on not finding volume, since I don't believe that info is written to the GBRAHMS .geo output file, hence tha module would always fail.
//
//      Revision 1.1  1999/11/30 22:22:17  videbaek
//      dded new class to digitize Multiplicity tiles.
//
//
//	
//

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