//-------------------------------------------------------------------------- // File and Version Information: // $Id: FsmDrcDisc.cc,v 1.7 2007/05/24 08:07:40 klausg Exp $ // // Description: // Class FsmDrcBarrel // // Implementation of the Disc DIRC for the FastSim // // This software was developed for the PANDA collaboration. If you // use all or part of it, please give an appropriate acknowledgement. // // Author List: // Klaus Goetzen Original Author // // Copyright Information: // Copyright (C) 2006 GSI // //------------------------------------------------------------------------ //----------------------- // This Class's Header -- //----------------------- #include "PndFsmDrcDisc.h" //------------- // C Headers -- //------------- //--------------- // C++ Headers -- //--------------- #include #include using std::cout; using std::endl; using std::ostream; using std::string; //------------------------------- // Collaborating Class Headers -- //------------------------------- #include "ArgList.h" #include "PndFsmResponse.h" #include "PndFsmTrack.h" #include "TH2F.h" #include "TFile.h" //----------------------------------------------------------------------- // Local Macros, Typedefs, Structures, Unions and Forward Declarations -- //----------------------------------------------------------------------- //---------------- // Constructors -- //---------------- PndFsmDrcDisc::PndFsmDrcDisc() { initParameters(); _thtMin=_thtMin*M_PI/180.0; _thtMax=_thtMax*M_PI/180.0; readParameters(); print(std::cout); } PndFsmDrcDisc::PndFsmDrcDisc(ArgList &par) { initParameters(); //set default parameter values and parses a parameter list //i.e. std::list of the form //"a=1" "b=2" "c=3" parseParameterList(par); _thtMin=_thtMin*M_PI/180.0; _thtMax=_thtMax*M_PI/180.0; readParameters(); print(std::cout); } //-------------- // Destructor -- //-------------- PndFsmDrcDisc::~PndFsmDrcDisc() { } //-------------- // Operations -- //-------------- PndFsmResponse* PndFsmDrcDisc::respond(PndFsmTrack *t) { PndFsmResponse *result=new PndFsmResponse(); result->setDetector(this); bool wasDetected=detected(t); result->setDetected(wasDetected); if (wasDetected && fabs(t->charge())>1e-8) { TParticlePDG* part = _fdbPDG->GetParticle(t->pdt()); double mass = (part) ? part->Mass() : t->p4().M(); double theta = t->p4().Theta(); double p=t->p4().Vect().Mag(); double ctht=cos(theta); //double p_t=p*sin(theta); int pid=abs(t->pdt()); int npid=-1; if (pid==11) npid=0; else if (pid==13) npid=1; else if (pid==211) npid=2; else if (pid==321) npid=3; else if (pid==2212) npid=4; double thtdeg=theta/M_PI*180.; if ( p==0 || ctht==0 ) {result->setDetected(false);return result;} // floating point check ******************** double lambda1 = 280e-9; //range of wavelength, which is seen by the PMT/PD double lambda2 = 330e-9; double alpha=7.2974e-3; //finestructure constant // curvature of track due to magnet field //double r = 3.3356 * p_t / _Bfield; // dip angle in phi direction (due to curvature of track in magnet field) //double psi = acos(_rBarrel/(2*r)); // path length in radiator double l = _dDisc/cos(theta);// _dSlab*sqrt( 1/(sin(theta)*sin(theta)) + 1/(tan(psi)*tan(psi)) ); // deteremine trapping fraction double trapped = _trap; if (npid>=0 && trapfrac[npid]) trapped = npid<0 ? 0.0 : trapfrac[npid]->GetBinContent(trapfrac[npid]->FindBin(p<6.0?p:6.0,thtdeg)); // estimate the number of initially produced cherenkov photons double nPhot = 2*M_PI*alpha*l*(1./lambda1 - 1./lambda2)*(1 - (mass*mass+p*p)/(p*p*_nRefrac*_nRefrac)); // dice a poisson value nPhot = _rand->Poisson(nPhot); // determine the number of photons hitting the sensor nPhot *= trapped*_effNPhotons; // ************** reset detected and quit due to low numbers of photons if (nPhot<=_nPhotMin) { result->setDetected(false); return result; } if (nPhot>100) nPhot=100; // overall resolution for the tht_c measurement double sig = _dthtc/sqrt(nPhot); double thtC = compThetaC(p,mass); double m_e = _fdbPDG->GetParticle(11)->Mass(); double m_mu = _fdbPDG->GetParticle(13)->Mass(); double m_pi = _fdbPDG->GetParticle(211)->Mass(); double m_K = _fdbPDG->GetParticle(321)->Mass(); double m_p = _fdbPDG->GetParticle(2212)->Mass(); // compute the expected cherenkov angles for all particle types // we need these to determine the pdf's (gaussian around nominal tht_c with res sig) // this Likelihood function has to be evaluated for the _measured_ momentum, which // is smeared with dp! double measp=_rand->Gaus(p,_dp*p); if (measp<0) measp=0; double thtc_e = compThetaC(measp,m_e); double thtc_mu = compThetaC(measp,m_mu); double thtc_pi = compThetaC(measp,m_pi); double thtc_K = compThetaC(measp,m_K); double thtc_p = compThetaC(measp,m_p); double measThetaC = _rand->Gaus(thtC,sig); if (measThetaC<0) measThetaC=0; result->setDrcDiscThtc(measThetaC,sig); if (thtc_e) result->setLHElectron( gauss(measThetaC,thtc_e,sig) ); if (thtc_mu) result->setLHMuon( gauss(measThetaC,thtc_mu,sig) ); if (thtc_pi) result->setLHPion( gauss(measThetaC,thtc_pi,sig) ); if (thtc_K) result->setLHKaon( gauss(measThetaC,thtc_K,sig) ); if (thtc_p) result->setLHProton(gauss(measThetaC,thtc_p,sig) ); } return result; } double PndFsmDrcDisc::gauss(double x, double x0, double s) { return (1.0/(sqrt(2.0*M_PI)*s))* exp(-(x-x0)*(x-x0)/(2.0*s*s)); } double PndFsmDrcDisc::compThetaC(double p, double m) { double val=0.0; if (p==0) return 0.0; if ( (val = (p*p+m*m)/(p*p*_nRefrac*_nRefrac)) <= 1.0 ) return acos(sqrt(val)); else return 0.0; } bool PndFsmDrcDisc::detected(PndFsmTrack *t) const { if (t->hitMapValid()) { return t->hitMapResponse(FsmDetEnum::Drc); } else { int lundId=abs(t->pdt()); TParticlePDG* part = _fdbPDG->GetParticle(lundId); double mass = (part) ? part->Mass() : t->p4().M(); double p_cerenkov_min=mass/sqrt(_nRefrac*_nRefrac - 1.0); double theta = t->p4().Theta(); double p=t->p4().Vect().Mag(); bool correctPidType=(lundId==11 || lundId==13 || lundId==211 || lundId==321 || lundId==2212); return ( p>p_cerenkov_min && theta>=_thtMin && theta<=_thtMax && correctPidType && _rand->Uniform()<=_efficiency); } } void PndFsmDrcDisc::print(ostream &o) { o <<"Parameters for detector <"<"<IsZombie()) { cout <<" -W- (PndFsmDrcDisc::readParameters) - file "<<_parFileName.c_str() <<" doesn't exist. Using constant trapping fraction _trap="<<_trap<Get("hacc0"); trapfrac[1]=(TH2F*)f->Get("hacc1"); trapfrac[2]=(TH2F*)f->Get("hacc2"); trapfrac[3]=(TH2F*)f->Get("hacc3"); trapfrac[4]=(TH2F*)f->Get("hacc4"); for (int i=0;i<5;i++) trapfrac[i]->SetDirectory(0); f->Close(); } delete f; return true; }