FairRoot/PandaRoot
PndFsmEmcFwCap.cxx
Go to the documentation of this file.
1 //--------------------------------------------------------------------------
2 // File and Version Information:
3 // $Id: FsmEmcFwCap.cc,v 1.9 2007/05/24 08:07:40 klausg Exp $
4 //
5 // Description:
6 // Class FsmEmcFwCap
7 //
8 // Implementation of the EMC Barrel part 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 "PndFsmEmcFwCap.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 #include "ArgList.h"
45 #include "PndFsmResponse.h"
46 #include "PndFsmTrack.h"
47 
48 
49 //-----------------------------------------------------------------------
50 // Local Macros, Typedefs, Structures, Unions and Forward Declarations --
51 //-----------------------------------------------------------------------
52 
53 //----------------
54 // Constructors --
55 //----------------
56 
58 {
60 
61  _thtMin=_thtMin*M_PI/180.0;
62  _thtMax=_thtMax*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  _thtMin=_thtMin*M_PI/180.0;
75  _thtMax=_thtMax*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  result->setdE(dE(t));
103  result->setdphi(dphi(t));
104  result->setdtheta(dtheta(t));
105  double zscale = t->p4().Pz()/_dist;
106  double x = t->p4().Px()/zscale;
107  double y = t->p4().Py()/zscale;
108  TVector3 hitpos(x,y,_dist);
109  t->setStopVtx(hitpos);
110  }
111  else
112  {
113  result->setdE(0.);
114  result->setdphi(0.);
115  result->setdtheta(0.);
116  }
117 
118  return result;
119 }
120 
121 bool
123 {
124  if (t->hitMapValid()) {
126  } else {
127  double theta = t->p4().Theta();
128  double E=t->p4().E();
129  double lund = t->pdt();
130  return ( lund==22 && theta>=_thtMin && theta<=_thtMax && E>_Emin && _rand->Rndm()<=_efficiency);
131  }
132 }
133 
134 double
136 {
137  double E = t->p4().E();
138 
139  return (_aPar+_bPar/E+_cPar/sqrt(E) ) * E;
140 }
141 
142 double
144 {
145  double theta = t->p4().Vect().Theta();
146 
147  return (_resFactor*M_PI/int(2*M_PI*_dist*tan(theta)/_xtalDim) );
148 }
149 
150 double
152 {
153  double theta = t->p4().Vect().Theta();
154  return ( _resFactor*atan(_xtalDim*cos(theta)/(2*_dist)) );
155 }
156 
157 void
159 {
160  o <<"Detector <"<<_detName<<">"<<endl;
161  o <<" _aPar = "<<_aPar<<endl;
162  o <<" _bPar = "<<_bPar<<endl;
163  o <<" _cPar = "<<_cPar<<endl;
164  o <<" _xtalDim = "<<_xtalDim<<endl;
165  o <<" _Emin = "<<_Emin<<endl;
166  o <<" _dist = "<<_dist<<endl;
167  o <<" _resFactor = "<<_resFactor<<endl;
168  o <<" _thtMin = "<<_thtMin<<endl;
169  o <<" _thtMax = "<<_thtMax<<endl;
170  o <<" _radiationLength = "<<_radiationLength<<endl;
171  o <<" _efficiency = "<<_efficiency<<endl;
172 }
173 
174 void
176 {
177  _detName = "EmcFwCap";
179  _aPar=4.52495e-3;
180  _bPar=2.9539e-3;
181  _cPar=7.7596e-3;
182  //_aPar = 0.01;
183  //_bPar = 0.01;
184  _xtalDim = 0.02;
185  _Emin = 0.01;
186  _dist = 2.5;
187  _resFactor = 0.25;
188  _thtMin = 5.0;
189  _thtMax = 22.0;
190  _radiationLength = 0.0;
191  _efficiency =1.0;
192 }
193 
194 bool
195 PndFsmEmcFwCap::setParameter(std::string &name, double value)
196 {
197  // *****************
198  // include here all parameters which should be settable via tcl
199  // *****************
200 
201  bool knownName=true;
202 
203  if (name == "aPar")
204  _aPar=value;
205  else
206  if (name == "bPar")
207  _bPar=value;
208  else
209  if (name == "cPar")
210  _cPar=value;
211  else
212  if (name == "xtalDim")
213  _xtalDim=value;
214  else
215  if (name == "Emin")
216  _Emin=value;
217  else
218  if (name == "dist")
219  _dist=value;
220  else
221  if (name == "resFactor")
222  _resFactor=value;
223  else
224  if (name == "thtMin")
225  _thtMin=value;
226  else
227  if (name == "thtMax")
228  _thtMax=value;
229  else
230  if (name == "radiationLength")
231  _radiationLength=value;
232  else
233  if (name == "efficiency")
234  _efficiency=value;
235  else
236  knownName=false;
237 
238  return knownName;
239 }
double _radiationLength
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
double _efficiency
Definition: PndFsmAbsDet.h:93
void setdphi(double val)
bool detected(PndFsmTrack *t) const
bool setParameter(std::string &name, double value)
void setStopVtx(TVector3 v)
std::list< std::string > ArgList
Definition: ArgList.h:7
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
bool hitMapResponse(unsigned int)
Definition: PndFsmTrack.h:86
Double_t par[3]
double charge()
Definition: PndFsmTrack.h:75
void parseParameterList(ArgList &par)
TRandom3 * _rand
Definition: PndFsmAbsDet.h:94
virtual PndFsmResponse * respond(PndFsmTrack *t)
virtual ~PndFsmEmcFwCap()
static const std::string & name(unsigned int t)
Definition: FsmDetTypes.h:26
TLorentzVector p4()
Definition: PndFsmTrack.h:72
basic_ostream< char, char_traits< char > > ostream
void setDetector(PndFsmAbsDet *detector)
double dphi(PndFsmTrack *t) const
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
TString name
bool hitMapValid()
Definition: PndFsmTrack.h:85
void setdE(double val)
Double_t x
std::string _detName
Definition: PndFsmAbsDet.h:92
void setDetected(bool isdet)
TTree * t
Definition: bump_analys.C:13
void setdtheta(double val)
Double_t y
double dE(PndFsmTrack *t) const
void print(std::ostream &o)
double dtheta(PndFsmTrack *t) const