|
//____________________________________________________________________ // // Rdo module for the pp TMrsF trigger counters. // // // The bb information is only used in case of creating an Ntuple. // This Ntuple is very useful to check internal resolution, timing offfsets // and slewing corrections. // // // TODO: // Add calibrations for Position and timing offsets. // Add calibrations for eff speed of light. // //____________________________________________________________________ // // $Id: BrTMrsFRdoModule.cxx,v 1.2 2002/08/06 19:10:59 hagel Exp $ // $Author: hagel $ // $Date: 2002/08/06 19:10:59 $ // $Copyright: (C) 2002 BRAHMS Collaboration <brahmlib@rhic.bnl.gov> // #ifndef BRAT_BrTMrsFRdoModule #include "BrTMrsFRdoModule.h" #endif #ifndef ROOT_TDirectory #include "TDirectory.h" #endif #include <BrIostream.h> #ifndef ROOT_TH1 #include "TH1.h" #endif #ifndef ROOT_TH2 #include "TH2.h" #endif #ifndef ROOT_TNtuple #include "TNtuple.h" #endif #ifndef BRAT_BrEventNode #include "BrEventNode.h" #endif #ifndef BRAT_BrTMrsFRdo #include "BrTMrsFRdo.h" #endif #ifndef BRAT_BrTofDig #include "BrTofDig.h" #endif #ifndef BRAT_BrBbVertex #include "BrBbVertex.h" #endif #ifndef BRAT_BrCalibrationManager #include "BrCalibrationManager.h" #endif #ifndef BRAT_BrTofCalibration #include "BrTofCalibration.h" #endif //____________________________________________________________________ ClassImp(BrTMrsFRdoModule); //____________________________________________________________________ BrTMrsFRdoModule::BrTMrsFRdoModule() { // Default constructor. DO NOT USE SetState(kSetup); SetDefaultParameters(); } //____________________________________________________________________ BrTMrsFRdoModule::BrTMrsFRdoModule(const Char_t* name, const Char_t* title) : BrModule(name, title) { // Named Constructor SetState(kSetup); SetDefaultParameters(); } //___________________________________________________________ void BrTMrsFRdoModule::SetDefaultParameters() { // default parameters, called in constructors SetNtuple(); // default is kFALSE SetUseBbVertex(); // default is kTRUE fPositionOffset = 0.0; // cm fTimeOffset = 0.0; // nsec } //____________________________________________________________________ void BrTMrsFRdoModule::DefineHistograms() { // Define histograms. They are: // <fill in here> if (GetState() != kInit) { Stop("DefineHistograms", "Must be called after Init"); return; } TDirectory* saveDir = gDirectory; TDirectory* histDir = gDirectory->mkdir(GetName()); histDir->cd(); // Make histograms here fhTimeSlat = new TH1F("timeSlat","Time vs slats", 1024, 0, 128); fhTimeUpSlat = new TH1F("timeUpSlat","TimeUp vs slats", 1024, 0, 128); fhTimeDownSlat = new TH1F("timeDownSlat","TimeDown vs slats", 1024, 0, 128); fhAdcUpSlat = new TH1F("adcUpSlat","Normalized Adc Up vs Slat",1024, 0, 8196); fhAdcDownSlat = new TH1F("adcDownSlat","Normalized Adc Down vs Slat",512, 0, 4096); fhPosSlat = new TH1F("posSlat","Average Pos vs Slat",100,-25, 25); Info(20,"Init","fUseNtuple = %dn",fUseNtuple); if(fUseNtuple) { fTMrsFTree = new TNtuple(Form("%sRdoTree", GetName()), Form("%s rdo hits tree", GetName()), "tof:tdcup:tdcdown:deup:dedown:y:" "bbt0:bbtl:bbtr:z"); Info(30,"Init","fTMrsFTree created = %d",Int_t(fTMrsFTree)); } gDirectory = saveDir; } //____________________________________________________________________ void BrTMrsFRdoModule::Init() { // Job-level initialisation SetState(kInit); if (DebugLevel() > 1) cout << "Entering Init of BrTMrsFRdoModule" << endl; // calibration BrCalibrationManager *parMan = BrCalibrationManager::Instance(); if (!parMan) { Abort("Init", "Couldn't instantiate calibration manager"); return; } // get the calibration element (created in main) fCalibration = (BrTofCalibration*)parMan->Register("BrTofCalibration", GetName()); if (!fCalibration) { Abort("Init", "Couldn't get calibration parameters for %s", GetName()); return; } // check calibration mode BrCalibration::EAccessMode mode = BrCalibration::kRead; const Int_t nSlats = 1; // pedestals if (!fCalibration->GetAccessMode("topPedestal") || !fCalibration->GetAccessMode("botPedestal")) { fCalibration->Use("topPedestal", mode, nSlats); fCalibration->Use("botPedestal", mode, nSlats); } if (!fCalibration->GetAccessMode("topTdcGain") || !fCalibration->GetAccessMode("botTdcGain")) { fCalibration->Use("topTdcGain", mode, nSlats); fCalibration->Use("botTdcGain", mode, nSlats); } if (!fCalibration->GetAccessMode("topGap") || !fCalibration->GetAccessMode("botGap")) { fCalibration->Use("topGap", mode, nSlats); fCalibration->Use("botGap", mode, nSlats); } if (!fCalibration->GetAccessMode("topGapStart") || !fCalibration->GetAccessMode("botGapStart")) { fCalibration->Use("topGapStart", mode, nSlats); fCalibration->Use("botGapStart", mode, nSlats); } if(!fCalibration->GetAccessMode("topAdcGain")){ fCalibration->Use("topAdcGain", mode, nSlats); fCalibration->Use("botAdcGain", mode, nSlats); } } //____________________________________________________________________ void BrTMrsFRdoModule::Begin() { // Run-level initialisation SetState(kBegin); } //____________________________________________________________________ void BrTMrsFRdoModule::Event(BrEventNode* inNode, BrEventNode* outNode) { // Per event method SetState(kEvent); Info(1,"Event","Starting Event method"); BrDataTable* table = inNode->GetDataTable ("DigTof TMrsF"); if(!table) return; BrBbVertex* bb=0; if(fUseBbVertex){ bb = (BrBbVertex*) inNode->GetObject("BB Vertex"); if (!bb) { if (DebugLevel() > 2) Warning("Event", "No Beam Beam vertex"); } } Float_t tdcUp, tdcDown; Float_t adcUp, adcDown; int nh= table->GetEntries(); Info(5,"Event","Table entries = %d",nh); for(int k=0; k < nh; k++){ BrTofDig* p = (BrTofDig*) table->At(k); if(p) { int slatc = 1; Float_t tPed = fCalibration->GetTopPedestal(slatc); Float_t bPed = fCalibration->GetBotPedestal(slatc); tdcUp= p->GetTdcUp(); tdcDown= p->GetTdcDown(); Float_t topGapStart = fCalibration->GetTopAdcGapStart(slatc); Float_t topGap = fCalibration->GetTopAdcGap(slatc); Float_t botGap = fCalibration->GetBotAdcGap(slatc); Float_t botGapStart = fCalibration->GetBotAdcGapStart(slatc); if(p->GetAdcUp() < topGapStart) adcUp= p->GetAdcUp(); else adcUp= p->GetAdcUp() - topGap; if(p->GetAdcDown() < botGapStart) adcDown= p->GetAdcDown(); else adcDown= p->GetAdcDown() -botGap; adcDown -=bPed; adcUp -=tPed; Info(20,"Event","adcDown = %f, adcUp = %f",adcDown,adcUp); } } BrTMrsFRdo* TMrsFRdo = new BrTMrsFRdo("TMrsFRdo","rdo"); outNode->AddObject(TMrsFRdo); for(int k=0;k<1;k++){ if(tdcUp<4000 && tdcDown<4000){ int slat = 1; Float_t tTGain = fCalibration->GetTopTdcGain(slat); Float_t bTGain = fCalibration->GetBotTdcGain(slat); Float_t tAdcGain = fCalibration->GetTopAdcGain(slat); Float_t bAdcGain = fCalibration->GetBotAdcGain(slat); Float_t tdcup =tdcUp*tTGain; Float_t tdcdown= tdcDown*bTGain; float deup = adcUp/tAdcGain*1000.; float dedown = adcDown/bAdcGain*1000.; if(deup<0.3 || dedown < 0.3) break; // // Poor mans slewing - but looks good when adjusting to Ntuple's (deup:tdcup-bbr) // float dtup = -1.2/TMath::Sqrt(deup/1000)+1.2; float dtdown = -1.2/TMath::Sqrt(dedown/1000.)+1.2; dtup=0; dtdown=0; Float_t t0 = (tdcup+dtup+tdcdown+dtdown)/2 - fTimeOffset; // nsec Float_t x = (tdcup+dtup-tdcdown-dtdown)*15. * 0.467 - fPositionOffset; Float_t y = 0.0; Float_t de = (deup+dedown)/2.0; TMrsFRdo->AddHit(t0, de, y, x); Info(20,"Event","t0 = %f, de = %f, y = %f, x = %f",t0,de,y,x); if(HistOn()){ fhTimeUpSlat->Fill(tdcup); fhTimeDownSlat->Fill(tdcdown); fhAdcUpSlat->Fill( deup); fhAdcDownSlat->Fill(dedown); fhTimeSlat->Fill(t0); fhPosSlat->Fill(x); if (fUseNtuple) { Float_t bbt0 = 999; Float_t bbz = 999; Float_t bbtl = 999; Float_t bbtr = -999; if(bb){ bbt0 = bb->GetTime0(); bbz = bb->GetZ0(); bbtl = bb->GetLeftTime(); bbtr = bb->GetRightTime(); } const Float_t ntVar[] = { t0, tdcup, tdcdown, deup, dedown, x, bbt0, bbtl, bbtr, bbz }; Info(30,"Event","Incrementing fTMrsFTree = %d",Int_t(fTMrsFTree)); fTMrsFTree->Fill(ntVar); } } } } // // if(Verbose()>10){ TMrsFRdo->Print("ad"); } } //____________________________________________________________________ void BrTMrsFRdoModule::End() { // Run-level finalisation SetState(kEnd); } //____________________________________________________________________ void BrTMrsFRdoModule::Finish() { // Job-level finalisation SetState(kFinish); } //____________________________________________________________________ void BrTMrsFRdoModule::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: Flemming Videbaek" << endl << " Last Modifications: " << endl << " $Author: hagel $" << endl << " $Date: 2002/08/06 19:10:59 $" << endl << " $Revision: 1.2 $ " << endl << endl << "-------------------------------------------------" << endl; } //____________________________________________________________________ // // $Log: BrTMrsFRdoModule.cxx,v $ // Revision 1.2 2002/08/06 19:10:59 hagel // Add verbosity for debugging; Change when ntuple is incremented so that it doesn't crash when HistOn() is false and fUseNtuple is true // // Revision 1.1 2002/06/07 16:30:11 videbaek // Add BrTMrsFRdoModule.h to the rdo directory. This hits are used in the tof // determination for Mrs data in the pp run. // |
||||||
This page automatically generated by script docBrat by Christian Holm |
Copyright ; 2002 BRAHMS Collaboration
<brahmlib@rcf.rhic.bnl.gov>
|