FairRoot/PandaRoot
PndFsmDrcDisc.cxx
Go to the documentation of this file.
1 //--------------------------------------------------------------------------
2 // File and Version Information:
3 // $Id: FsmDrcDisc.cc,v 1.7 2007/05/24 08:07:40 klausg Exp $
4 //
5 // Description:
6 // Class FsmDrcBarrel
7 //
8 // Implementation of the Disc DIRC 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) 2006 GSI
18 //
19 //------------------------------------------------------------------------
20 
21 //-----------------------
22 // This Class's Header --
23 //-----------------------
24 #include "PndFsmDrcDisc.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 #include "TH2F.h"
50 #include "TFile.h"
51 
52 //-----------------------------------------------------------------------
53 // Local Macros, Typedefs, Structures, Unions and Forward Declarations --
54 //-----------------------------------------------------------------------
55 
56 //----------------
57 // Constructors --
58 //----------------
59 
61 {
63 
64  _thtMin=_thtMin*M_PI/180.0;
65  _thtMax=_thtMax*M_PI/180.0;
66 
68 
69  //print(std::cout);
70 }
71 
73 {
75  //set default parameter values and parses a parameter list
76  //i.e. std::list<std::string> of the form
77  //"a=1" "b=2" "c=3"
78  parseParameterList(par);
79 
80  _thtMin=_thtMin*M_PI/180.0;
81  _thtMax=_thtMax*M_PI/180.0;
82 
84 
85  //print(std::cout);
86 }
87 
88 //--------------
89 // Destructor --
90 //--------------
91 
93 {
94 }
95 
96 //--------------
97 // Operations --
98 //--------------
99 
102 {
103  PndFsmResponse *result=new PndFsmResponse();
104 
105  result->setDetector(this);
106  bool wasDetected=detected(t);
107  result->setDetected(wasDetected);
108 
109  if (wasDetected && fabs(t->charge())>1e-8)
110  {
111  TParticlePDG* part = _fdbPDG->GetParticle(t->pdt());
112  double mass = (part) ? part->Mass() : t->p4().M();
113  double theta = t->p4().Theta();
114  double p=t->p4().Vect().Mag();
115  double ctht=cos(theta);
116  //double p_t=p*sin(theta);
117 
118  int pid=abs(t->pdt());
119  int npid=-1;
120 
121  if (pid==11) npid=0;
122  else if (pid==13) npid=1;
123  else if (pid==211) npid=2;
124  else if (pid==321) npid=3;
125  else if (pid==2212) npid=4;
126 
127  double thtdeg=theta/M_PI*180.;
128 
129  if ( p==0 || ctht==0 ) {result->setDetected(false);return result;} // floating point check ********************
130 
131  double lambda1 = 280e-9; //range of wavelength, which is seen by the PMT/PD
132  double lambda2 = 330e-9;
133 
134  double alpha=7.2974e-3; //finestructure constant
135 
136  // curvature of track due to magnet field
137  //double r = 3.3356 * p_t / _Bfield;
138 
139  // dip angle in phi direction (due to curvature of track in magnet field)
140  //double psi = acos(_rBarrel/(2*r));
141 
142  // path length in radiator
143  double l = _dDisc/cos(theta);// _dSlab*sqrt( 1/(sin(theta)*sin(theta)) + 1/(tan(psi)*tan(psi)) );
144 
145  // deteremine trapping fraction
146  double trapped = _trap;
147 
148  if (npid>=0 && trapfrac[npid])
149  trapped = npid<0 ? 0.0 : trapfrac[npid]->GetBinContent(trapfrac[npid]->FindBin(p<6.0?p:6.0,thtdeg));
150 
151  // estimate the number of initially produced cherenkov photons
152  double nPhot = 2*M_PI*alpha*l*(1./lambda1 - 1./lambda2)*(1 - (mass*mass+p*p)/(p*p*_nRefrac*_nRefrac));
153 
154  // dice a poisson value
155  nPhot = _rand->Poisson(nPhot);
156 
157  // determine the number of photons hitting the sensor
158  nPhot *= trapped*_effNPhotons;
159 
160  // ************** reset detected and quit due to low numbers of photons
161 
162  if (nPhot<=_nPhotMin) {
163  result->setDetected(false);
164  return result;
165  }
166 
167  if (nPhot>100) nPhot=100;
168 
169  // overall resolution for the tht_c measurement
170  double sig = _dthtc/sqrt(nPhot);
171  double thtC = compThetaC(p,mass);
172 
173  double m_e = _fdbPDG->GetParticle(11)->Mass();
174  double m_mu = _fdbPDG->GetParticle(13)->Mass();
175  double m_pi = _fdbPDG->GetParticle(211)->Mass();
176  double m_K = _fdbPDG->GetParticle(321)->Mass();
177  double m_p = _fdbPDG->GetParticle(2212)->Mass();
178 
179  // compute the expected cherenkov angles for all particle types
180  // we need these to determine the pdf's (gaussian around nominal tht_c with res sig)
181  // this Likelihood function has to be evaluated for the _measured_ momentum, which
182  // is smeared with dp!
183 
184  double measp=_rand->Gaus(p,_dp*p);
185  if (measp<0) measp=0;
186 
187  double thtc_e = compThetaC(measp,m_e);
188  double thtc_mu = compThetaC(measp,m_mu);
189  double thtc_pi = compThetaC(measp,m_pi);
190  double thtc_K = compThetaC(measp,m_K);
191  double thtc_p = compThetaC(measp,m_p);
192 
193  double measThetaC = _rand->Gaus(thtC,sig);
194  if (measThetaC<0) measThetaC=0;
195 
196  result->setDrcDiscThtc(measThetaC,sig);
197 
198  if (thtc_e) result->setLHElectron( gauss(measThetaC,thtc_e,sig) );
199  if (thtc_mu) result->setLHMuon( gauss(measThetaC,thtc_mu,sig) );
200  if (thtc_pi) result->setLHPion( gauss(measThetaC,thtc_pi,sig) );
201  if (thtc_K) result->setLHKaon( gauss(measThetaC,thtc_K,sig) );
202  if (thtc_p) result->setLHProton(gauss(measThetaC,thtc_p,sig) );
203 
204  }
205 
206  return result;
207 }
208 
209 double
210 PndFsmDrcDisc::gauss(double x, double x0, double s)
211 {
212  return (1.0/(sqrt(2.0*M_PI)*s))*
213  exp(-(x-x0)*(x-x0)/(2.0*s*s));
214 }
215 
216 double
218 {
219  double val=0.0;
220  if (p==0) return 0.0;
221  if ( (val = (p*p+m*m)/(p*p*_nRefrac*_nRefrac)) <= 1.0 ) return acos(sqrt(val));
222  else return 0.0;
223 }
224 
225 
226 bool
228 {
229  if (t->hitMapValid()) {
230  return t->hitMapResponse(FsmDetEnum::Drc);
231  } else {
232  int lundId=abs(t->pdt());
233  TParticlePDG* part = _fdbPDG->GetParticle(lundId);
234  double mass = (part) ? part->Mass() : t->p4().M();
235  double p_cerenkov_min=mass/sqrt(_nRefrac*_nRefrac - 1.0);
236  double theta = t->p4().Theta();
237  double p=t->p4().Vect().Mag();
238 
239  bool correctPidType=(lundId==11 || lundId==13 || lundId==211 || lundId==321 || lundId==2212);
240 
241  return ( p>p_cerenkov_min && theta>=_thtMin && theta<=_thtMax && correctPidType && _rand->Uniform()<=_efficiency);
242  }
243 }
244 
245 void
247 {
248  o <<"Parameters for detector <"<<detName()<<">"<<endl;
249  o <<" _thtMin = "<<_thtMin<<endl;
250  o <<" _thtMax = "<<_thtMax<<endl;
251  o <<" _radiationLength = "<<_radiationLength<<endl;
252  o <<" _pmin = "<<_pmin<<endl;
253  o <<" _dthtc = "<<_dthtc<<endl;
254  o <<" _nPhotMin = "<<_nPhotMin<<endl;
255  o <<" _nRefrac = "<<_nRefrac<<endl;
256  o <<" _Bfield = "<<_Bfield<<endl;
257  o <<" _effNPhotons = "<<_effNPhotons<<endl;
258  o <<" _dDisc = "<<_dDisc<<endl;
259  o <<" _dp = "<<_dp<<endl;
260  o <<" _trap = "<<_trap<<endl;
261  o <<" _efficiency = "<<_efficiency<<endl;
262  o <<" _parFileName = "<<_parFileName<<endl;
263 }
264 
265 void
267 {
268  _detName = "DrcDisc";
269  _thtMin = 5.0;
270  _thtMax = 22.0;
271  _radiationLength = 0.0;
272  _pmin = 0.0;
273  _dthtc = 0.01;
274  _nPhotMin = 5;
275  _nRefrac = 1.472;
276  _Bfield = 2.;
277  _effNPhotons = 0.1;
278  _dDisc = 0.017;
279  _dp = 0.01;
280  _trap = 0.7;
281  _efficiency = 1.0;
282  _parFileName = "$VMCWORKDIR/fastsim/trapfrac_disc.root";
283 }
284 
285 bool
286 PndFsmDrcDisc::setParameter(std::string &name,std::string &value)
287 {
288  // *****************
289  // include here all string parameters which should be settable
290  // *****************
291 
292  bool knownName=true;
293 
294  if (name == "parFileName")
295  _parFileName=value;
296  else
297  knownName=false;
298 
299  return knownName;
300 }
301 
302 bool
303 PndFsmDrcDisc::setParameter(std::string &name, double value)
304 {
305  // *****************
306  // include here all parameters which should be settable via tcl
307  // *****************
308 
309  bool knownName=true;
310 
311  if (name == "thtMin")
312  _thtMin=value;
313  else
314  if (name == "thtMax")
315  _thtMax=value;
316  else
317  if (name == "radiationLength")
318  _radiationLength=value;
319  else
320  if (name == "pmin")
321  _pmin=value;
322  else
323  if (name == "dthtc")
324  _dthtc=value;
325  else
326  if (name == "nPhotMin")
327  _nPhotMin=value;
328  else
329  if (name == "nRefrac")
330  _nRefrac=value;
331  else
332  if (name == "Bfield")
333  _Bfield=value;
334  else
335  if (name == "effNPhotons")
336  _effNPhotons=value;
337  else
338  if (name == "dDisc")
339  _dDisc=value;
340  else
341  if (name == "dp")
342  _dp=value;
343  else
344  if (name == "trap")
345  _trap=value;
346  else
347  if (name == "efficiency")
348  _efficiency=value;
349  else
350  knownName=false;
351 
352  return knownName;
353 }
354 
356 {
357  TFile *f=new TFile(_parFileName.c_str());
358 
359  for (int i=0;i<5;i++)
360  {
361  trapfrac[i]=0;
362  }
363 
364  if (f->IsZombie())
365  {
366  cout <<" -W- (PndFsmDrcDisc::readParameters) - file "<<_parFileName.c_str()
367  <<" doesn't exist. Using constant trapping fraction _trap="<<_trap<<endl;
368  }
369  else
370  {
371  trapfrac[0]=(TH2F*)f->Get("hacc0");
372  trapfrac[1]=(TH2F*)f->Get("hacc1");
373  trapfrac[2]=(TH2F*)f->Get("hacc2");
374  trapfrac[3]=(TH2F*)f->Get("hacc3");
375  trapfrac[4]=(TH2F*)f->Get("hacc4");
376 
377  for (int i=0;i<5;i++) trapfrac[i]->SetDirectory(0);
378 
379  f->Close();
380  }
381  delete f;
382 
383  return true;
384 }
385 
friend F32vec4 acos(const F32vec4 &a)
Definition: P4_F32vec4.h:113
Double_t x0
Definition: checkhelixhit.C:70
Double_t p
Definition: anasim.C:58
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
double _efficiency
Definition: PndFsmAbsDet.h:93
void setLHElectron(double val)
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
Int_t i
Definition: run_full.C:25
__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
int pid()
Double_t par[3]
double charge()
Definition: PndFsmTrack.h:75
void parseParameterList(ArgList &par)
bool detected(PndFsmTrack *t) const
virtual ~PndFsmDrcDisc()
TRandom3 * _rand
Definition: PndFsmAbsDet.h:94
TLorentzVector p4()
Definition: PndFsmTrack.h:72
basic_ostream< char, char_traits< char > > ostream
void setDetector(PndFsmAbsDet *detector)
void setDrcDiscThtc(double val, double err=0)
const std::string & detName()
Definition: PndFsmAbsDet.h:74
double _effNPhotons
Definition: PndFsmDrcDisc.h:95
TH2F * trapfrac[5]
Definition: PndFsmDrcDisc.h:85
std::string _parFileName
Definition: PndFsmDrcDisc.h:99
TFile * f
Definition: bump_analys.C:12
void setLHKaon(double val)
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
virtual PndFsmResponse * respond(PndFsmTrack *t)
TString name
bool hitMapValid()
Definition: PndFsmTrack.h:85
double gauss(double x, double x0, double s)
void print(std::ostream &o)
Double_t x
bool setParameter(std::string &name, double value)
std::string _detName
Definition: PndFsmAbsDet.h:92
double _nPhotMin
Definition: PndFsmDrcDisc.h:92
void setDetected(bool isdet)
TTree * t
Definition: bump_analys.C:13
double compThetaC(double p, double m)
double alpha
Definition: f_Init.h:9
double _radiationLength
Definition: PndFsmDrcDisc.h:89
void setLHPion(double val)