|
//____________________________________________________________________ // // Creates new BrDataTable of reduced hits in any given Drift Chamber. // It requires RawDC data in the inNode and the created table is // inserted into the outNode as a BrDataTable ( "RdoDC T3", "RdoDC T4" // or "RdoDC T3" ) of BrDetectorHits. // Using this tables eliminates the necessity of using both BrDcCalibration // and BrDriverDC since it contains 'z' position (Pos[2]) and distance // from center ( xwire -/+ xdrift -> Pos[0]/Pos[1] -> ambiguity included ) // of all the hits. // //____________________________________________________________________ // // $Id: BrDcRdoModule.cxx,v 1.2 2002/05/09 17:58:38 ouerdane Exp $ // $Author: ouerdane $ // $Date: 2002/05/09 17:58:38 $ // $Copyright: (C) 2002 BRAHMS Collaboration <brahmlib@rhic.bnl.gov> // //----------------------------------------------------------------// // #ifndef BRAT_BrDcRdoModule #include "BrDcRdoModule.h" #endif #ifndef BRAT_BrDcDig #include "BrDcDig.h" #endif #ifndef BRAT_BrCalibrationManager #include "BrCalibrationManager.h" #endif #ifndef BRAT_BrDataTable #include "BrDataTable.h" #endif #ifndef BRAT_BrHit #include "BrHit.h" #endif #ifndef BRAT_BrDetectorHit #include "BrDetectorHit.h" #endif #ifndef ROOT_TString #include "TString.h" #endif #ifndef ROOT_TDirectory #include "TDirectory.h" #endif #ifndef __IOSTREAM__ #include <iostream> #endif #ifndef __IOMANIP__ #include <iomanip> #endif #ifndef WIN32 #include <iostream> #include <iomanip> #else #include <iostream.h> #include <iomanip.h> #endif //____________________________________________________________________ ClassImp(BrDcRdoModule); //____________________________________________________________________ BrDcRdoModule::BrDcRdoModule() : BrModule() { // Default constructor. DO NOT USE SetState(kSetup); fDCDriver = 0; fMinTdcPP = 0; fPromptDelayPP = 0; fMaxTdcPP = 0; fMinWidth = 0; } //____________________________________________________________________ BrDcRdoModule::BrDcRdoModule(const Char_t* name, const Char_t* title) : BrModule(name, title) { // ------------------------------------ // Normal constructor. // This Name & title is informational. // ------------------------------------ SetState(kSetup); fDCDriver = 0; if (!fName.CompareTo("T3", TString::kIgnoreCase)) { fMinTdcPP = 3970; fPromptDelayPP = 0; fMaxTdcPP = 4420; fMinWidth = 16; } if (!fName.CompareTo("T4", TString::kIgnoreCase)) { fMinTdcPP = 2250; fPromptDelayPP = 3520; fMaxTdcPP = 4470; } if (!fName.CompareTo("T5", TString::kIgnoreCase)) { fMinTdcPP = 2250; fPromptDelayPP = 3320; fMaxTdcPP = 4370; } } //____________________________________________________________________ void BrDcRdoModule::DefineHistograms() { // Define histograms. They are: // Tdc of hits // Multiplicity of hits // Drift distance per plane/view if (GetState() != kInit) { Stop("DefineHistograms", "Must be called after Init"); return; } TDirectory* saveDir = gDirectory; TDirectory* histDir = gDirectory->mkdir(Form("Rdo%s", GetName())); histDir->cd(); Int_t index[10]; Int_t detectorId = 0; if(!fName.CompareTo("T3", TString::kIgnoreCase)) { detectorId = 3; index[0] = 1; index[1] = 2; index[2] = 3; index[3] = 1; index[4] = 2; index[5] = 3; index[6] = 1; index[7] = 2; index[8] = 1; index[9] = 2; } else if(!fName.CompareTo("T4", TString::kIgnoreCase)) { detectorId = 4; index[0] = 1; index[1] = 2; index[2] = 1; index[3] = 2; index[4] = 1; index[5] = 2; index[6] = 1; index[7] = 2; } else if(!fName.CompareTo("T5", TString::kIgnoreCase)) { detectorId = 5; index[0] = 1; index[1] = 2; index[2] = 1; index[3] = 2; index[4] = 1; index[5] = 2; index[6] = 1; index[7] = 2; } else { Abort("DefineHistograms", "unknown drift chamber '%s'", GetName()); return; } TString histName(Form("t%dHitsAll", detectorId)); TString histTitle(Form("Multiplicity in T%d",detectorId)); fHitsMultiplicity = new TH1F (histName.Data(), histTitle.Data(), 301, -0.5, 300.5 ); histName = Form("t%dTime",detectorId); histTitle = Form("Tdc of T%d hits",detectorId); fTime = new TH1F (histName.Data(), histTitle.Data(), 400, 1500, 5500); Char_t view[4]; view[0]='X'; view[1]='Y'; view[2]='U'; view[3]='V'; // Drift distance in different planes for (Int_t i = 0 ; i < 3 ; i++) { for (Int_t j = 0 ; j < fDCDriver->GetMaxPlaneNum() ; j++) { histName = Form("t%d%d%c%dHitsPos", detectorId, i+1, view[fDCDriver->GetView(j)], index[j]); histTitle = Form("T%d, module %d, plane %c%d, position of hits", detectorId, i+1, view[fDCDriver->GetView(j)], index[j]); fHitsPositions[i][j] = new TH1F(histName.Data(), histTitle.Data(), 100, -0.5, 1.5 ); } } // Drift distance in different views for (Int_t k = 0 ; k < 4 ; k++) { histName = Form("t%d%cHitsPos",detectorId, view[k]); histTitle = Form("Position in %c of T%d hits",view[k],detectorId); fHitsPos[k] = new TH1F (histName.Data(), histTitle.Data(), 100, -0.5, 1.5 ); } gDirectory = saveDir; } //____________________________________________________________________ void BrDcRdoModule::Init() { // Job-level initialisation SetState(kInit); // ----------------------------------------- // Initialze the matching module and setup BrCalibrationManager* calMan = BrCalibrationManager::Instance(); fCalibration = (BrDcCalibration*) calMan->Register("BrDcCalibration", GetName()); if (!fCalibration) { Abort("Init", "Couldn't get a pointer to the DC Calibration parameter"); return; } Debug(5, "Init", "using SQL"); fCalibration->Use("tdcOffset", BrCalibration::kRead, 1); fCalibration->Use("driftTime", BrCalibration::kRead, 1); // DO: Driver setting should be done at Begin time only if(!fDCDriver) fDCDriver = new BrDriverDC(GetName(),GetName()); } //____________________________________________________________________ void BrDcRdoModule::Begin() { // Run-level initialisation SetState(kBegin); if (!fCalibration->RevisionExists("*")) { Abort("Begin", "Could not find any revision of tdcOffset for %s", GetName()); return; } if (!fCalibration->ValidCalibration()) { Abort("Begin", "One calibration is missing, either drift time or tdc offset, " " check revisions with the calib tool on the DAQ homepage", GetName()); return; } fGlobalTdcOffset = fCalibration->GetTdcOffset(); fCalibChan = fCalibration->GetDriftTime(); fMinTdc = fMinTdcPP + fGlobalTdcOffset; fPromptDelay = fPromptDelayPP + fGlobalTdcOffset; fMaxTdc = fMaxTdcPP + fGlobalTdcOffset; if (!fOffsetFile.IsNull()) fDCDriver->ReadOffsets(fOffsetFile, fGlobalTdcOffset); if (!fCalibFile.IsNull()) fDCDriver->ReadCalibData(fCalibFile, fCalibChan); fDCDriver->SetMinTdc(fMinTdc); fDCDriver->SetMaxTdc(fMaxTdc); Debug(5, "Begin", "driver settings: %d %d %d", fMinTdc, fMaxTdc, fMinWidth); } //____________________________________________________________________ void BrDcRdoModule::Event(BrEventNode* inNode, BrEventNode* outNode) { // Per event method SetState(kEvent); BrDataTable* rdoTable = new BrDataTable(Form("RdoDC %s", GetName()), "Rdo data"); //BrDataTable* newRdoTable = //new BrDataTable(Form("NewRdoDC %s", GetName()), "Rdo data"); outNode->AddDataTable(rdoTable); //outNode->AddDataTable(newRdoTable); BrDataTable *rawTable= 0; if (!(rawTable = inNode->GetDataTable(Form("DigDC %s", GetName())))) { Info(10, "Event", "couldn't get the inputtable 'DigDc %s'", GetName()); return; } Int_t nHits = rawTable->GetEntries(); if(nHits > 1999) { Info(10, "Event", "found %d (>1999) hits, givin up", nHits); return; } Int_t cleanHits = 0; for(Int_t i = 0; i < nHits; i++) { BrDcDig* digDc = (BrDcDig*)rawTable->At(i); Int_t module = digDc->GetSubModule() - 1; Int_t plane = digDc->GetPlane() - 1; Int_t wire = digDc->GetWire() - 1; Int_t time = digDc->GetTime(); Int_t width = digDc->GetWidth(); if(!fDCDriver->HitIsOK(time, width)) { Debug(10, "Event", "hit at %d, %d not vaild", time, width); continue; } // We're working on either T4 or T5 // In T4 and T5 there is rather complicated spectrum in time, // and that is why it looks so strange, but is OK if (fName.CompareTo("T3", TString::kIgnoreCase)) { if (!fPromptDelay) wire = wire % 100; else { if (wire < 100) { if (time <= fPromptDelay) // Is within the delay length, so continue continue; } else { if (time > fPromptDelay) // time is outside of the delay, so continue continue; else wire = wire % 100; } } } // There can be two ways of writing RdoDc data, as far // the one with BrDetectorHit is used // BrHit* newHit = new BrHit; BrDetectorHit* hit = new BrDetectorHit(); fDCDriver->AtPlane(plane); Float_t pos[] = { (fDCDriver->GetWirePosition(wire) - fDCDriver->GetDriftDist(module, plane, wire, time)), (fDCDriver->GetWirePosition(wire) + fDCDriver->GetDriftDist(module, plane, wire, time)), fDCDriver->GetPlanePosition(module, plane) }; hit->SetPos(pos); hit->SetID(wire); hit->SetImod(100 * module + plane); cleanHits++; //newHit->SetX(pos[0]); //newHit->SetY(pos[1]); //newHit->SetZ(pos[2]); //newHit->SetUniqueID((100 * module + plane) * 1000 + wire); //newRdoTable->Add(newHit); rdoTable->Add(hit); if (!HistOn()) continue; fTime->Fill(time); fHitsPositions[module][plane] ->Fill(fDCDriver->GetDriftDist(module, plane, wire, time)); fHitsPos[fDCDriver->GetView(plane)] ->Fill(fDCDriver->GetDriftDist(module, plane, wire, time)); } if (HistOn()) fHitsMultiplicity->Fill(cleanHits); } //____________________________________________________________________ void BrDcRdoModule::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: Radek Karabowicz" << endl << " Last Modifications: " << endl << " $Author: ouerdane $" << endl << " $Date: 2002/05/09 17:58:38 $" << endl << " $Revision: 1.2 $ " << endl << endl << "-------------------------------------------------" << endl; } //____________________________________________________________________ // // $Log: BrDcRdoModule.cxx,v $ // Revision 1.2 2002/05/09 17:58:38 ouerdane // updated init and begin to get the driftTime from the DB. Fixed wrong code structure in init // // Revision 1.1 2002/04/11 18:11:15 radziu // Module for creating the reduced data for the Drift Chambers. // A new BrDataTable is created ( named "RdoDC Tx", x=3,4or5 ), // consisting of BrDetectorHits. // // |
||||||
This page automatically generated by script docBrat by Christian Holm |
Copyright ; 2002 BRAHMS Collaboration
<brahmlib@rcf.rhic.bnl.gov>
|