|
//____________________________________________________________________ // // Class BrFfsTofMatchingModule // // Module scanning the global track table and picking up those // matching valid tof hits. // Selected tracks saved into a table of BrTofTrackMatch named // "FFS Matching" // // Algorithm (Event method): // ------------------------- // // 1- pick up tables of H1 hits and FFS tracks // if one of them missing, return // // 2- Call method MatchTof(...) // a- Select tof hits by checking the TDCs and ADCs of both tubes // // b- loop over tracks and sub-loop over hits // i - evaluate the dX between track projection and hit position // ii - minimize this distance and store the "winner" hit // slat number in an array that has the dimension of // the track table. It might be that the minimal dX is // beyond the matching criterium, then, the hit id is -1 // // c- check "multiple hits" when more than 1 track share the same // hit. In that case, ignore these tracks (it doesn't remove // completely multiple hits since there might be one track // only pointing to such hits - the other tracks having not // been reconstructed) // // d- loop over remaining matchings and store result into table // of BrTofTrackMatch object (cf class implementation file) // // 3- If FFS and BFS are aligned and if there are some matchings in // H1 (a bit strict but ...), try the same game with H2 hits by // swimming tracks through D3 and D4 and projecting on H2 // //____________________________________________________________________ // // $Id: BrFfsTofMatchingModule.cxx,v 1.7 2002/04/16 11:48:17 ekman Exp $ // $Author: ekman $ // $Date: 2002/04/16 11:48:17 $ // $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov> // #ifndef BRAT_BrFfsTofMatchingModule #include "BrFfsTofMatchingModule.h" #endif #ifndef WIN32 #include <iostream> #include <iomanip> #include <fstream> #else #include <iomanip.h> #include <fstream.h> #include <iostream.h> #endif #ifndef BRAT_BrDataTable #include "BrDataTable.h" #endif #ifndef BRAT_BrDetectorTrack #include "BrDetectorTrack.h" #endif #ifndef BRAT_BrRunInfoManager #include "BrRunInfoManager.h" #endif #ifndef BRAT_BrParameterDbManager #include "BrParameterDbManager.h" #endif #ifndef BRAT_BrCalibrationManager #include "BrCalibrationManager.h" #endif #ifndef BRAT_BrGeometryDbManager #include "BrGeometryDbManager.h" #endif #ifndef BRAT_BrSpectrometerTracks #include "BrSpectrometerTracks.h" #endif #ifndef BRAT_BrLine3D #include "BrLine3D.h" #endif #ifndef BRAT_BrPlane3D #include "BrPlane3D.h" #endif //____________________________________________________________________ ClassImp(BrFfsTofMatchingModule); //____________________________________________________________________ BrFfsTofMatchingModule::BrFfsTofMatchingModule() : BrModule() { // Default constructor. DO NOT USE SetState(kSetup); SetDefaultParameters(); fCalib[kH1] = 0; fCalib[kH2] = 0; fT2Vol = 0; fTofVol[kH1] = 0; fTofVol[kH2] = 0; fTofPar[kH1] = 0; fTofPar[kH2] = 0; fD3Vol = 0; fD4Vol = 0; fUseH2 = kFALSE; } //____________________________________________________________________ BrFfsTofMatchingModule::BrFfsTofMatchingModule(const Char_t* name, const Char_t* title) : BrModule(name, title) { // Named Constructor SetState(kSetup); SetDefaultParameters(); fCalib[kH1] = 0; fCalib[kH2] = 0; fT2Vol = 0; fTofVol[kH1] = 0; fTofVol[kH2] = 0; fTofPar[kH1] = 0; fTofPar[kH2] = 0; fD3Vol = 0; fD4Vol = 0; fUseH2 = kFALSE; } //____________________________________________________________________ void BrFfsTofMatchingModule::SetDefaultParameters() { // Set default parameters SetMaxH1DX(); // default is +/- 0.5 cm SetMaxH2DX(); // default is +/- 0.5 cm SetT2H1XOffset(); // default is 0 SetT2H2XOffset(); // default is 0 SetD3PScale(); // default is 0.95 SetD4PScale(); // default is 0.9 SetUseH2(); // default is false } //____________________________________________________________________ void BrFfsTofMatchingModule::DefineHistograms() { // Define histograms. They are: // <fill in here> if (GetState() != kInit) { Stop("DefineHistograms", "Must be called after Init"); return; } TDirectory* saveDir = gDirectory; fHistDir[kH1] = gDirectory->mkdir("FFS_TofMatch"); fHistDir[kH2] = 0; fHistDir[kH1]->cd(); // list of histograms if (fNtuple) fNtMatch[kH1] = new TNtuple("H1Match", "Matched FFS Tracks - H1 Hits", "nTrack:yTrack:xTrack:tTime:bTime:slat:tEner:bEner" ":xHit:momentum"); fSlatCheck = new TH1F("slatCheck", "Pointed slat - hit Slat number for pairs matching", 10, -4.5, 5.5); fDxAll[kH1] = new TH1F("H1dXAll", "H1: #DeltaX FFS Track - H1 hit, " "all combinations", 400, -20, 20); fDxMatch[kH1] = new TH1F("H1dXMatch", "H1: #DeltaX FFS Track - H1 hit " "within matching condition", 400, -20, 20); fHistDir[kH1]->mkdir("Efficiency"); fHistDir[kH1]->cd("Efficiency"); fTopEff = new TH1F* [fTofPar[kH1]->GetNoSlats()]; fBotEff = new TH1F* [fTofPar[kH1]->GetNoSlats()]; for (Int_t s = 0; s < fTofPar[kH1]->GetNoSlats(); s++) { fTopEff[s] = new TH1F(Form("topEff%02d", s+1), Form("Efficiency, top tube %d", s+1), 30, -0.5, 2.5); fBotEff[s] = new TH1F(Form("botEff%02d", s+1), Form("Efficiency, bot tube %d", s+1), 30, -0.5, 2.5); } fTrkEff = new TH1F("trkEff", "Track-Hit eff - 0: at least a bad sig., 1: all signals OK", 20, -0.5, 1.5); fHistDir[kH1]->cd(); fMatchStat[kH1] = new TH1F("H1Stat", "Number of matchings per event in H1", 21, -0.5, 20.5); f2dXAll[kH1] = new TH2F("T2H1XAll", "X correlation T2 vs H1, all combinations", 51, -25.5, 25.5, 510, -25.5, 25.5); f2dXMatch[kH1] = new TH2F("T2H1XMatch", "X correlation T2 vs H1, matching pairs", 51, -25.5, 25.5, 510, -25.5, 25.5); if (fUseH2) { fHistDir[kH2] = fHistDir[kH1]->mkdir("H2"); fHistDir[kH2]->cd(); fDxAll[kH2] = new TH1F("H2dXAll", "H2: #DeltaX FFS Track - H2 hit, " "all combinations", 400, -100, 100); fDxMatch[kH2] = new TH1F("H2dXMatch", "H2: #DeltaX FFS Track - H2 hit " "within matching condition", 400, -100, 100); f2dXAll[kH2] = new TH2F("T2H2XAll", "X correlation T2 vs H2, all combinations", 51, -25.5, 25.5, 400, -100, 100); f2dXMatch[kH2] = new TH2F("T2H2XMatch", "X correlation T2 vs H2, matching pairs", 51, -25.5, 25.5, 400, -100, 100); fMatchStat[kH2] = new TH1F("H2Stat", "Number of matchings per event in H2", 21, -0.5, 20.5); if (fNtuple) fNtMatch[kH2] = new TNtuple("H2Match", "Matched FFS Tracks - H2 Hits", "nTrack:yTrack:xTrack:tTime:bTime:slat:tEner:bEner" ":xHit:momentum"); } gDirectory = saveDir; } //____________________________________________________________________ void BrFfsTofMatchingModule::Init() { // Job-level initialisation SetState(kInit); //--------------- // first check if fs is aligned if (fUseH2) { BrRunInfoManager* runMan = BrRunInfoManager::Instance(); if (!runMan) { Abort("Init", "Could not get an instance of BrRunInfoManager"); return; } const BrRunInfo* run = runMan->GetCurrentRun(); if (run->GetRunNo() == -1) { Abort("Init", "Run number is -1, check it out"); return; } if (run->GetFFSAngle() != run->GetBFSAngle()) { Warning("Init", "FFS and BFS are not aligned, disabling the use of H2"); fUseH2 = kFALSE; } } // Detector parameters and calibrations InitParams(); // geometry InitGeo(); // check tof pedestal calibration CheckTofCal(); // initialize all slat to good (let's be optimistic :) fValidSlat[kH1] = new Bool_t [fTofPar[kH1]->GetNoSlats()]; if (fUseH2) fValidSlat[kH2] = new Bool_t [fTofPar[kH2]->GetNoSlats()]; } //____________________________________________________________________ void BrFfsTofMatchingModule::InitParams() { // private // ------------------------------- // initialize detector parameters // ------------------------------- BrParameterDbManager* parDb = BrParameterDbManager::Instance(); BrCalibrationManager* calDb = BrCalibrationManager::Instance(); fCalib[kH1] = (BrTofCalibration*)calDb-> Register("BrTofCalibration", "TOF1"); fTofPar[kH1] = (BrDetectorParamsTof*)parDb-> GetDetectorParameters("BrDetectorParamsTof", "TOF1"); if (!fUseH2) return; fCalib[kH2] = (BrTofCalibration*)calDb-> Register("BrTofCalibration", "TOF2"); fTofPar[kH2] = (BrDetectorParamsTof*)parDb-> GetDetectorParameters("BrDetectorParamsTof", "TOF2"); } //____________________________________________________________________ void BrFfsTofMatchingModule::InitGeo() { // private // ------------------------------- // initialize geometry // ------------------------------- BrGeometryDbManager* geoDb = BrGeometryDbManager::Instance(); fT2Vol = (BrDetectorVolume*)geoDb-> GetDetectorVolume("BrDetectorVolume", "T2"); fTofVol[kH1] = (BrDetectorVolume*)geoDb-> GetDetectorVolume("BrDetectorVolume", "TOF1"); if (!fUseH2) return; fD3Vol = (BrMagnetVolume*)geoDb-> GetDetectorVolume("BrMagnetVolume", "D3"); fD4Vol = (BrMagnetVolume*)geoDb-> GetDetectorVolume("BrMagnetVolume", "D4"); fTofVol[kH2] = (BrDetectorVolume*)geoDb-> GetDetectorVolume("BrDetectorVolume", "TOF2"); } //____________________________________________________________________ void BrFfsTofMatchingModule::CheckTofCal() { // private // --------------------------- // check pedestal calibration // --------------------------- Int_t nSlats = fTofPar[kH1]->GetNoSlats(); // use pedestals and widths BrCalibration::EAccessMode mode = BrCalibration::kRead; // check needed calibration, if already loaded (ascii file) // by the tof pedestal module (cf brat guide) // if not, try to get them from DB if (!fCalib[kH1]->GetAccessMode("topPedestal") || !fCalib[kH1]->GetAccessMode("botPedestal") || !fCalib[kH1]->GetAccessMode("topPedestalWidth") || !fCalib[kH1]->GetAccessMode("botPedestalWidth")) { mode = BrCalibration::kRead; fCalib[kH1]->Use("topPedestal" ,mode, nSlats); fCalib[kH1]->Use("botPedestal", mode, nSlats); fCalib[kH1]->Use("topPedestalWidth", mode, nSlats); fCalib[kH1]->Use("botPedestalWidth", mode, nSlats); } // ---- if (!fUseH2) return; nSlats = fTofPar[kH2]->GetNoSlats(); // use pedestals and widths mode = BrCalibration::kRead; // check needed calibration, if already loaded (ascii file) // by the tof pedestal module (cf brat guide) // if not, try to get them from DB if (!fCalib[kH2]->GetAccessMode("topPedestal") || !fCalib[kH2]->GetAccessMode("botPedestal") || !fCalib[kH2]->GetAccessMode("topPedestalWidth") || !fCalib[kH2]->GetAccessMode("botPedestalWidth")) { mode = BrCalibration::kRead; fCalib[kH2]->Use("topPedestal" ,mode, nSlats); fCalib[kH2]->Use("botPedestal", mode, nSlats); fCalib[kH2]->Use("topPedestalWidth", mode, nSlats); fCalib[kH2]->Use("botPedestalWidth", mode, nSlats); } } //____________________________________________________________________ void BrFfsTofMatchingModule::Begin() { // Run-level initialisation SetState(kBegin); // --------------------------- // check pedestal revisions // --------------------------- // -------- check H1 if (!fCalib[kH1]->RevisionExists("*")) { Abort("Begin", "Could not find H1 pedestal revision. Aborting..."); return; } // check calibrated slats for(Int_t s = 0; s < fTofPar[kH1]->GetNoSlats(); s++) { fValidSlat[kH1][s] = kTRUE; if (!IsValid(fCalib[kH1]->GetBotPedestal(s+1)) || !IsValid(fCalib[kH1]->GetTopPedestal(s+1)) || !IsValid(fCalib[kH1]->GetBotPedestalWidth(s+1)) || !IsValid(fCalib[kH1]->GetTopPedestalWidth(s+1))) fValidSlat[kH1][s] = kFALSE; } // -------- check H2 if (!fUseH2) return; if (!fCalib[kH2]->RevisionExists("*")) { Warning("Begin", "Could not find H2 pedestal revision. " "Proceeding without..."); fUseH2 = kFALSE; } // check calibrated slats if (fUseH2) for(Int_t s = 0; s < fTofPar[kH2]->GetNoSlats(); s++) { fValidSlat[kH2][s] = kTRUE; if (!IsValid(fCalib[kH2]->GetBotPedestal(s+1)) || !IsValid(fCalib[kH2]->GetTopPedestal(s+1)) || !IsValid(fCalib[kH2]->GetBotPedestalWidth(s+1)) || !IsValid(fCalib[kH2]->GetTopPedestalWidth(s+1))) fValidSlat[kH2][s] = kFALSE; } } //____________________________________________________________________ void BrFfsTofMatchingModule::Event(BrEventNode* inNode, BrEventNode* outNode) { // Per event method SetState(kEvent); // ---------------------------------------------------------------- // Event method: // 1- pick up tables of H1 hits and FFS tracks // if one of them missing, return // // 2- Call method MatchTof(...) // a- Select tof hits by checking the TDCs and ADCs of both tubes // // b- loop over tracks and sub-loop over hits // i - evaluate the dX between track projection and hit position // ii - minimize this distance and store the "winner" hit // slat number in an array that has the dimension of // the track table. It might be that the minimal dX is // beyond the matching criterium, then, the hit id is -1 // // c- check "multiple hits" when more than 1 track share the same // hit. In that case, ignore these tracks (it doesn't remove // completely multiple hits since there might be one track // only pointing to such hits - the other tracks having not // been reconstructed) // // d- loop over remaining matchings and store result into table // of BrTofTrackMatch object (cf class implementation file) // // 3- If FFS and BFS are aligned and if there are some matchings in // H1 (a bit strict but ...), try the same game with H2 hits by // swimming tracks through D3 and D4 and projecting on H2 // ------------------------------------------------------------------- // ------- pick up ffs tracks BrDataTable* glbTracks = inNode->GetDataTable("FfsTracks"); if (!glbTracks) glbTracks = outNode->GetDataTable("FfsTracks"); if (!glbTracks) { if (Verbose() > 10) Warning("Event", "No FFS tracks in this event"); return; } // ------- pick up tof digits BrDataTable* h1Hits = inNode->GetDataTable("DigTof TOF1"); if (!h1Hits) h1Hits = outNode->GetDataTable("DigTof TOF1"); if (!h1Hits) { if (Verbose() > 10) Warning("Event", "No TOF1 hits in this event"); return; } // check efficiency: if (HistOn()) CheckEfficiency(glbTracks, h1Hits); // ------- call matching method, check if H2 can be used Int_t ntrk = glbTracks->GetEntries(); Int_t matchId1[ntrk]; // array of track indexes matching tof1 hits BrTofDig* accHit1[h1Hits->GetEntries()]; Int_t nH1Match = MatchTof(kH1, h1Hits, glbTracks, matchId1, accHit1); if (!nH1Match) return; // ------- pick up tof digits Bool_t h2Ok = fUseH2; BrDataTable* h2Hits = inNode->GetDataTable("DigTof TOF2"); if (!h2Hits && fUseH2) h2Hits = outNode->GetDataTable("DigTof TOF2"); if (!h2Hits && fUseH2) { if (Verbose() > 10) Warning("Event", "No TOF2 hits in this event"); h2Ok = kFALSE; } // try matching with H2 Int_t nh = 0; Int_t nt = 0; if (h2Ok) { nt = ntrk; nh = h2Hits->GetEntries(); } BrTofDig* accHit2[nh]; Int_t matchId2[nt]; // array of track indexes matching tof1 hits Int_t nH2Match = 0; if (h2Ok) nH2Match = MatchTof(kH2, h2Hits, glbTracks, matchId2, accHit2, matchId1); // --------------------------------------- // prepare output BrDataTable* matchTab = new BrDataTable("FFS Matching"); // ------- loop over pairs for (Int_t t = 0; t < ntrk; t++) { if (matchId1[t] == -1) continue; BrTofDig* hit1 = accHit1[matchId1[t]]; Int_t slat1 = hit1->GetSlatno(); BrFfsTrack* trk = (BrFfsTrack*)glbTracks->At(t); // H2 stuff BrTofDig* hit2 = 0; Int_t slat2 = 0; if (nH2Match) if (matchId2[t] > -1) { hit2 = accHit2[matchId2[t]]; slat2 = hit2->GetSlatno(); } if (HistOn()) fSlatCheck->Fill(trk->GetPointedSlat() - slat1); // create object if (trk->GetPointedSlat() - slat1 == 0) { BrTofTrackMatch* tt = new BrTofTrackMatch(); tt->SetMatching(trk->GetTrackId(), slat1, 0, slat2); // 0 is for H1 panel id matchTab->Add(tt); } // ------------------------------------------------- if (!HistOn()) continue; Double_t ttdc = hit1->GetTdcUp(); Double_t tadc = hit1->GetAdcUp() - fCalib[kH1]->GetTopPedestal(slat1); Double_t btdc = hit1->GetTdcDown(); Double_t badc = hit1->GetAdcDown() - fCalib[kH1]->GetBotPedestal(slat1); Double_t xHit = fTofPar[kH1]->GetSlatPos(slat1)[0] + fXOffset[kH1]; BrVector3D proj = trk->GetProjOnTof(); fDxMatch[kH1]->Fill(proj(0) - xHit); f2dXMatch[kH1]->Fill(xHit, proj(0)); if (fNtuple) fNtMatch[kH1]->Fill(nH1Match, proj(1), proj(0), ttdc, btdc, slat1, tadc, badc, xHit, trk->GetMomentum()); if (hit2) { ttdc = hit2->GetTdcUp(); tadc = hit2->GetAdcUp() - fCalib[kH2]->GetTopPedestal(slat2); btdc = hit2->GetTdcDown(); badc = hit2->GetAdcDown() - fCalib[kH2]->GetBotPedestal(slat2); xHit = fTofPar[kH2]->GetSlatPos(slat2)[0] + fXOffset[kH2]; proj = SwimToH2(trk->GetBackTrack(), trk->GetMomentum()); fDxMatch[kH2]->Fill(proj(0) - xHit); f2dXMatch[kH2]->Fill(xHit, proj(0)); if (fNtuple) fNtMatch[kH2]->Fill(nH2Match, proj(1), proj(0), ttdc, btdc, slat2, tadc, badc, xHit, trk->GetMomentum()); } } // check match table if (matchTab->GetEntries()) outNode->AddDataTable(matchTab); else delete matchTab; } //____________________________________________________________________ void BrFfsTofMatchingModule::CheckEfficiency(BrDataTable* trks, BrDataTable* hits) { // private // ---------------------------------------------------------------- // check that when a track points a slat, there's a valid signal // in top and bot tubes // ---------------------------------------------------------------- for (Int_t t = 0; t < trks->GetEntries(); t++) { BrVector3D proj = ((BrFfsTrack*)trks->At(t))->GetProjOnTof(); BrVector3D off(fXOffset[kH1], 0, 0); proj -= off; Int_t psl = fTofPar[kH1]->GetSlatNo(proj(0)); if (psl < 1 || psl > fTofPar[kH1]->GetNoSlats()) continue; Bool_t good = kFALSE; for (Int_t h = 0; h < hits->GetEntries(); h++) { BrTofDig* hit = (BrTofDig*)hits->At(h); if (psl != hit->GetSlatno()) continue; Int_t slat = hit->GetSlatno(); // pedestals Float_t tPed = fCalib[kH1]->GetTopPedestal(slat); Float_t bPed = fCalib[kH1]->GetBotPedestal(slat); Float_t tWid = fCalib[kH1]->GetTopPedestalWidth(slat); Float_t bWid = fCalib[kH1]->GetBotPedestalWidth(slat); // hit info Double_t tadc = hit->GetAdcUp() - tPed; Double_t badc = hit->GetAdcDown() - bPed; Double_t ttdc = hit->GetTdcUp(); Double_t btdc = hit->GetTdcDown(); // top adc and tdc checks if (tadc > fNoPedWidth*tWid) { fTopEff[slat-1]->Fill(1.); if (ttdc < 4000 && ttdc > 10) { fTopEff[slat-1]->Fill(2.); good = kTRUE; } } else fTopEff[slat-1]->Fill(0.); // bot adc and tdc checks if (badc > fNoPedWidth*tWid) { fBotEff[slat-1]->Fill(1.); if (btdc < 4000 && btdc > 10) fBotEff[slat-1]->Fill(2.); else good = kFALSE; } else { fBotEff[slat-1]->Fill(0.); good = kFALSE; } } if (good) fTrkEff->Fill(1.); else fTrkEff->Fill(0.); } } //____________________________________________________________________ Int_t BrFfsTofMatchingModule::MatchTof(ETofId id, BrDataTable* hitTab, BrDataTable* trkTab, Int_t* matchId, BrTofDig** hits, Int_t* h1Id) { Int_t ntrk = trkTab->GetEntries(); // number of global tracks in event Int_t nhit = hitTab->GetEntries(); // number of tof hits in event // ------------------------------------------------------- // ------ 1st: selection of tof hits ------ // array of accepted hits for (Int_t i = 0; i < nhit; i++) hits[i] = 0; Int_t nacc = SelectHits(id, hitTab, hits); if (Verbose() > 30) cout << " Number of selected hits in TOF" << id+1 << " : " << nacc << endl; if (!nacc) return 0; // if H2, we have to set the D3 and D4 fields if (id == kH2) if (!SetMagnetFields()) return 0; // ------ 2nd match tracks to hits ------ // ---------------------------------------------------------------- // loop over tracks and tof hits to find best hit for given track // ---------------------------------------------------------------- Bool_t goOn = kFALSE; for (Int_t t = 0; t < ntrk; t++) { // loop over tracks matchId[t] = -1; // initialize matching index // if H2, check if this track already matches with H1: if (id == kH2 && h1Id) if (h1Id[t] == -1) continue; BrFfsTrack* trk = (BrFfsTrack*)trkTab->At(t); Double_t dX = 999999.; // dX between hit and back track proj BrVector3D proj = trk->GetProjOnTof(); // local to Tof plane // check that the track status is OK if (trk->GetStatus() != BrMatchedTrack::kOk) continue; if (id == kH2) { proj = SwimToH2(trk->GetBackTrack(), trk->GetMomentum()); if (proj(0) == -999) continue; } // ---- loop over hits for(Int_t h = 0; h < nacc; h++) { // determine x pos. of tof hit (middle of a slat) Int_t slat = hits[h]->GetSlatno(); BrVector3D hitPos = fTofPar[id]->GetSlatPos(slat); // correct for the offset in X between T2 and H1 hitPos(0) += fXOffset[id]; // have now X hit position on tof plane // gonna minimize the difference with track proj X cut Double_t diffX = proj(0) - hitPos(0); if (HistOn()) { fDxAll[id]->Fill(diffX); f2dXAll[id]->Fill(hitPos(0), proj(0)); } // ignore cases where we are too far off if (TMath::Abs(diffX) > fMaxDX[id]) continue; if (TMath::Abs(diffX) < TMath::Abs(dX)) { dX = diffX; matchId[t] = h; } } if (matchId[t] > -1) goOn = kTRUE; if (Verbose() > 10) { cout << "FFS track " << t << " : dX = " << dX; if (matchId[t] > -1) cout << " with TOF" << id+1 << " hit at slat " << hits[matchId[t]]->GetSlatno() << endl; else cout << " ---> no hits in TOF" << id+1 << " within matching conditions " << endl; } } if (!goOn) return 0; // ------- check multiple hits for (Int_t t = 0; t < ntrk - 1; t++) for (Int_t i = t + 1; i < ntrk; i++) { if (matchId[i] == -1 || matchId[t] == -1) continue; if (matchId[i] == matchId[t]) { // multiple hit (?) matchId[i] = -1; matchId[t] = -1; } } // ------- evaluate the number of matchings and return it Int_t nMatch = 0; for (Int_t t = 0; t < ntrk; t++) { if (matchId[t] == -1) continue; nMatch++; } if (HistOn()) fMatchStat[id]->Fill(nMatch); if (Verbose() > 10) cout << " Found " << nMatch << " FFS tracks matching TOF" << id+1 << " hits" << " within a matching limit of " << fMaxDX[id] << " cm " << endl; return nMatch; } //___________________________________________________________________ Bool_t BrFfsTofMatchingModule::SetMagnetFields() { // private // ---------------------------------------------------- // get field information from database // update info each time a field is needed (event loop) // ---------------------------------------------------- BrRunInfoManager* runMan = BrRunInfoManager::Instance(); if (!runMan) { Warning("SetMagnetFields", "Could not get an instance of BrRunInfoManager"); return kFALSE; } const BrRunInfo* run = runMan->GetCurrentRun(); if (run->GetRunNo() == -1) { Warning("SetMagnetFields", "Run number is -1, check it out"); return kFALSE; } // get magnet settings fD3Vol->SetCurrent(run->GetD3Set(), run->GetD3Pol()); fD4Vol->SetCurrent(run->GetD4Set(), run->GetD4Pol()); if (Verbose() > 10) cout << " D3 field: " << fD3Vol->GetField() << " D4 field: " << fD4Vol->GetField() << endl; return kTRUE; } //_____________________________________________________________ Int_t BrFfsTofMatchingModule::SelectHits(ETofId id, BrDataTable* hits, BrTofDig** sel) { // private // -------------------------------- // tof hit selection: // 1- check ADCs and TDCs // 2- Fill array of good hits // 3- return their number // -------------------------------- Int_t hitOk = 0; for(Int_t h = 0; h < hits->GetEntries(); h++) { BrTofDig* hit = (BrTofDig*)hits->At(h); Int_t slat = hit->GetSlatno(); // check if slat ok if (!fValidSlat[id][slat-1]) continue; // pedestals Float_t tPed = fCalib[id]->GetTopPedestal(slat); Float_t bPed = fCalib[id]->GetBotPedestal(slat); Float_t tWid = fCalib[id]->GetTopPedestalWidth(slat); Float_t bWid = fCalib[id]->GetBotPedestalWidth(slat); // hit info Double_t tadc = hit->GetAdcUp() - tPed; Double_t badc = hit->GetAdcDown() - bPed; Double_t ttdc = hit->GetTdcUp(); Double_t btdc = hit->GetTdcDown(); // adc and tdc checks if (tadc < fNoPedWidth*tWid || badc < fNoPedWidth*bWid || ttdc > fMaxTdc || btdc > fMaxTdc || ttdc < fMinTdc || btdc < fMinTdc) continue; sel[hitOk] = hit; hitOk++; } return hitOk; } //_____________________________________________________________ BrVector3D BrFfsTofMatchingModule::SwimToH2(BrDetectorTrack* trk, Double_t p) { // private // swim FFS back track to H2 through D3 and D4 BrVector3D bad(-999, -999, -999); BrLine3D gT2Line(fT2Vol->LocalToGlobal(trk->GetTrackLine())); BrLine3D t2LineN(fD3Vol->GlobalToLocal(gT2Line)); BrLine3D t2LineX(fD3Vol->SwimForward(t2LineN, fD3PScale*p)); if (!fD3Vol->GetSwimStatus(t2LineN, t2LineX, 0.01)) return bad; gT2Line = fD3Vol->LocalToGlobal(t2LineX); // D4 t2LineN = fD4Vol->GlobalToLocal(gT2Line); t2LineX = fD4Vol->SwimForward(t2LineN, fD4PScale*p); if (!fD4Vol->GetSwimStatus(t2LineN, t2LineX, 0.01)) return bad; gT2Line = fD4Vol->LocalToGlobal(t2LineX); // H2 BrLine3D h2Line(fTofVol[kH2]->GlobalToLocal(gT2Line)); BrPlane3D h2Plane(0, 0, 0, 0, 1, 0, 1, 0, 0); return h2Plane.GetIntersectionWithLine(h2Line); } //____________________________________________________________________ void BrFfsTofMatchingModule::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: Djamel Ouerdane" << endl << " Last Modifications: " << endl << " $Author: ekman $" << endl << " $Date: 2002/04/16 11:48:17 $" << endl << " $Revision: 1.7 $ " << endl << endl << "-------------------------------------------------" << endl; } //____________________________________________________________________ // // $Log: BrFfsTofMatchingModule.cxx,v $ // Revision 1.7 2002/04/16 11:48:17 ekman // Check that the status of the ffs track is OK before matching with tof hits. // // Revision 1.6 2002/04/09 02:09:31 ouerdane // removed name restriction from constructor (useless) // // Revision 1.5 2001/11/05 07:11:39 ouerdane // changed FFS to Ffs, BFS to Bfs, MRS to Mrs. Fixed bugs in BFS module // // Revision 1.4 2001/10/19 15:54:37 ouerdane // added histos for TOF1 efficiency, fixed a couple of things in the slat pos // // Revision 1.3 2001/10/08 11:29:57 cholm // Changed to use new DB access classes // // Revision 1.2 2001/10/02 20:20:39 ouerdane // remove member fNoPanels and setter method since it is now part of the // tof detector parameter class. // Corrected errors due to fast rewriting of the modules. should be ok now. // // Revision 1.1 2001/10/02 02:05:25 ouerdane // Added modules combining and saving Tracks and Tof hits // // |
||||||
This page automatically generated by script docBrat by Christian Holm |
Copyright ; 2002 BRAHMS Collaboration
<brahmlib@rcf.rhic.bnl.gov>
|