Newer
Older
/** \class THcHodoscope
\ingroup Detectors
Class for a generic hodoscope consisting of multiple
planes with multiple paddles with phototubes on both ends.
This differs from Hall A scintillator class in that it is the whole
hodoscope array, not just one plane.
*/
#include "THcCherenkov.h"
#include "THcHallCSpectrometer.h"
#include "THcHitList.h"
#include "THcRawShowerHit.h"
#include "TClass.h"
#include "math.h"
#include "THaSubDetector.h"
#include "THcHodoscope.h"
#include "THaEvData.h"
#include "THaDetMap.h"
#include "THcDetectorMap.h"
#include "THcGlobals.h"
#include "THcParmList.h"
#include "VarDef.h"
#include "VarType.h"
#include "THaTrack.h"
#include "TClonesArray.h"
#include "TMath.h"
#include <cstring>
#include <cstdio>
#include <cstdlib>
#include <iostream>
#include <iomanip>
Gabriel Niculescu
committed
#include <fstream>
using namespace std;
//_____________________________________________________________________________
THcHodoscope::THcHodoscope( const char* name, const char* description,
THaApparatus* apparatus ) :
THaNonTrackingDetector(name,description,apparatus)
{
// Constructor
//fTrackProj = new TClonesArray( "THaTrackProj", 5 );
// Construct the planes
Stephen A. Wood
committed
fNPlanes = 0; // No planes until we make them
//_____________________________________________________________________________
THcHodoscope::THcHodoscope( ) :
THaNonTrackingDetector()
{
// Constructor
}
//_____________________________________________________________________________
void THcHodoscope::Setup(const char* name, const char* description)
{
// static const char* const here = "Setup()";
// static const char* const message =
// "Must construct %s detector with valid name! Object construction failed.";
Stephen A. Wood
committed
cout << "In THcHodoscope::Setup()" << endl;
// Base class constructor failed?
if( IsZombie()) return;
fDebug = 1; // Keep this at one while we're working on the code
char prefix[2];
prefix[0]=tolower(GetApparatus()->GetName()[0]);
prefix[1]='\0';
TString temp(prefix[0]);
fSHMS=kFALSE;
if (temp == "p" ) fSHMS=kTRUE;
cout << " fSHMS = " << fSHMS << endl;
DBRequest listextra[]={
{"hodo_num_planes", &fNPlanes, kInt},
{"hodo_plane_names",&planenamelist, kString},
//fNPlanes = 4; // Default if not defined
gHcParms->LoadParmValues((DBRequest*)&listextra,prefix);
cout << "Plane Name List : " << planenamelist << endl;
vector<string> plane_names = vsplit(planenamelist);
if(plane_names.size() != (UInt_t) fNPlanes) {
cout << "ERROR: Number of planes " << fNPlanes << " doesn't agree with number of plane names " << plane_names.size() << endl;
// Should quit. Is there an official way to quit?
}
fPlaneNames = new char* [fNPlanes];
for(Int_t i=0;i<fNPlanes;i++) {
fPlaneNames[i] = new char[plane_names[i].length()+1];
strcpy(fPlaneNames[i], plane_names[i].c_str());
}
// Probably shouldn't assume that description is defined
char* desc = new char[strlen(description)+100];
fPlanes = new THcScintillatorPlane* [fNPlanes];
for(Int_t i=0;i < fNPlanes;i++) {
strcpy(desc, description);
strcat(desc, " Plane ");
strcat(desc, fPlaneNames[i]);
fPlanes[i] = new THcScintillatorPlane(fPlaneNames[i], desc, i+1, this); // Number planes starting from zero!!
cout << "Created Scintillator Plane " << fPlaneNames[i] << ", " << desc << endl;
// --------------- To get energy from THcShower ----------------------
const char* shower_detector_name = "cal";
// THaApparatus* app;
THcHallCSpectrometer *app = dynamic_cast<THcHallCSpectrometer*>(GetApparatus());
THaDetector* det = app->GetDetector( shower_detector_name );
if( dynamic_cast<THcShower*>(det) ) {
fShower = dynamic_cast<THcShower*>(det);
}
else if( !dynamic_cast<THcShower*>(det) ) {
cout << "Warining: calorimeter analysis module "
<< shower_detector_name << " not loaded for spectrometer "
<< prefix << endl;
// --------------- To get energy from THcShower ----------------------
// --------------- To get NPEs from THcCherenkov -------------------
const char* chern_detector_name = "cer";
THaDetector* detc = app->GetDetector( chern_detector_name );
if( dynamic_cast<THcCherenkov*>(detc) ) {
fChern = dynamic_cast<THcCherenkov*>(detc);
}
else if( !dynamic_cast<THcCherenkov*>(detc) ) {
cout << "Warining: Cherenkov detector analysis module "
<< chern_detector_name << " not loaded for spectrometer "
<< prefix << endl;
// --------------- To get NPEs from THcCherenkov -------------------
fScinShould = 0;
fScinDid = 0;
gHcParms->Define(Form("%shodo_did",prefix),"Total hodo tracks",fScinDid);
gHcParms->Define(Form("%shodo_should",prefix),"Total hodo triggers",fScinShould);
// Save the nominal particle mass
fPartMass = app->GetParticleMass();
fBetaNominal = app->GetBetaAtPcentral();
}
//_____________________________________________________________________________
THaAnalysisObject::EStatus THcHodoscope::Init( const TDatime& date )
{
cout << "In THcHodoscope::Init()" << endl;
Setup(GetName(), GetTitle());
char EngineDID[] = "xSCIN";
EngineDID[0] = toupper(GetApparatus()->GetName()[0]);
if( gHcDetectorMap->FillMap(fDetMap, EngineDID) < 0 ) {
static const char* const here = "Init()";
Error( Here(here), "Error filling detectormap for %s.", EngineDID );
return kInitError;
}
// Should probably put this in ReadDatabase as we will know the
// maximum number of hits after setting up the detector map
// But it needs to happen before the sub detectors are initialized
// so that they can get the pointer to the hitlist.
InitHitList(fDetMap, "THcRawHodoHit", fDetMap->GetTotNumChan()+1);
EStatus status;
// This triggers call of ReadDatabase and DefineVariables
if( (status = THaNonTrackingDetector::Init( date )) )
return fStatus=status;
for(Int_t ip=0;ip<fNPlanes;ip++) {
if((status = fPlanes[ip]->Init( date ))) {
return fStatus=status;
}
}
// Replace with what we need for Hall C
// const DataDest tmp[NDEST] = {
// { &fRTNhit, &fRANhit, fRT, fRT_c, fRA, fRA_p, fRA_c, fROff, fRPed, fRGain },
// { &fLTNhit, &fLANhit, fLT, fLT_c, fLA, fLA_p, fLA_c, fLOff, fLPed, fLGain }
// };
// memcpy( fDataDest, tmp, NDEST*sizeof(DataDest) );
fNScinHits = new Int_t [fNPlanes];
fGoodPlaneTime = new Bool_t [fNPlanes];
fNPlaneTime = new Int_t [fNPlanes];
fSumPlaneTime = new Double_t [fNPlanes];
// Double_t fHitCnt4 = 0., fHitCnt3 = 0.;
// fScinHit = new Double_t*[fNPlanes];
// for ( m = 0; m < fNPlanes; m++ ){
// fScinHit[m] = new Double_t[fNPaddle[0]];
// }
for (int ip=0; ip<fNPlanes; ++ip) {
fScinHitPaddle.push_back(std::vector<Int_t>(fNPaddle[ip], 0));
}
return fStatus = kOK;
}
//_____________________________________________________________________________
Int_t THcHodoscope::ReadDatabase( const TDatime& date )
{
// Read this detector's parameters from the database file 'fi'.
// This function is called by THaDetectorBase::Init() once at the
// beginning of the analysis.
// 'date' contains the date/time of the run being analyzed.
// static const char* const here = "ReadDatabase()";
char prefix[2];
char parname[100];
// Read data from database
// Pull values from the THcParmList instead of reading a database
// file like Hall A does.
// Will need to determine which spectrometer in order to construct
// the parameter names (e.g. hscin_1x_nr vs. sscin_1x_nr)
prefix[0]=tolower(GetApparatus()->GetName()[0]);
prefix[1]='\0';
strcpy(parname,prefix);
strcat(parname,"scin_");
// Int_t plen=strlen(parname);
cout << " readdatabse hodo fnplanes = " << fNPlanes << endl;
fFPTime = new Double_t [fNPlanes];
fPlaneCenter = new Double_t[fNPlanes];
fPlaneSpacing = new Double_t[fNPlanes];
prefix[0]=tolower(GetApparatus()->GetName()[0]);
//
prefix[1]='\0';
for(Int_t i=0;i<fNPlanes;i++) {
DBRequest list[]={
{Form("scin_%s_nr",fPlaneNames[i]), &fNPaddle[i], kInt},
{0}
};
gHcParms->LoadParmValues((DBRequest*)&list,prefix);
// GN added
// reading variables from *hodo.param
Gabriel Niculescu
committed
fMaxScinPerPlane=fNPaddle[0];
for (Int_t i=1;i<fNPlanes;i++) {
Gabriel Niculescu
committed
fMaxScinPerPlane=(fMaxScinPerPlane > fNPaddle[i])? fMaxScinPerPlane : fNPaddle[i];
}
// need this for "padded arrays" i.e. 4x16 lists of parameters (GN)
fMaxHodoScin=fMaxScinPerPlane*fNPlanes;
Gabriel Niculescu
committed
if (fDebug>=1) cout <<"fMaxScinPerPlane = "<<fMaxScinPerPlane<<" fMaxHodoScin = "<<fMaxHodoScin<<endl;
fHodoVelLight=new Double_t [fMaxHodoScin];
fHodoPosSigma=new Double_t [fMaxHodoScin];
fHodoNegSigma=new Double_t [fMaxHodoScin];
fHodoPosMinPh=new Double_t [fMaxHodoScin];
fHodoNegMinPh=new Double_t [fMaxHodoScin];
fHodoPosPhcCoeff=new Double_t [fMaxHodoScin];
fHodoNegPhcCoeff=new Double_t [fMaxHodoScin];
fHodoPosTimeOffset=new Double_t [fMaxHodoScin];
fHodoNegTimeOffset=new Double_t [fMaxHodoScin];
fHodoPosPedLimit=new Int_t [fMaxHodoScin];
fHodoNegPedLimit=new Int_t [fMaxHodoScin];
fHodoPosInvAdcOffset=new Double_t [fMaxHodoScin];
fHodoNegInvAdcOffset=new Double_t [fMaxHodoScin];
fHodoPosInvAdcLinear=new Double_t [fMaxHodoScin];
fHodoNegInvAdcLinear=new Double_t [fMaxHodoScin];
fHodoPosInvAdcAdc=new Double_t [fMaxHodoScin];
fHodoNegInvAdcAdc=new Double_t [fMaxHodoScin];
fxLoScin = new Int_t [fNHodoscopes];
fxHiScin = new Int_t [fNHodoscopes];
fyLoScin = new Int_t [fNHodoscopes];
fyHiScin = new Int_t [fNHodoscopes];
fHodoSlop = new Double_t [fNPlanes];
fTdcOffset = new Int_t [fNPlanes];
for(Int_t ip=0;ip<fNPlanes;ip++) { // Set a large default window
fTdcOffset[ip] = 0 ;
}
DBRequest list[]={
{"cosmicflag", &fCosmicFlag, kInt, 0, 1},
{"start_time_center", &fStartTimeCenter, kDouble},
{"start_time_slop", &fStartTimeSlop, kDouble},
{"scin_tdc_to_time", &fScinTdcToTime, kDouble},
{"scin_tdc_min", &fScinTdcMin, kDouble},
{"scin_tdc_max", &fScinTdcMax, kDouble},
{"tof_tolerance", &fTofTolerance, kDouble, 0, 1},
{"pathlength_central", &fPathLengthCentral, kDouble},
{"hodo_pos_sigma", &fHodoPosSigma[0], kDouble, fMaxHodoScin},
{"hodo_neg_sigma", &fHodoNegSigma[0], kDouble, fMaxHodoScin},
{"hodo_pos_ped_limit", &fHodoPosPedLimit[0], kInt, fMaxHodoScin},
{"hodo_neg_ped_limit", &fHodoNegPedLimit[0], kInt, fMaxHodoScin},
{"tofusinginvadc", &fTofUsingInvAdc, kInt, 0, 1},
{"xloscin", &fxLoScin[0], kInt, (UInt_t) fNHodoscopes},
{"xhiscin", &fxHiScin[0], kInt, (UInt_t) fNHodoscopes},
{"yloscin", &fyLoScin[0], kInt, (UInt_t) fNHodoscopes},
{"yhiscin", &fyHiScin[0], kInt, (UInt_t) fNHodoscopes},
{"track_eff_test_num_scin_planes", &fTrackEffTestNScinPlanes, kInt},
{"cer_npe", &fNCerNPE, kDouble, 0, 1},
{"normalized_energy_tot", &fNormETot, kDouble, 0, 1},
{"hodo_slop", fHodoSlop, kDouble, (UInt_t) fNPlanes},
{"debugprintscinraw", &fdebugprintscinraw, kInt, 0,1},
{"hodo_tdc_offset", fTdcOffset, kInt, (UInt_t) fNPlanes, 1},
{"dumptof", &fDumpTOF, kInt, 0, 1},
{"dumptof_filename", &fTOFDumpFile, kString, 0, 1},
{0}
};
// Defaults if not defined in parameter file
fdebugprintscinraw=0;
fDumpTOF = 0;
fTOFDumpFile="";
fTofUsingInvAdc = 1;
fTofTolerance = 3.0;
fNCerNPE = 2.0;
fNormETot = 0.7;
// Gets added to each reference time corrected raw TDC value
// to make sure valid range is all positive.
gHcParms->LoadParmValues((DBRequest*)&list,prefix);
if (fCosmicFlag==1) cout << " setup for cosmics in TOF"<< endl;
cout << " cosmic flag = " << fCosmicFlag << endl;
fDumpOut.open(fTOFDumpFile.c_str());
if(fDumpOut.is_open()) {
//fDumpOut << "Hodoscope Time of Flight calibration data" << endl;
} else {
fDumpTOF = 0;
cout << "WARNING: Unable to open TOF Dump file " << fTOFDumpFile << endl;
cout << "Data for TOF calibration not being written." << endl;
}
}
cout << " x1 lo = " << fxLoScin[0]
<< " x2 lo = " << fxLoScin[1]
<< " x1 hi = " << fxHiScin[0]
<< " x2 hi = " << fxHiScin[1]
cout << " y1 lo = " << fyLoScin[0]
<< " y2 lo = " << fyLoScin[1]
<< " y1 hi = " << fyHiScin[0]
<< " y2 hi = " << fyHiScin[1]
cout << "Hdososcope planes hits for trigger = " << fTrackEffTestNScinPlanes
<< " normalized energy min = " << fNormETot
<< " number of photo electrons = " << fNCerNPE
<< endl;
if (fTofUsingInvAdc) {
DBRequest list2[]={
{"hodo_vel_light", &fHodoVelLight[0], kDouble, fMaxHodoScin},
{"hodo_pos_invadc_offset",&fHodoPosInvAdcOffset[0],kDouble,fMaxHodoScin},
{"hodo_neg_invadc_offset",&fHodoNegInvAdcOffset[0],kDouble,fMaxHodoScin},
{"hodo_pos_invadc_linear",&fHodoPosInvAdcLinear[0],kDouble,fMaxHodoScin},
{"hodo_neg_invadc_linear",&fHodoNegInvAdcLinear[0],kDouble,fMaxHodoScin},
{"hodo_pos_invadc_adc",&fHodoPosInvAdcAdc[0],kDouble,fMaxHodoScin},
{"hodo_neg_invadc_adc",&fHodoNegInvAdcAdc[0],kDouble,fMaxHodoScin},
{0}
};
gHcParms->LoadParmValues((DBRequest*)&list2,prefix);
};
if (!fTofUsingInvAdc) {
DBRequest list3[]={
{"hodo_vel_light", &fHodoVelLight[0], kDouble, fMaxHodoScin},
{"hodo_pos_minph", &fHodoPosMinPh[0], kDouble, fMaxHodoScin},
{"hodo_neg_minph", &fHodoNegMinPh[0], kDouble, fMaxHodoScin},
{"hodo_pos_phc_coeff", &fHodoPosPhcCoeff[0], kDouble, fMaxHodoScin},
{"hodo_neg_phc_coeff", &fHodoNegPhcCoeff[0], kDouble, fMaxHodoScin},
{"hodo_pos_time_offset", &fHodoPosTimeOffset[0], kDouble, fMaxHodoScin},
{"hodo_neg_time_offset", &fHodoNegTimeOffset[0], kDouble, fMaxHodoScin},
{0}
};
gHcParms->LoadParmValues((DBRequest*)&list3,prefix);
};
Gabriel Niculescu
committed
if (fDebug >=1) {
cout <<"******* Testing Hodoscope Parameter Reading ***\n";
cout<<"StarTimeCenter = "<<fStartTimeCenter<<endl;
cout<<"StartTimeSlop = "<<fStartTimeSlop<<endl;
cout <<"ScintTdcToTime = "<<fScinTdcToTime<<endl;
cout <<"TdcMin = "<<fScinTdcMin<<" TdcMax = "<<fScinTdcMax<<endl;
cout <<"TofTolerance = "<<fTofTolerance<<endl;
cout <<"*** VelLight ***\n";
Gabriel Niculescu
committed
cout<<"Plane "<<i1<<endl;
Gabriel Niculescu
committed
cout<<fHodoVelLight[GetScinIndex(i1,i2)]<<" ";
Gabriel Niculescu
committed
cout <<endl;
Gabriel Niculescu
committed
cout <<endl<<endl;
// check fHodoPosPhcCoeff
/*
cout <<"fHodoPosPhcCoeff = ";
for (int i1=0;i1<fMaxHodoScin;i1++) {
cout<<this->GetHodoPosPhcCoeff(i1)<<" ";
}
cout<<endl;
*/
}
//
if ((fTofTolerance > 0.5) && (fTofTolerance < 10000.)) {
cout << "USING "<<fTofTolerance<<" NSEC WINDOW FOR FP NO_TRACK CALCULATIONS.\n";
}
else {
fTofTolerance= 3.0;
cout << "*** USING DEFAULT 3 NSEC WINDOW FOR FP NO_TRACK CALCULATIONS!! ***\n";
fIsInit = true;
return kOK;
}
//_____________________________________________________________________________
Int_t THcHodoscope::DefineVariables( EMode mode )
{
// Initialize global variables and lookup table for decoder
cout << "THcHodoscope::DefineVariables called " << GetName() << endl;
if( mode == kDefine && fIsSetup ) return kOK;
fIsSetup = ( mode == kDefine );
// Register variables in global list
// Move these into THcHallCSpectrometer using track fTracks
{"beta", "Beta including track info", "fBeta"},
{"betanotrack", "Beta from scintillator hits", "fBetaNoTrk"},
{"betachisqnotrack", "Chi square of beta from scintillator hits", "fBetaNoTrkChiSq"},
{"fpHitsTime", "Time at focal plane from all hits", "fFPTimeAll"},
{"goodstarttime", "Hodoscope Good Start Time (logical flag)", "fGoodStartTime"},
{"goodscinhit", "Hit in fid area", "fGoodScinHits"},
// {"goodscinhitx", "Hit in fid x range", "fGoodScinHitsX"},
{"scinshould", "Total scin Hits in fid area", "fScinShould"},
{"scindid", "Total scin Hits in fid area with a track", "fScinDid"},
return DefineVarsFromList( vars, mode );
// return kOK;
}
//_____________________________________________________________________________
THcHodoscope::~THcHodoscope()
{
// Destructor. Remove variables from global list.
delete [] fPlaneCenter;
delete [] fPlaneSpacing;
if( fIsSetup )
RemoveVariables();
if( fIsInit )
DeleteArrays();
if (fTrackProj) {
fTrackProj->Clear();
delete fTrackProj; fTrackProj = 0;
}
}
//_____________________________________________________________________________
void THcHodoscope::DeleteArrays()
{
// Delete member arrays. Used by destructor.
// for( k = 0; k < fNPlanes; k++){
// delete [] fScinHit[k];
// }
// delete [] fScinHit;
delete [] fxLoScin; fxLoScin = NULL;
delete [] fxHiScin; fxHiScin = NULL;
delete [] fHodoSlop; fHodoSlop = NULL;
delete [] fNPaddle; fNPaddle = NULL;
delete [] fHodoVelLight; fHodoVelLight = NULL;
delete [] fHodoPosSigma; fHodoPosSigma = NULL;
delete [] fHodoNegSigma; fHodoNegSigma = NULL;
delete [] fHodoPosMinPh; fHodoPosMinPh = NULL;
delete [] fHodoNegMinPh; fHodoNegMinPh = NULL;
delete [] fHodoPosPhcCoeff; fHodoPosPhcCoeff = NULL;
delete [] fHodoNegPhcCoeff; fHodoNegPhcCoeff = NULL;
delete [] fHodoPosTimeOffset; fHodoPosTimeOffset = NULL;
delete [] fHodoNegTimeOffset; fHodoNegTimeOffset = NULL;
delete [] fHodoPosPedLimit; fHodoPosPedLimit = NULL;
delete [] fHodoNegPedLimit; fHodoNegPedLimit = NULL;
delete [] fHodoPosInvAdcOffset; fHodoPosInvAdcOffset = NULL;
delete [] fHodoNegInvAdcOffset; fHodoNegInvAdcOffset = NULL;
delete [] fHodoPosInvAdcLinear; fHodoPosInvAdcLinear = NULL;
delete [] fHodoNegInvAdcLinear; fHodoNegInvAdcLinear = NULL;
delete [] fHodoPosInvAdcAdc; fHodoPosInvAdcAdc = NULL;
delete [] fGoodPlaneTime; fGoodPlaneTime = NULL;
delete [] fNPlaneTime; fNPlaneTime = NULL;
delete [] fSumPlaneTime; fSumPlaneTime = NULL;
delete [] fNScinHits; fNScinHits = NULL;
delete [] fTdcOffset; fTdcOffset = NULL;
}
//_____________________________________________________________________________
void THcHodoscope::ClearEvent()
/*! \brief Clears variables
*
* Called by THcHodoscope::Decode
*
*/
fStartTime = 0.0;
fFPTimeAll= -1000.;
fGoodStartTime = kFALSE;
fGoodScinHits = 0;
fScinShould = 0;
fScinDid = 0;
for(Int_t ip=0;ip<fNPlanes;ip++) {
fPlanes[ip]->Clear();
fPlaneCenter[ip]=0.;
fPlaneSpacing[ip]=0.;
for(UInt_t iPaddle=0;iPaddle<fNPaddle[ip]; ++iPaddle) {
fdEdX.clear();
fNScinHit.clear();
fNClust.clear();
fThreeScin.clear();
fGoodScinHitsX.clear();
}
//_____________________________________________________________________________
Int_t THcHodoscope::Decode( const THaEvData& evdata )
{
/*! \brief Decodes raw data and processes raw data into hits for each instance of THcScintillatorPlane
*
* - Calls THcHodoscope::ClearEvent
* - Reads raw data using THcHitList::DecodeToHitList
* - If one wants to subtract pedestals (assumed to be a set of data at beginning of run)
* + Must define "Pedestal_event" cut in the cuts definition file
* + For each "Pedestal_event" calls THcScintillatorPlane::AccumulatePedestals and returns
* + After First event which is not a "Pedestal_event" calls THcScintillatorPlane::CalculatePedestals
* - For each scintillator plane THcScintillatorPlane::ProcessHits
* - Calls THcHodoscope::EstimateFocalPlaneTime
ClearEvent();
// Get the Hall C style hitlist (fRawHitList) for this event
fNHits = DecodeToHitList(evdata);
//
// GN: print event number so we can cross-check with engine
// if (evdata.GetEvNum()>1000)
fCheckEvent = evdata.GetEvNum();
if(gHaCuts->Result("Pedestal_event")) {
Int_t nexthit = 0;
for(Int_t ip=0;ip<fNPlanes;ip++) {
nexthit = fPlanes[ip]->AccumulatePedestals(fRawHitList, nexthit);
}
fAnalyzePedestals = 1; // Analyze pedestals first normal events
return(0);
}
if(fAnalyzePedestals) {
for(Int_t ip=0;ip<fNPlanes;ip++) {
fPlanes[ip]->CalculatePedestals();
}
fAnalyzePedestals = 0; // Don't analyze pedestals next event
}
// Let each plane get its hits
Int_t nexthit = 0;
fNfptimes=0;
for(Int_t ip=0;ip<fNPlanes;ip++) {
fPlaneCenter[ip] = fPlanes[ip]->GetPosCenter(0) + fPlanes[ip]->GetPosOffset();
fPlaneSpacing[ip] = fPlanes[ip]->GetSpacing();
// nexthit = fPlanes[ip]->ProcessHits(fRawHitList, nexthit);
// GN: select only events that have reasonable TDC values to start with
// as per the Engine h_strip_scin.f
nexthit = fPlanes[ip]->ProcessHits(fRawHitList,nexthit);
EstimateFocalPlaneTime();
if (fdebugprintscinraw == 1) {
Mark Jones
committed
cout << " Event number = " << evdata.GetEvNum()<<endl;
for(UInt_t ihit = 0; ihit < fNRawHits ; ihit++) {
// THcRawHodoHit* hit = (THcRawHodoHit *) fRawHitList->At(ihit);
// cout << ihit << " : " << hit->fPlane << ":" << hit->fCounter << " : "
// << hit->fADC_pos << " " << hit->fADC_neg << " " << hit->fTDC_pos
// << " " << hit->fTDC_neg << endl;
}
cout << endl;
Stephen A. Wood
committed
return fNHits;
//_____________________________________________________________________________
void THcHodoscope::EstimateFocalPlaneTime( void )
{
/*! \brief Calculates the Drift Chamber start time and fBetaNoTrk (velocity determined without track info)
*
* - Called by THcHodoscope::Decode
* + loops through hits in each scintillator plane and fills histogram array, "timehist", with corrected times for positive
* + Determines the peak of "timehist"
Int_t timehist[200];
for (Int_t i=0;i<200;i++) {
timehist[i] = 0;
}
Int_t ihit=0;
Int_t nscinhits=0; // Total # hits with at least one good tdc
for(Int_t ip=0;ip<fNPlanes;ip++) {
Int_t nphits=fPlanes[ip]->GetNScinHits();
nscinhits += nphits;
TClonesArray* hodoHits = fPlanes[ip]->GetHits();
for(Int_t i=0;i<nphits;i++) {
THcHodoHit *hit = (THcHodoHit*)hodoHits->At(i);
if(hit->GetHasCorrectedTimes()) {
Double_t postime=hit->GetPosTOFCorrectedTime();
Double_t negtime=hit->GetNegTOFCorrectedTime();
for (Int_t k=0;k<200;k++) {
Double_t tmin=0.5*(k+1);
if ((postime> tmin) && (postime < tmin+fTofTolerance)) {
timehist[k]++;
}
if ((negtime> tmin) && (negtime < tmin+fTofTolerance)) {
timehist[k]++;
}
}
}
}
}
// Find the bin with most hits
Int_t binmax=0;
Int_t maxhit=0;
for(Int_t i=0;i<200;i++) {
if(timehist[i]>maxhit) {
binmax = i+1;
maxhit = timehist[i];
}
}
ihit = 0;
Double_t fpTimeSum = 0.0;
fNfptimes=0;
Mark Jones
committed
Int_t Ngood_hits_plane=0;
Double_t Plane_fptime_sum=0.0;
Bool_t goodplanetime[fNPlanes];
Bool_t twogoodtimes[nscinhits];
Int_t temp_planes=fNPlanes;
if (fSHMS) temp_planes=3;
for(Int_t ip=0;ip<temp_planes;ip++) {
goodplanetime[ip] = kFALSE;
Int_t nphits=fPlanes[ip]->GetNScinHits();
TClonesArray* hodoHits = fPlanes[ip]->GetHits();
Mark Jones
committed
Ngood_hits_plane=0;
Plane_fptime_sum=0.0;
for(Int_t i=0;i<nphits;i++) {
THcHodoHit *hit = (THcHodoHit*)hodoHits->At(i);
twogoodtimes[ihit] = kFALSE;
if(hit->GetHasCorrectedTimes()) {
Double_t tmin = 0.5*binmax;
Double_t postime=hit->GetPosTOFCorrectedTime();
Double_t negtime=hit->GetNegTOFCorrectedTime();
if ((postime>tmin) && (postime<tmin+fTofTolerance) &&
(negtime>tmin) && (negtime<tmin+fTofTolerance)) {
hit->SetTwoGoodTimes(kTRUE);
twogoodtimes[ihit] = kTRUE; // Both tubes fired
Int_t index=hit->GetPaddleNumber()-1; //
Double_t fptime;
if(fCosmicFlag==1) {
fptime = hit->GetScinCorrectedTime()
+ (fPlanes[ip]->GetZpos()+(index%2)*fPlanes[ip]->GetDzpos())
/ (29.979 * fBetaNominal);
}else{
fptime = hit->GetScinCorrectedTime()
- (fPlanes[ip]->GetZpos()+(index%2)*fPlanes[ip]->GetDzpos())
/ (29.979 * fBetaNominal);
if(TMath::Abs(fptime-fStartTimeCenter)<=fStartTimeSlop) {
Mark Jones
committed
Ngood_hits_plane++;
Plane_fptime_sum+=fptime;
fpTimeSum += fptime;
fNfptimes++;
goodplanetime[ip] = kTRUE;
} else {
hit->SetTwoGoodTimes(kFALSE);
}
ihit++;
fPlanes[ip]->SetFpTime(Plane_fptime_sum/float(Ngood_hits_plane));
fPlanes[ip]->SetNGoodHits(Ngood_hits_plane);
if(fNfptimes>0) {
fStartTime = fpTimeSum/fNfptimes;
fGoodStartTime=kTRUE;
} else {
fStartTime = fStartTimeCenter;
fGoodStartTime=kFALSE;
if((goodplanetime[0]||goodplanetime[1])
&&(goodplanetime[2]||goodplanetime[3])) {
Double_t sumW = 0.;
Double_t sumT = 0.;
Double_t sumZ = 0.;
Double_t sumZZ = 0.;
Double_t sumTZ = 0.;
Int_t ihhit = 0;
temp_planes=fNPlanes;
if (fSHMS) temp_planes=3;
for(Int_t ip=0;ip<temp_planes;ip++) {
Int_t nphits=fPlanes[ip]->GetNScinHits();
TClonesArray* hodoHits = fPlanes[ip]->GetHits();
for(Int_t i=0;i<nphits;i++) {
Int_t index=((THcHodoHit*)hodoHits->At(i))->GetPaddleNumber()-1;
if(twogoodtimes[ihhit]){
Double_t sigma = 0.5 * ( TMath::Sqrt( TMath::Power( fHodoPosSigma[GetScinIndex(ip,index)],2) +
TMath::Power( fHodoNegSigma[GetScinIndex(ip,index)],2) ) );
Double_t scinWeight = 1 / TMath::Power(sigma,2);
Double_t zPosition = fPlanes[ip]->GetZpos() + (index%2)*fPlanes[ip]->GetDzpos();
// cout << "hit = " << ihhit + 1 << " zpos = " << zPosition << " sigma = " << sigma << endl;
sumW += scinWeight;
sumT += scinWeight * ((THcHodoHit*)hodoHits->At(i))->GetScinCorrectedTime();
sumZ += scinWeight * zPosition;
sumZZ += scinWeight * ( zPosition * zPosition );
sumTZ += scinWeight * zPosition * ((THcHodoHit*)hodoHits->At(i))->GetScinCorrectedTime();
} // condition of good scin time
ihhit ++;
} // loop over hits of plane
} // loop over planes
Double_t tmp = sumW * sumZZ - sumZ * sumZ ;
Double_t t0 = ( sumT * sumZZ - sumZ * sumTZ ) / tmp ;
Double_t tmpDenom = sumW * sumTZ - sumZ * sumT;
if ( TMath::Abs( tmpDenom ) > ( 1 / 10000000000.0 ) ) {
fBetaNoTrkChiSq = 0.;
for (Int_t ip = 0; ip < temp_planes; ip++ ){ // Loop over planes
Int_t nphits=fPlanes[ip]->GetNScinHits();
TClonesArray* hodoHits = fPlanes[ip]->GetHits();
for(Int_t i=0;i<nphits;i++) {
Int_t index=((THcHodoHit*)hodoHits->At(i))->GetPaddleNumber()-1;
if(twogoodtimes[ihhit]) {
Double_t zPosition = fPlanes[ip]->GetZpos() + (index%2)*fPlanes[ip]->GetDzpos();
Double_t timeDif = ( ((THcHodoHit*)hodoHits->At(i))->GetScinCorrectedTime() - t0 );
Double_t sigma = 0.5 * ( TMath::Sqrt( TMath::Power( fHodoPosSigma[GetScinIndex(ip,index)],2) +
TMath::Power( fHodoNegSigma[GetScinIndex(ip,index)],2) ) );
fBetaNoTrkChiSq += ( ( zPosition / fBetaNoTrk - timeDif ) *
( zPosition / fBetaNoTrk - timeDif ) ) / ( sigma * sigma );
} // condition for good scin time
ihhit++;
} // loop over hits of a plane
} // loop over planes
fBetaNoTrk = fBetaNoTrk / 29.979; // velocity / c
} // condition for fTmpDenom
else {
fBetaNoTrk = 0.;
fBetaNoTrkChiSq = -2.;
} // else condition for fTmpDenom
fGoodEventTOFCalib=kFALSE;
if (!fSHMS&&goodplanetime[0]&&goodplanetime[1]&&goodplanetime[2]&&goodplanetime[3]&&fPlanes[0]->GetNGoodHits()==1&&fPlanes[1]->GetNGoodHits()==1&&fPlanes[2]->GetNGoodHits()==1&&fPlanes[3]->GetNGoodHits()==1) fGoodEventTOFCalib=kTRUE;
if (fSHMS&&goodplanetime[0]&&goodplanetime[1]&&goodplanetime[2]&&fPlanes[0]->GetNGoodHits()==1&&fPlanes[1]->GetNGoodHits()==1&&fPlanes[2]->GetNGoodHits()==1) fGoodEventTOFCalib=kTRUE;
//_____________________________________________________________________________
Int_t THcHodoscope::ApplyCorrections( void )
{
return(0);
}
//_____________________________________________________________________________
Double_t THcHodoscope::TimeWalkCorrection(const Int_t& paddle,
const ESide side)
{
return(0.0);
}
//_____________________________________________________________________________
Int_t THcHodoscope::CoarseProcess( TClonesArray& tracks )
Int_t ntracks = tracks.GetLast()+1; // Number of reconstructed tracks
Int_t timehist[200];
// -------------------------------------------------
// fDumpOut << " ntrack = " << ntracks << endl;
Int_t temp_planes=fNPlanes;
if (fSHMS) temp_planes=3;
// **MAIN LOOP: Loop over all tracks and get corrected time, tof, beta...
Double_t* nPmtHit = new Double_t [ntracks];
Double_t* timeAtFP = new Double_t [ntracks];
for ( Int_t itrack = 0; itrack < ntracks; itrack++ ) { // Line 133
nPmtHit[itrack]=0;
timeAtFP[itrack]=0;
THaTrack* theTrack = dynamic_cast<THaTrack*>( tracks.At(itrack) );
if (!theTrack) return -1;
for (Int_t ip = 0; ip < temp_planes; ip++ ){
fGoodPlaneTime[ip] = kFALSE;
fNPlaneTime[ip] = 0;
fSumPlaneTime[ip] = 0.;
std::vector<Double_t> dedx_temp;
fdEdX.push_back(dedx_temp); // Create array of dedx per hit
std::vector<std::vector<GoodFlags> > goodflagstmp1;
fGoodFlags.push_back(goodflagstmp1);
// timeAtFP[itrack] = 0.;
Double_t sumFPTime = 0.; // Line 138
fNScinHit.push_back(0);
//! Calculate all corrected hit times and histogram
//! This uses a copy of code below. Results are save in time_pos,neg
//! including the z-pos. correction assuming nominal value of betap
//! Code is currently hard-wired to look for a peak in the
//! range of 0 to 100 nsec, with a group of times that all
//! agree withing a time_tolerance of time_tolerance nsec. The normal
//! peak position appears to be around 35 nsec.
//! NOTE: if want to find particles with beta different than
//! reference particle, need to make sure this is big enough
//! to accomodate difference in TOF for other particles
//! Default value in case user hasnt defined something reasonable
for (Int_t j=0; j<200; j++) { timehist[j]=0; } // Line 176
// Loop over scintillator planes.
// In ENGINE, its loop over good scintillator hits.
fTOFCalc.clear(); // SAW - Can we
fTOFPInfo.clear(); // SAW - combine these two?
Int_t ihhit = 0; // Hit # overall
for(Int_t ip = 0; ip < temp_planes; ip++ ) {
std::vector<GoodFlags> goodflagstmp2;
fGoodFlags[itrack].push_back(goodflagstmp2);
TClonesArray* hodoHits = fPlanes[ip]->GetHits();
Double_t zPos = fPlanes[ip]->GetZpos();
Double_t dzPos = fPlanes[ip]->GetDzpos();
// first loop over hits with in a single plane
for (Int_t iphit = 0; iphit < fNScinHits[ip]; iphit++ ){
// iphit is hit # within a plane
THcHodoHit *hit = (THcHodoHit*)hodoHits->At(iphit);
fTOFPInfo.push_back(TOFPInfo());
// Can remove these as we will initialize in the constructor
// fTOFPInfo[ihhit].time_pos = -99.0;
// fTOFPInfo[ihhit].time_neg = -99.0;
// fTOFPInfo[ihhit].keep_pos = kFALSE;
// fTOFPInfo[ihhit].keep_neg = kFALSE;
fTOFPInfo[ihhit].scin_pos_time = 0.0;
fTOFPInfo[ihhit].scin_neg_time = 0.0;
fTOFPInfo[ihhit].hit = hit;
fTOFPInfo[ihhit].planeIndex = ip;
fTOFPInfo[ihhit].hitNumInPlane = iphit;
fTOFPInfo[ihhit].onTrack = kFALSE;
Int_t paddle = hit->GetPaddleNumber()-1;
Double_t zposition = zPos + (paddle%2)*dzPos;
Double_t xHitCoord = theTrack->GetX() + theTrack->GetTheta() *
( zposition ); // Line 183
Double_t yHitCoord = theTrack->GetY() + theTrack->GetPhi() *
( zposition ); // Line 184
if ( ( ip == 0 ) || ( ip == 2 ) ){ // !x plane. Line 185
scinTrnsCoord = xHitCoord;
scinLongCoord = yHitCoord;
else if ( ( ip == 1 ) || ( ip == 3 ) ){ // !y plane. Line 188
scinTrnsCoord = yHitCoord;
scinLongCoord = xHitCoord;
}
else { return -1; } // Line 195
fTOFPInfo[ihhit].scinTrnsCoord = scinTrnsCoord;
fTOFPInfo[ihhit].scinLongCoord = scinLongCoord;
Double_t scinCenter = fPlanes[ip]->GetPosCenter(paddle) + fPlanes[ip]->GetPosOffset();
// Index to access the 2d arrays of paddle/scintillator properties
Int_t fPIndex = GetScinIndex(ip,paddle);
( fPlanes[ip]->GetSize() * 0.5 + fPlanes[ip]->GetHodoSlop() ) ){ // Line 293
fTOFPInfo[ihhit].onTrack = kTRUE;
Double_t zcor = zposition/(29.979*fBetaNominal)*
TMath::Sqrt(1. + theTrack->GetTheta()*theTrack->GetTheta()
+ theTrack->GetPhi()*theTrack->GetPhi());
fTOFPInfo[ihhit].zcor = zcor;
if (fCosmicFlag) {
Double_t zcor = -zposition/(29.979*1.0)*
TMath::Sqrt(1. + theTrack->GetTheta()*theTrack->GetTheta()
+ theTrack->GetPhi()*theTrack->GetPhi());
fTOFPInfo[ihhit].zcor = zcor;