|
//____________________________________________________________________ // // C1 pid module - makes particle identification from C1 data. // // Finds Ffs track table and C1Rdo object in inNode. Matches tracks // and hits in the sensitive plane of C1 and returns (called C1Pid) // data table with BrC1Pid objects to the outNode. // // The track is pointed to the backplane of C1. The blob (of cherenkov // light) has an overlap with a number of photo tubes (0<n<10). The // calibrated adc signal from these tubes are added - this is called // the blob energy. // // If two blobs (corresponding to two tracks) have overlapping tubes // they should be separated - this method is not implemented yet. // //____________________________________________________________________ // // $Id: BrC1PidModule.cxx,v 1.13 2002/06/10 13:29:50 ekman Exp $ // $Author: ekman $ // $Date: 2002/06/10 13:29:50 $ // $Copyright: (C) 2001 BRAHMS Collaboration <brahmlib@rhic.bnl.gov> // #ifndef BRAT_BrC1PidModule #include "BrC1PidModule.h" #endif #ifndef BRAT_BrGeometryDbManager #include "BrGeometryDbManager.h" #endif #ifndef BRAT_BrParameterDbManager #include "BrParameterDbManager.h" #endif #ifndef BRAT_BrChkvRdo #include "BrChkvRdo.h" #endif #ifndef BRAT_BrSpectrometerTracks #include "BrSpectrometerTracks.h" #endif #ifndef BRAT_BrTofPid #include "BrTofPid.h" #endif #ifndef BRAT_BrC1Pid #include "BrC1Pid.h" #endif #ifndef ROOT_TMath #include "TMath.h" #endif #ifndef ROOT_TDirectory #include "TDirectory.h" #endif #ifndef ROOT_TRandom #include "TRandom.h" #endif #ifndef WIN32 #include <iostream> #else #include <iostream.h> #endif //____________________________________________________________________ ClassImp(BrC1PidModule); //____________________________________________________________________ BrC1PidModule::BrC1PidModule() { // Default constructor. DO NOT USE SetState(kSetup); } //____________________________________________________________________ BrC1PidModule::BrC1PidModule(const Char_t* name, const Char_t* title) : BrModule(name, title) { // Named Constructor SetState(kSetup); SetBackPlaneOffset(); SetBlobRadius(); SetCheckGeometry(); // default is kFALSE SetEstimateBackground(); // default is kFALSE SetBetaCut(); SetRefractiveIndex(); SetUseH1Beta(); } //____________________________________________________________________ void BrC1PidModule::DefineHistograms() { // Define histograms. They are: if (GetState() != kInit) { Stop("DefineHistograms", "Must be called after Init"); return; } TDirectory* saveDir = gDirectory; fHistDir = gDirectory->mkdir(GetName()); fHistDir->cd(); // Make histograms here Char_t* formula_c1pion = Form("%f*75*sin(acos((1/%f)*sqrt(1+(0.139**2)/(x**2))))",3.6,fRefractiveIndex); Char_t* formula_c1kaon = Form("%f*75*sin(acos((1/%f)*sqrt(1+(0.494**2)/(x**2))))",3.6,fRefractiveIndex); fPionEVsP = new TF1("photonsVsP_pions", formula_c1pion, -20, 20); fKaonEVsP = new TF1("photonsVsP_kaons", formula_c1kaon, -20, 20); fHistDir->Add(fPionEVsP); fHistDir->Add(fKaonEVsP); fBlobEnergy = new TH1F("blobEnergy", "Summed adc for tubes hit by blob",255,-1,50); fBlobEnergyVsP = new TH2F("blobEnergyVsP", "Blob energy versus track momentum", 200,-15,15,255,-1,50); fNTubesInBlob = new TH1F("nTubesHitByBlob", "Number of tubes hit by each blob", 10,-0.5,9.5); fBlobEnergyVsNTubes = new TH2F("blobEnergyVsNTubes", "Summed adc vs number of tubes", 10,-0.5,9.5,100,0,200); fBackground = new TH1F("background", "Background estimate",255,-1,50); fBackgroundToSignal = new TH1F("backgroundToSignal", "background/signal",255,-1,50); fOverlaps = new TH1F("numberOfOverlaps", "Number of overlapping tubes in events with nTrack>1", 10,-0.5, 9.5); fOverlapEnergy = new TH1F("overlapEnergy", "Calibrated adc signal in overlapping tubes", 100,0,10); fOverlapEnergyR = new TH1F("overlapEnergyRelativ", "Overlap energy divided by total blob energy", 500,-0.1,1.1); fConflicts = new TH1F("numberOfConflicts", "Number of not-solved conflicts in events with nTrack>1", 10,-0.5, 9.5); fConfLevel = new TH1F("confidenceLevel", "Confidence Level of C1 pids", 110,-1.1, 1.1); fStatus = new TH1F("status","GhostTrack:-2, Conflict:-1, Ok:1",6,-3.5,2.5); fHistDir->mkdir("Geometry"); fHistDir->cd("Geometry"); fTubeToHit = new TH2F("dXdY_TubeToHit","dXdY",80,-40,40,50,-25,25); fTubeToHitX = new TH1F("dX_TubeToHit","dX",80,-40,40); fTubeToHitY = new TH1F("dY_TubeToHit","dY",50,-25,25); fIntWithBackPlane = new TH2F("intersectionWithBackPlane", "Intersection With Back Plane", 100,-50,50,100,-20,20); fUsedBlobRadius = new TH1F("blobRadius","Radius of light blob", 100, -0.5, 4); fSnapshots = new TList(); fHistDir->cd(); fHistDir->mkdir("Snapshots"); fHistDir->cd("Snapshots"); for (Int_t i = 0; i <= fMaxNSnapshots; i++) fSnapshots->Add(new TH2F(Form("Snapshot%03d", i), Form("Snapshot of C1 plane, no.%d", i), 10,-31.75,31.75,6,-19.05,19.05)); fBlobEnergy ->SetXTitle("Number of photons detected"); fBlobEnergyVsP ->SetXTitle("#frac{1}{q} reconstructed track momentum [GeV/c]"); fBlobEnergyVsP ->SetYTitle("Number of photons detected"); fNTubesInBlob ->SetXTitle("Number of tubes in blob"); fTubeToHit ->SetXTitle("hitX - hitTubeX [cm]"); fTubeToHit ->SetYTitle("hitY - hitTubeY [cm]"); fIntWithBackPlane->SetXTitle("track hit x [cm]"); fIntWithBackPlane->SetYTitle("track hit y [cm]"); fBackground ->SetXTitle("Number of photons detected"); fBackgroundToSignal->SetXTitle("Number of photons detected"); gDirectory = saveDir; } //____________________________________________________________________ void BrC1PidModule::Init() { // Job-level initialisation SetState(kInit); // get geometry and parameters BrGeometryDbManager *geoDb = BrGeometryDbManager::Instance(); BrParameterDbManager *parDb = BrParameterDbManager::Instance(); if (!geoDb) { Failure("Init ", "Couldn't instantiate geometry manager"); return; } if (!parDb) { Failure("Init ", "Couldn't instantiate detector parameter manager"); return; } fC1Params = (BrChkvParameters*)parDb->GetDetectorParameters("BrChkvParameters", "C1"); fC1Volume = (BrDetectorVolume*)geoDb->GetDetectorVolume("BrDetectorVolume","C1"); if(!fC1Volume) Abort("Init","No C1 Detector Volume : %s", GetName()); fT2Volume = (BrDetectorVolume*)geoDb->GetDetectorVolume("BrDetectorVolume","T2"); fC1BackPlane = fC1Volume->GetBackPlane(); fPionThreshold = 0.13957 / TMath::Sqrt(fRefractiveIndex*fRefractiveIndex-1); } //____________________________________________________________________ void BrC1PidModule::Begin() { // Run-level initialisation SetState(kBegin); } //____________________________________________________________________ void BrC1PidModule::Event(BrEventNode* inNode, BrEventNode* outNode) { // Requires: ChkvRdo C1 and FfsTracks objects in inNode or outNode // Generates: BrDataTable with BrC1Pid objects // Per event method SetState(kEvent); // get ffs track table BrDataTable* ffsTrackTable = (BrDataTable*) inNode->GetDataTable("FfsTracks"); if (!ffsTrackTable) { ffsTrackTable = (BrDataTable*) outNode->GetDataTable("FfsTracks"); if (!ffsTrackTable) { if (Verbose()>10) Warning("Event","Could not find Ffs track table."); return; } } // return if there is no tracks Int_t nTracks = ffsTrackTable->GetEntries(); if (nTracks==0) return; // get c1 rdo BrChkvRdo* c1Rdo = (BrChkvRdo*) inNode->GetObject("ChkvRdo C1"); if (!c1Rdo) { c1Rdo = (BrChkvRdo*) outNode->GetObject("ChkvRdo C1"); if (!c1Rdo) { if (Verbose()>10) Warning("Event","Could not find ChkvRdo C1."); return; } } // get tof pid (for beta) BrDataTable* tofPidTable = (BrDataTable*) inNode->GetDataTable("TofPid TOF1"); if (!tofPidTable) tofPidTable = (BrDataTable*) outNode->GetDataTable("TofPid TOF1"); // create data table with BrC1Pid objects BrDataTable* c1PidTable = new BrDataTable("C1Pid","C1 pid object table"); outNode->AddDataTable(c1PidTable); // if there is only one global track if (nTracks==1) { BrFfsTrack* track = (BrFfsTrack*)ffsTrackTable->At(0); Float_t p = track->GetMomentum(); Float_t blobR = fBlobRadius; Float_t beta = 99; if (fUseH1Beta) { //get matching h1 pid BrTofPid* h1Pid; if (tofPidTable) { h1Pid = (BrTofPid*) tofPidTable->At(0); if (h1Pid) if (h1Pid->GetTrackId()==track->GetId()) beta = h1Pid->GetBeta(); } } blobR = EstimateBlobRadius(p,beta); BrVector3D c1BPHit = PointToC1BackPlane(track->GetBackTrack()); BLOB blob; FindTubesHitByBlob(c1Rdo,c1BPHit, blobR, blob); Float_t energy = blob.eSum; if (track->GetStatus()==1) c1PidTable->Add(CreatePid(energy,1,track->GetId())); else { if (HistOn()) fStatus->Fill(-2); return; } if (HistOn()) { fBlobEnergy ->Fill(energy); fBlobEnergyVsP ->Fill(p,energy); fBlobEnergyVsNTubes->Fill(blob.n,energy); fUsedBlobRadius->Fill(blobR); fStatus->Fill(1); Snapshot(c1Rdo, c1BPHit, p); //check geometry if p>3 if (fCheckGeometry && TMath::Abs(p)>3) CheckGeometry(c1Rdo, c1BPHit); if (fEstimateBackground) EstimateBackground(c1Rdo, blob); } } // if there is more than one global track if (nTracks>1) { //make usefull arrays: Float_t p[nTracks]; Int_t trackId[nTracks]; Int_t trackStatus[nTracks]; BrVector3D c1BPHit[nTracks]; BLOB blob[nTracks]; Bool_t conflict[nTracks]; Float_t blobR[nTracks]; for (Int_t i=0; i<nTracks; i++){ BrFfsTrack* track = (BrFfsTrack*)ffsTrackTable->At(i); p[i] = track->GetMomentum(); trackStatus[i] = track->GetStatus(); blobR[i] = fBlobRadius; Float_t beta = 99; if (fUseH1Beta) { //get matching h1 pid BrTofPid* h1Pid; if (tofPidTable) for (Int_t h=0; h<tofPidTable->GetEntries(); h++) { h1Pid = (BrTofPid*) tofPidTable->At(h); if (h1Pid) if (h1Pid->GetTrackId()==track->GetId()) beta = h1Pid->GetBeta(); } } blobR[i] = EstimateBlobRadius(p[i],beta); trackId[i] = track->GetId(); c1BPHit[i] = PointToC1BackPlane(track->GetBackTrack()); conflict[i] = kFALSE; FindTubesHitByBlob(c1Rdo,c1BPHit[i], blobR[i], blob[i]); } Int_t nConflicts =0; //loop over ffs tracks for (Int_t i=0; i<nTracks-1; i++){ BrFfsTrack* track1 = (BrFfsTrack*)ffsTrackTable->At(i); if ((TMath::Abs(p[i])<fPionThreshold)|| (trackStatus[i]!=1)) continue; for(Int_t j=i+1; j<nTracks; j++){ BrFfsTrack* track2 = (BrFfsTrack*)ffsTrackTable->At(i); if ((TMath::Abs(p[j])<fPionThreshold)|| (trackStatus[i]!=1)) continue; //check if there are overlapping tubes if (Overlap(blob[i],blob[j])) { //try to seperate blobs if (SeparateBlobs(blob[i],blob[j])<0) { conflict[i] = kTRUE; conflict[j] = kTRUE; blob[i].eSum = -blob[i].eSum;; nConflicts++; } } } } if (HistOn()) fConflicts->Fill(nConflicts); for (Int_t i=0; i<nTracks; i++) { // if no conflicts and good track if (!conflict[i] && trackStatus[i]==1) c1PidTable->Add(CreatePid(blob[i].eSum,1,trackId[i])); if (conflict[i] && trackStatus[i]==1) c1PidTable->Add(CreatePid(blob[i].eSum,1,trackId[i])); if (conflict[i] && trackStatus[i]!=1) c1PidTable->Add(CreatePid(-1,-1,trackId[i])); if (HistOn()) { fBlobEnergy ->Fill(blob[i].eSum); fBlobEnergyVsP ->Fill(p[i],blob[i].eSum); fUsedBlobRadius->Fill(blobR[i]); Snapshot(c1Rdo, c1BPHit[i], p[i]); if (conflict[i]) fStatus->Fill(-1); else fStatus->Fill(1); } if (DebugLevel()>10) { cout << " BrC1PidModule: " << nTracks << " tracks: "; cout << " Track " << i << " in conflict with " << nConflicts << " other tracks." << endl; } } } // endif more than one track } //____________________________________________________________________ void BrC1PidModule::FindTubesHitByBlob(BrChkvRdo* allHits, const BrVector3D& xy, Float_t r, BLOB& blob) { blob.n = 0; blob.eSum = 0; blob.xy = xy; blob.r = r; if (blob.r <= 0) return; for (Int_t h = 0; h<allHits->GetEntries(); h++){ BrChkvRdo::BrChkvHit* c1Hit = (BrChkvRdo::BrChkvHit*) allHits->GetHit(h); if (IsTubeHitByBlob(c1Hit->GetTubeNum(), xy, r)) { blob.t[blob.n] = c1Hit->GetTubeNum(); blob.e[blob.n] = c1Hit->GetEnergy(); blob.eSum = blob.eSum + blob.e[blob.n]; blob.n++; } } if (HistOn()) fNTubesInBlob->Fill(blob.n); return; } //____________________________________________________________________ Bool_t BrC1PidModule::IsTubeHitByBlob(Int_t tubeNo, const BrVector3D& xy, Float_t r) { // finds the blob circle (center = (cX,cY), radius = fBlobRadius) // and the tube square of tube with tubeNo. Is there and overlap? // returns kTRUE is yes, kFALSE if no. if (tubeNo<0) return kFALSE; Float_t cX = xy.X(); Float_t cY = xy.Y(); Float_t tubeX = fC1Params->GetTubePosition(tubeNo).GetX(); Float_t tubeY = fC1Params->GetTubePosition(tubeNo).GetY(); Float_t tSize = fC1Params->GetTubeSize(); Float_t tubeX1 = tubeX - tSize/2; Float_t tubeX2 = tubeX + tSize/2; Float_t tubeY1 = tubeY - tSize/2; Float_t tubeY2 = tubeY + tSize/2; /* // Tracks pointing inside the outer area can create photons that hit // the sensitive tube area: // // | | // | s | // c |_____| c // / // ----| +-----+ |---- // | |+++++| | // s | |++T++| | s // | |+++++| | // ----| +-----+ |---- // _______/ // c | | c // | s | // | | // // - first check if one the corners is inside the blob. // - then check if the blob has an overlap at the sides. */ // check distance to tube corners if ((TMath::Sqrt((cX-tubeX1)*(cX-tubeX1)+(cY-tubeY1)*(cY-tubeY1)) < r)|| (TMath::Sqrt((cX-tubeX1)*(cX-tubeX1)+(cY-tubeY2)*(cY-tubeY2)) < r)|| (TMath::Sqrt((cX-tubeX2)*(cX-tubeX2)+(cY-tubeY1)*(cY-tubeY1)) < r)|| (TMath::Sqrt((cX-tubeX2)*(cX-tubeX2)+(cY-tubeY2)*(cY-tubeY2)) < r)) return kTRUE; // check sides if ((cX>tubeX1) && (cX<tubeX2)) if (((cY-r)<tubeY2) && ((cY+r)>tubeY1)) return kTRUE; if ((cY>tubeY1) && (cY<tubeY2)) if (((cX-r)<tubeX2) && ((cX+r)>tubeX1)) return kTRUE; return kFALSE; } //____________________________________________________________________ Bool_t BrC1PidModule::Overlap(BLOB& blob1, BLOB& blob2) { // check if there are overlapping tubes in the two arrays. Int_t nOverlaps=0; for(Int_t i=0; i<blob1.n; i++) { for(Int_t j=0; j<blob2.n; j++) { if (blob1.t[i] == blob2.t[j]) { nOverlaps++; } } } if (HistOn()) fOverlaps->Fill(nOverlaps); if (nOverlaps>0) return kTRUE; return kFALSE; } //____________________________________________________________________ Int_t BrC1PidModule::SeparateBlobs(BLOB& blob1, BLOB& blob2) { // For now this method always returns -1!!!!! // try to separate blobs (arrays of hits) // // 1) Check if the blobs have energy above threshold // if one/both of them does/do not, there's no idea in trying // to separate them. // 2) Check if they are to close to separate // 3) ...??? // // return 1 if separation succeeded, -1 if not. /* if (HistOn()) { fOverlapEnergy->Fill(energyInOverlap); fOverlapEnergyR->Fill(energyInOverlap/energyOfBlob1); fOverlapEnergyR->Fill(energyInOverlap/energyOfBlob2); } */ return -1; } //____________________________________________________________________ Float_t BrC1PidModule::EstimateBlobRadius(Float_t p, Float_t beta) { // Estimate the blob radius from the momentum (asume pion) // r = L * tan(theta) Float_t blobR = 0; if ((beta > fBetaCut)||(!fUseH1Beta)) blobR = 75*tan(TMath::ACos((1/fRefractiveIndex)*TMath::Sqrt(1+(0.139*0.139)/(p*p)))); if (blobR > fBlobRadius) blobR = fBlobRadius; return blobR; } //____________________________________________________________________ void BrC1PidModule::CheckGeometry(BrChkvRdo* allHits, const BrVector3D& xy) { //fills "dist to hit" histogram - used for geometry alignment //loop over C1 hits for (Int_t i = 0; i<allHits->GetEntries(); i++){ BrChkvRdo::BrChkvHit* c1Hit = (BrChkvRdo::BrChkvHit*) allHits->GetHit(i); // only look at "good" hits... Float_t energy = c1Hit->GetEnergy(); if (energy>3 && energy<25) { Int_t tubeNo = c1Hit->GetTubeNum(); Float_t xDiff = xy.X() - fC1Params->GetTubePosition(tubeNo).GetX(); Float_t yDiff = xy.Y() - fC1Params->GetTubePosition(tubeNo).GetY(); fTubeToHit ->Fill(xDiff,yDiff,energy); fTubeToHitX->Fill(xDiff, energy); fTubeToHitY->Fill(yDiff, energy); } } } //____________________________________________________________________ void BrC1PidModule::EstimateBackground(BrChkvRdo* allHits, BLOB& blob) { // find n tubes that are not in blobHits // n = number of hits in blobHits // sum the adc's TRandom ran = TRandom(846528); Bool_t findNew; Float_t sumE =0; // find blobHits->GetEntries() that are not in blobHits for (Int_t n=0; n<blob.n; n++) { Int_t t; findNew = kTRUE; while(findNew) { t = ran.Integer(32)+1; findNew = kFALSE; // check if tube is in hitsArray for(Int_t i=0; i<blob.n; i++) if (t==blob.t[i]) findNew = kTRUE; } // check if it's in all Hits for(Int_t i=0; i<allHits->GetEntries(); i++) if (t==((BrChkvRdo::BrChkvHit*)allHits->GetHit(i))->GetTubeNum()) sumE = sumE + ((BrChkvRdo::BrChkvHit*)allHits->GetHit(i))->GetEnergy(); } fBackground->Fill(sumE); } //____________________________________________________________________ void BrC1PidModule::Snapshot(BrChkvRdo* allHits, const BrVector3D& xy, Float_t p=99) { if (fNSnapshots<fMaxNSnapshots) { //make the snapshot TH2F* hist = (TH2F*)fSnapshots->FindObject(Form("Snapshot%03d",fNSnapshots)); fNSnapshots++; for(Int_t i=0; i<allHits->GetEntries(); i++) { BrChkvRdo::BrChkvHit* hit = (BrChkvRdo::BrChkvHit*)allHits->GetHit(i); BrVector3D tubePos = fC1Params->GetTubePosition(hit->GetTubeNum()); Float_t energy = hit->GetEnergy(); hist->Fill(tubePos.GetX(),tubePos.GetY(),energy); const Char_t* title = Form("Snapshot%d (track x=%03f,y%06f,p=%06f)",i,xy.X(),xy.Y(),p); hist->SetTitle(title); } } else { if (Verbose()>10) Warning("Snapshot","NSnapshots > MaxNSnapshot. Image not created!"); } } //____________________________________________________________________ BrVector3D BrC1PidModule::PointToC1BackPlane(BrDetectorTrack* t2track) { // Returns the position (in local C1 coordinates) of the // t2 track intersection with C1 backplane. // If outside backplane returns -999,-999,-999 BrLine3D trackLineG = fT2Volume->LocalToGlobal(t2track->GetTrackLine()); BrVector3D c1HitG = fC1BackPlane.GetIntersectionWithLine(trackLineG); BrVector3D c1HitL; fC1Volume->GlobalToLocal(c1HitG,c1HitL,0); // add the offsets c1HitL.SetX(c1HitL.X()-fBackPlaneOffset.X()); c1HitL.SetY(c1HitL.Y()-fBackPlaneOffset.Y()); if (HistOn()) fIntWithBackPlane->Fill(c1HitL.GetX(),c1HitL.GetY()); //check if track hits sensitive plane Float_t dimX = fC1Params->GetTubeSize() * fC1Params->GetNoTubesX(); Float_t dimY = fC1Params->GetTubeSize() * fC1Params->GetNoTubesY(); if((TMath::Abs(c1HitL.GetX())>(0.5*dimX))|| (TMath::Abs(c1HitL.GetY())>(0.5*dimY))) { if(Verbose()>20) Warning("Event","Track outside C1 sensitive plane"); return BrVector3D::BrVector3D(-999, -999, -999); } return c1HitL; } //____________________________________________________________________ BrC1Pid* BrC1PidModule::CreatePid(Float_t e, Float_t cl, Int_t trackId) { BrC1Pid* c1Pid = new BrC1Pid(); c1Pid->SetBlobEnergy(e); c1Pid->SetConfidenceLevel(cl); c1Pid->SetTrackId(trackId); if (HistOn()) fConfLevel->Fill(cl); return c1Pid; } //____________________________________________________________________ void BrC1PidModule::End() { // Run-level finalisation SetState(kEnd); } //____________________________________________________________________ void BrC1PidModule::Finish() { // Job-level finalisation SetState(kFinish); if (!HistOn()) return; TDirectory* saveDir = gDirectory; fHistDir->cd(); const Char_t* name = "blobEnergyVsPProfileX"; Int_t firstYBin = fBlobEnergyVsP->GetYaxis()->FindBin(1); fBlobEnergyVsPProfile = fBlobEnergyVsP->ProfileX(name,firstYBin); TH1F* bg = (TH1F*)fBackground->Clone(); TH1F* s = (TH1F*)fBlobEnergy->Clone(); bg->SetNormFactor(bg->GetEntries()); s->SetNormFactor(s->GetEntries()); // s->Add(bg); fBackgroundToSignal->Divide(bg,s); fBackgroundToSignal->Rebin(2); delete bg; delete s; if (fCheckGeometry) { if (fTubeToHitX->GetEntries()>5000) { Float_t meanX = fTubeToHitX->GetBinCenter(fTubeToHitX->GetMaximumBin()); Float_t meanY = fTubeToHitY->GetBinCenter(fTubeToHitY->GetMaximumBin()); TF1* fitX = new TF1("fitX","gaus"); TF1* fitY = new TF1("fitY","gaus"); fTubeToHitX->Fit(fitX,"Q0","Q0",meanX-5, meanX+5); fTubeToHitY->Fit(fitY,"Q0","Q0",meanY-5, meanY+5); cout << "BrC1PidModule summary:" << endl; if ((TMath::Abs(fitX->GetParameter(1))>0.005) || (TMath::Abs(fitY->GetParameter(1))>0.005)) { cout << "Try SetBackPlaneOffset(" << setprecision(3) << fBackPlaneOffset.X()+fitX->GetParameter(1) << "," << setprecision(3) << fBackPlaneOffset.Y()+fitY->GetParameter(1) << ");" << endl; } } } gDirectory = saveDir; } //____________________________________________________________________ void BrC1PidModule::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: Claus O. E. Jorgensen" << endl << " Last Modifications: " << endl << " $Author: ekman $" << endl << " $Date: 2002/06/10 13:29:50 $" << endl << " $Revision: 1.13 $ " << endl << endl << "-------------------------------------------------" << endl; } //____________________________________________________________________ // // $Log: BrC1PidModule.cxx,v $ // Revision 1.13 2002/06/10 13:29:50 ekman // Now the nphotons is also calculated for tracks that overlaps with other tracks - // the value is set to -nphotons. // // Revision 1.12 2002/04/08 10:44:38 ekman // Added check of ffs track status. If trackStatus!=1 no C1 Pid is made. // // Revision 1.11 2002/04/05 16:01:24 ekman // Introduced fBackPlaneOffset. Fixed a mistake in overlap loop (we // shouldn't check for overlaps if the momentum is below pion threshold) // and a couple of other small things. // // Revision 1.10 2002/02/14 15:13:45 videbaek // Added protection in Init for no volume // // Revision 1.9 2001/11/19 15:03:22 ouerdane // check in Finish if histograms are booked before any actions on histograms, added TDirectory member in the C1 module // // Revision 1.8 2001/11/13 23:05:26 ouerdane // bumped version number to 2.1.32 // // Revision 1.7 2001/11/05 17:36:49 ekman // Major changes - new handling of rdo data, new histograms, new setter functions. // // Revision 1.5 2001/10/25 04:48:43 ekman // Blob radius was set in Init() - moved it to the constructor. // // Revision 1.4 2001/10/19 16:42:44 ekman // Fixed a bug in the check geometry. Added energy vs n hits histogram. // // Revision 1.3 2001/10/17 23:44:56 ekman // Small adjustments in the code. Deleted some histograms, and changed a few // methods. // // Revision 1.2 2001/10/02 01:54:49 ouerdane // removed named contructor calls to Pid data classes // // Revision 1.1 2001/09/03 17:54:13 ekman // Added the BrC1PidModule, BrPid and BrC1Pid classes. // // |
||||||
This page automatically generated by script docBrat by Christian Holm |
Copyright ; 2002 BRAHMS Collaboration
<brahmlib@rcf.rhic.bnl.gov>
|