FairRoot/PandaRoot
PndFsmRich.cxx
Go to the documentation of this file.
1 //--------------------------------------------------------------------------
2 // File and Version Information:
3 // $Id: FsmRich.cc,v 1.1 2007/05/24 08:07:40 klausg Exp $
4 //
5 // Description:
6 // Class FsmRich
7 //
8 // Implementation of the Far Forward Rich for the FastSim
9 //
10 // This software was developed for the PANDA collaboration. If you
11 // use all or part of it, please give an appropriate acknowledgement.
12 //
13 // Author List:
14 // Klaus Goetzen Original Author
15 //
16 // Copyright Information:
17 // Copyright (C) 2007 GSI
18 //
19 //------------------------------------------------------------------------
20 
21 //-----------------------
22 // This Class's Header --
23 //-----------------------
24 #include "PndFsmRich.h"
25 
26 //-------------
27 // C Headers --
28 //-------------
29 
30 //---------------
31 // C++ Headers --
32 //---------------
33 #include <math.h>
34 #include <iostream>
35 
36 using std::cout;
37 using std::endl;
38 using std::ostream;
39 using std::string;
40 
41 //-------------------------------
42 // Collaborating Class Headers --
43 //-------------------------------
44 
45 #include "ArgList.h"
46 #include "PndFsmResponse.h"
47 #include "PndFsmTrack.h"
48 
49 //-----------------------------------------------------------------------
50 // Local Macros, Typedefs, Structures, Unions and Forward Declarations --
51 //-----------------------------------------------------------------------
52 
53 //----------------
54 // Constructors --
55 //----------------
56 
58 {
60 
61  _angleXMax *= M_PI/180.0;
62  _angleYMax *= M_PI/180.0;
63  //print(std::cout);
64 }
65 
67 {
69  //set default parameter values and parses a parameter list
70  //i.e. std::list<std::string> of the form
71  //"a=1" "b=2" "c=3"
72  parseParameterList(par);
73 
74  _angleXMax *= M_PI/180.0;
75  _angleYMax *= M_PI/180.0;
76  //print(std::cout);
77 }
78 
79 //--------------
80 // Destructor --
81 //--------------
82 
84 {
85 }
86 
87 //--------------
88 // Operations --
89 //--------------
90 
93 {
94  PndFsmResponse *result=new PndFsmResponse();
95 
96  result->setDetector(this);
97  bool wasDetected=detected(t);
98  result->setDetected(wasDetected);
99 
100  if (wasDetected && fabs(t->charge())>1e-8)
101  {
102  TParticlePDG* part = _fdbPDG->GetParticle(t->pdt());
103  double mass = (part) ? part->Mass() : t->p4().M();
104  double theta = t->p4().Theta();
105  double p=t->p4().Vect().Mag();
106  double ctht=cos(theta);
107  //double p_t=p*sin(theta);
108 
109  if ( p==0 || ctht==0 ) {result->setDetected(false);return result;} // floating point check ********************
110 
111  double lambda1 = 280e-9; //range of wavelength, which is seen by the PMT/PD
112  double lambda2 = 330e-9;
113 
114  double alpha=7.2974e-3; //finestructure constant
115 
116  // curvature of track due to magnet field
117  //double r = 3.3356 * p_t / _Bfield;
118 
119  // dip angle in phi direction (due to curvature of track in magnet field)
120  //double psi = acos(_rBarrel/(2*r));
121 
122  // path length in radiator
123  double l = _dRich/ctht;// _dSlab*sqrt( 1/(sin(theta)*sin(theta)) + 1/(tan(psi)*tan(psi)) );
124 
125  // estimate the number of cherenkov photons hitting our photosensor
126  double nPhot = _effNPhotons*2*M_PI*alpha*l*(1./lambda1 - 1./lambda2)*(1 - (mass*mass+p*p)/(p*p*_nRefrac*_nRefrac));
127 
128  if (nPhot==0) {result->setDetected(false);return result;} // floating point check ********************
129 
130  nPhot = _rand->Poisson(nPhot);
131 
132  // ************** reset detected and quit due to low numbers of photons
133 
134  if (nPhot<=_nPhotMin) {
135  result->setDetected(false);
136  return result;
137  }
138 
139  if (nPhot>100) nPhot=100;
140 
141  // overall resolution for the tht_c measurement
142  double sig = _dthtc/sqrt(nPhot);
143  double thtC = compThetaC(p,mass);
144 
145  double m_e = _fdbPDG->GetParticle(11)->Mass();
146  double m_mu = _fdbPDG->GetParticle(13)->Mass();
147  double m_pi = _fdbPDG->GetParticle(211)->Mass();
148  double m_K = _fdbPDG->GetParticle(321)->Mass();
149  double m_p = _fdbPDG->GetParticle(2212)->Mass();
150 
151  // compute the expected cherenkov angles for all particle types
152  // we need these to determine the pdf's (gaussian around nominal tht_c with res sig)
153  // this Likelihood function has to be evaluated for the _measured_ momentum, which
154  // is smeared with dp!
155 
156  double measp=_rand->Gaus(p,_dp*p);
157  if (measp<0) measp=0;
158 
159  double thtc_e = compThetaC(measp,m_e);
160  double thtc_mu = compThetaC(measp,m_mu);
161  double thtc_pi = compThetaC(measp,m_pi);
162  double thtc_K = compThetaC(measp,m_K);
163  double thtc_p = compThetaC(measp,m_p);
164 
165  double measThetaC = _rand->Gaus(thtC,sig);
166  if (measThetaC<0) measThetaC=0;
167 
168  result->setRichThtc(measThetaC,sig);
169 
170  if (thtc_e) result->setLHElectron( gauss(measThetaC,thtc_e,sig) );
171  if (thtc_mu) result->setLHMuon( gauss(measThetaC,thtc_mu,sig) );
172  if (thtc_pi) result->setLHPion( gauss(measThetaC,thtc_pi,sig) );
173  if (thtc_K) result->setLHKaon( gauss(measThetaC,thtc_K,sig) );
174  if (thtc_p) result->setLHProton(gauss(measThetaC,thtc_p,sig) );
175 
176  }
177 
178  return result;
179 }
180 
181 double
182 PndFsmRich::gauss(double x, double x0, double s)
183 {
184  return (1.0/(sqrt(2.0*M_PI)*s))*
185  exp(-(x-x0)*(x-x0)/(2.0*s*s));
186 }
187 
188 double
189 PndFsmRich::compThetaC(double p, double m)
190 {
191  double val=0.0;
192  if (p==0) return 0.0;
193  if ( (val = (p*p+m*m)/(p*p*_nRefrac*_nRefrac)) <= 1.0 ) return acos(sqrt(val));
194  else return 0.0;
195 }
196 
197 
198 bool
200 {
201  if (t->hitMapValid()) {
202  return t->hitMapResponse(FsmDetEnum::Drc);
203  } else {
204  int lundId=abs(t->pdt());
205  TParticlePDG* part = _fdbPDG->GetParticle(lundId);
206  double mass = (part) ? part->Mass() : t->p4().M();
207  // double mass = t->pdt()->mass();
208  double p_cerenkov_min=mass/sqrt(_nRefrac*_nRefrac - 1.0);
209 
210  if (t->p4().Z()==0) return false;
211 
212  double angleX = fabs(atan(t->p4().X()/t->p4().Z()));
213  double angleY = fabs(atan(t->p4().Y()/t->p4().Z()));
214  double p=t->p4().Vect().Mag();
215 
216  bool correctPidType=false;
217 
218  if (lundId==11 || lundId==13 || lundId==211 || lundId==321 || lundId==2212) correctPidType=true;
219 
220  return ( p>p_cerenkov_min && angleX<=_angleXMax && angleY<=_angleYMax && correctPidType && _rand->Rndm()<=_efficiency);
221  }
222 }
223 
224 void
226 {
227  o <<"Parameters for detector <"<<detName()<<">"<<endl;
228  o <<" _angleXMax = "<<_angleXMax<<endl;
229  o <<" _angleYMax = "<<_angleYMax<<endl;
230  o <<" _radiationLength = "<<_radiationLength<<endl;
231  o <<" _pmin = "<<_pmin<<endl;
232  o <<" _dthtc = "<<_dthtc<<endl;
233  o <<" _nPhotMin = "<<_nPhotMin<<endl;
234  o <<" _nRefrac = "<<_nRefrac<<endl;
235  o <<" _Bfield = "<<_Bfield<<endl;
236  o <<" _effNPhotons = "<<_effNPhotons<<endl;
237  o <<" _dRich = "<<_dRich<<endl;
238  o <<" _dp = "<<_dp<<endl;
239  o <<" _efficiency = "<<_efficiency<<endl;
240 }
241 
242 void
244 {
245  _detName = "Rich";
246  _angleXMax = 5.0;
247  _angleYMax = 10.0;
248  _radiationLength = 0.0;
249  _pmin = 0.0;
250  _dthtc = 0.01;
251  _nPhotMin = 5;
252  _nRefrac = 1.05;
253  _Bfield = 2.;
254  _effNPhotons = 0.075;
255  _dRich = 1.0;
256  _dp = 0.01;
257  _efficiency =1.0;
258 }
259 
260 bool
261 PndFsmRich::setParameter(std::string &name, double value)
262 {
263  // *****************
264  // include here all parameters which should be settable via tcl
265  // *****************
266 
267  bool knownName=true;
268 
269  if (name == "angleXMax")
270  _angleXMax=value;
271  else
272  if (name == "angleYMax")
273  _angleYMax=value;
274  else
275  if (name == "radiationLength")
276  _radiationLength=value;
277  else
278  if (name == "pmin")
279  _pmin=value;
280  else
281  if (name == "dthtc")
282  _dthtc=value;
283  else
284  if (name == "nPhotMin")
285  _nPhotMin=value;
286  else
287  if (name == "nRefrac")
288  _nRefrac=value;
289  else
290  if (name == "Bfield")
291  _Bfield=value;
292  else
293  if (name == "effNPhotons")
294  _effNPhotons=value;
295  else
296  if (name == "dRich")
297  _dRich=value;
298  else
299  if (name == "dp")
300  _dp=value;
301  else
302  if (name == "efficiency")
303  _efficiency=value;
304  else
305  knownName=false;
306 
307  return knownName;
308 }
309 
void print(std::ostream &o)
Definition: PndFsmRich.cxx:225
friend F32vec4 acos(const F32vec4 &a)
Definition: P4_F32vec4.h:113
Double_t x0
Definition: checkhelixhit.C:70
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
double _efficiency
Definition: PndFsmAbsDet.h:93
double _dRich
Definition: PndFsmRich.h:88
void setLHElectron(double val)
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
__m128 m
Definition: P4_F32vec4.h:28
TDatabasePDG * _fdbPDG
Definition: PndFsmAbsDet.h:95
void setLHProton(double val)
std::list< std::string > ArgList
Definition: ArgList.h:7
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
void setLHMuon(double val)
bool hitMapResponse(unsigned int)
Definition: PndFsmTrack.h:86
TLorentzVector s
Definition: Pnd2DStar.C:50
TString ctht(TString pts, TString exts="px py pz")
Definition: invexp.C:179
Double_t par[3]
double _nRefrac
Definition: PndFsmRich.h:85
double _dp
Definition: PndFsmRich.h:89
double charge()
Definition: PndFsmTrack.h:75
void parseParameterList(ArgList &par)
TRandom3 * _rand
Definition: PndFsmAbsDet.h:94
Double_t p
Definition: anasim.C:58
double _pmin
Definition: PndFsmRich.h:82
double _angleYMax
Definition: PndFsmRich.h:80
TLorentzVector p4()
Definition: PndFsmTrack.h:72
basic_ostream< char, char_traits< char > > ostream
void setDetector(PndFsmAbsDet *detector)
virtual PndFsmResponse * respond(PndFsmTrack *t)
Definition: PndFsmRich.cxx:92
const std::string & detName()
Definition: PndFsmAbsDet.h:74
void setLHKaon(double val)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
double _effNPhotons
Definition: PndFsmRich.h:87
TString name
bool hitMapValid()
Definition: PndFsmTrack.h:85
double _dthtc
Definition: PndFsmRich.h:83
bool setParameter(std::string &name, double value)
Definition: PndFsmRich.cxx:261
Double_t x
virtual ~PndFsmRich()
Definition: PndFsmRich.cxx:83
void setRichThtc(double val, double err=0)
std::string _detName
Definition: PndFsmAbsDet.h:92
double compThetaC(double p, double m)
Definition: PndFsmRich.cxx:189
void setDetected(bool isdet)
TTree * t
Definition: bump_analys.C:13
double gauss(double x, double x0, double s)
Definition: PndFsmRich.cxx:182
double alpha
Definition: f_Init.h:9
bool detected(PndFsmTrack *t) const
Definition: PndFsmRich.cxx:199
void initParameters()
Definition: PndFsmRich.cxx:243
double _nPhotMin
Definition: PndFsmRich.h:84
double _radiationLength
Definition: PndFsmRich.h:81
void setLHPion(double val)
double _Bfield
Definition: PndFsmRich.h:86
double _angleXMax
Definition: PndFsmRich.h:79