EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CbmRichProjectionProducer.cxx
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file CbmRichProjectionProducer.cxx
1 
9 #include "CbmRichHitProducer.h"
10 
11 #include "FairRootManager.h"
12 #include "CbmMCTrack.h"
13 #include "FairTrackParam.h"
14 #include "FairRuntimeDb.h"
15 #include "FairRun.h"
16 #include "FairGeoNode.h"
17 #include "CbmGeoRichPar.h"
18 #include "FairGeoTransform.h"
19 #include "FairGeoVector.h"
20 #include "FairRunAna.h"
21 
22 #include "TVector3.h"
23 #include "TClonesArray.h"
24 #include "TMatrixFSym.h"
25 
26 #include <iostream>
27 
28 using std::cout;
29 using std::endl;
30 
31 
33  Int_t zflag):
35  fListRICHImPlanePoint(NULL),
36 
37  fNHits(0),
38  fEvent(0),
39 
40  fDetX(0.0),
41  fDetY(0.0),
42  fDetZ(0.0),
43  fDetWidthX(0.0),
44  fDetWidthY(0.0),
45  fThetaDet(0.0),
46  fPhiDet(0.0),
47 
48  fDetXTransf(0.0),
49  fDetYTransf(0.0),
50  fDetZTransf(0.0),
51 
52  fZm(0.0),
53  fYm(0.0),
54  fXm(0.0),
55  fR(0.0),
56 
57  fMaxXTrackExtr(0.0),
58  fMaxYTrackExtr(0.0),
59 
60  fSensNodes(NULL),
61  fPassNodes(NULL),
62  fPar(NULL)
63 {
64 }
65 
67 {
69  fManager->Write();
70 }
71 
73 {
75  FairRuntimeDb* rtdb=sim->GetRuntimeDb();
76  fPar = (CbmGeoRichPar*)(rtdb->getContainer("CbmGeoRichPar"));
77 }
78 
80 {
82 
85 
86  // get detector position:
87  FairGeoNode *det= (FairGeoNode *) fSensNodes->FindObject("rich1d#1");
88  FairGeoTransform* detTr=det->getLabTransform(); // detector position in labsystem
89  FairGeoVector detPosLab=detTr->getTranslation(); // ... in cm
90  FairGeoTransform detCen=det->getCenterPosition(); // center in Detector system
91  FairGeoVector detPosCen=detCen.getTranslation();
92  fDetZ = detPosLab.Z() + detPosCen.Z(); // z coordinate of photodetector (Labsystem, cm)
93  fDetY = detPosLab.Y() + detPosCen.Y(); // y coordinate of photodetector (Labsystem, cm)
94  fDetX = detPosLab.X() + detPosCen.X(); // x coordinate of photodetector (Labsystem, cm)
95 
96  TArrayD *fdetA=det->getParameters();
97  fDetWidthX = fdetA->At(0);
98  fDetWidthY = fdetA->At(1);
99 // for(Int_t i=0;i<fdetA->GetSize();i++) cout << "Array detector " << fdetA->At(i)<< endl;
100  // detector might be rotated by theta around x-axis:
101  FairGeoRotation fdetR=detTr->getRotMatrix();
102 
103  // possible tilting around x-axis (theta) and y-axis (phi)
104  // fdetR(0) = cos(phi)
105  // fdetR(1) = 0
106  // fdetR(2) = -sin(phi)
107  // fdetR(3) = -sin(theta)sin(phi)
108  // fdetR(4) = cos(theta)
109  // fdetR(5) = -sin(theta)cos(phi)
110  // fdetR(6) = cos(theta)sin(phi)
111  // fdetR(7) = sin(theta)
112  // fdetR(8) = cos(theta)cos(phi)
113 
114  // theta = tilting angle around x-axis
115  fThetaDet = TMath::ASin(fdetR(7));
116  // phi = tilting angle around y-axis
117  fPhiDet = -1.*TMath::ASin(fdetR(2));
118 
119  cout << "---------------------- RICH Projection Producer ---------------------------------------" << endl;
120  cout << " detector position in (x,y,z): " << fDetX << " " << fDetY << " " << fDetZ << endl;
121  cout << " detector size in x and y: " << fDetWidthX << " " << fDetWidthY << endl;
122  cout << " detector tilting angle (around x): " << fThetaDet*180./TMath::Pi() << " degrees" << endl;
123  cout << " detector tilting angle (around y): " << fPhiDet*180./TMath::Pi() << " degrees" << endl;
124 
125  // transform nominal detector position (for tilted photodetector):
126  // shift x back by fDetZ_org*TMath::Sin(phi) in order to avoid overlap
127  fDetXTransf = fDetX*TMath::Cos(fPhiDet)+fDetZ*TMath::Sin(fPhiDet)-fDetZ*TMath::Sin(fPhiDet);
128  fDetYTransf = -fDetX*TMath::Sin(fThetaDet)*TMath::Sin(fPhiDet) + fDetY*TMath::Cos(fThetaDet) + fDetZ*TMath::Sin(fThetaDet)*TMath::Cos(fPhiDet);
129  fDetZTransf = -fDetX*TMath::Cos(fThetaDet)*TMath::Sin(fPhiDet) - fDetY*TMath::Sin(fThetaDet) + fDetZ*TMath::Cos(fThetaDet)*TMath::Cos(fPhiDet);
130 
131  // get mirror position:
132  //FairGeoNode *mir= (FairGeoNode *) fPassNodes->FindObject("rich1mgl#1");
133  FairGeoNode *mir= (FairGeoNode *) fSensNodes->FindObject("rich1mgl#1");
134  FairGeoTransform* mirTr=mir->getLabTransform(); // position of mirror center in labsystem
135  FairGeoVector mirPosLab=mirTr->getTranslation(); // ... in cm
136  fZm = mirPosLab.Z();
137  fYm = mirPosLab.Y();
138  fXm = mirPosLab.X();
139 
140  TArrayD *fmirA=mir->getParameters(); // get other geometry parameters: radius,
141  fR = fmirA->At(0); // mirror radius
142  Double_t spheTheta = TMath::Abs(90. - fmirA->At(2)); // opening angle for SPHERE in theta (90 degree +- theta)
143  Double_t sphePhi = TMath::Abs(90. - fmirA->At(4)); // opening angle for SPHERE in phi (90 degree +- phi)
144  // from that calculate (with safety factor 1.3) maximum x-y positions for track extrapolation:
145  fMaxXTrackExtr = 1.3*(fR*TMath::Tan(sphePhi*TMath::Pi()/180.));
146  fMaxYTrackExtr = 1.3*(TMath::Abs(fYm) + fR*TMath::Tan(spheTheta*TMath::Pi()/180.));
147 
148  // mirror might be rotated by theta around x-axis:
149  FairGeoRotation fmirR=mirTr->getRotMatrix();
150  Double_t thetaM = -1.*TMath::ASin(fmirR(5)) - TMath::Pi()/2 ;
151 
152  // note that mirror is by default tilted by 90 degrees in order to get the necessary shape in GEANT
153  // the "extra" tilting angle is then: fThetaM = -1.*TMath::ASin(fmirR(5)) - TMath::Pi()/2.
154  cout << "Mirror center (x,y,z): " << fXm << " " << fYm << " " << fZm << endl;
155  cout << "Mirror radius: " << fR << endl;
156  cout << "Mirror tilting angle: " << thetaM*180./TMath::Pi() << " degrees" << endl;
157 
158  fEvent = 0;
159 
160  fListRICHImPlanePoint = (TClonesArray*)fManager->GetObject("RichTrackParamZ");
161  if (fZflag == 1) cout << " use tracks in imaginary plane for projection to photodetector plane" << endl;
162  if (fZflag == 2) cout << " use tracks in RICH mirror for projection to photodetector plane" << endl;
163  if ( NULL == fListRICHImPlanePoint) {
164  cout << "-W- CbmRichProjectionProducer::Init: No Rich Z-Point array!" << endl;
165  }
166 }
167 
169  TClonesArray* richProj)
170 {
171  fEvent++;
172  cout << "CbmRichProjectionProducer: event " << fEvent << endl;
173 
174  richProj->Clear();
175  TMatrixFSym covMat(5);
176  for(Int_t i = 0; i < 5; i++){
177  for(Int_t j=0; j<=i; j++){
178  covMat(i,j) = 0;
179  }
180  }
181  covMat(0,0) = covMat(1,1) = covMat(2,2) = covMat(3,3) = covMat(4,4) = 1.e-4;
182 
183  for(Int_t j = 0; j < fListRICHImPlanePoint->GetEntriesFast(); j++) {
185  new((*richProj)[j]) FairTrackParam(0., 0., 0., 0., 0., 0., covMat);
186 
187  // check if Array was filled
188  if (point->GetX() == 0 && point->GetY() == 0 && point->GetZ() == 0 &&
189  point->GetTx() == 0 && point->GetTy() ==0) continue;
190  if (point->GetQp()==0) continue;
191 
192  // check that x and y value make sense (sometimes strange extrapolations may appear)
193  //if (TMath::Abs(point->GetX()) > fMaxXTrackExtr || TMath::Abs(point->GetY()) > fMaxYTrackExtr){
194  // cout << " -W- RichProjectionProducer: strange (x,y) values for track extrapolation: " <<
195  // point->GetX() << " " << point->GetY() << endl;
196  // continue;
197  //}
198 
199  Double_t rho1 = 0.;
200  Double_t rho2 = 0.;
201  TVector3 startP, momP, crossP, centerP;
202 
203  // operate on ImPlane point
204  if (fZflag ==1) {
205  Double_t p = 1./TMath::Abs(point->GetQp());
206  Double_t pz;
207  if ((1+point->GetTx()*point->GetTx()+point->GetTy()*point->GetTy()) > 0. )
208  pz = p/TMath::Sqrt(1+point->GetTx()*point->GetTx()+point->GetTy()*point->GetTy());
209  else {
210  cout << " -E- RichProjectionProducer: strange value for calculating pz: " <<
211  (1+point->GetTx()*point->GetTx()+point->GetTy()*point->GetTy()) << endl;
212  pz = 0.;
213  }
214  Double_t px = pz*point->GetTx();
215  Double_t py = pz*point->GetTy();
216  momP.SetXYZ(px,py,pz);
217  point->Position(startP);
218  if ((fYm*startP.y())<0) fYm = -fYm; // check that mirror center and startP are in same hemisphere
219 
220  // calculation of intersection of track with selected mirror
221  // corresponds to calculation of intersection between a straight line and a sphere:
222  // vector: r = startP - mirrorCenter
223  // RxP = r*momP
224  // normP2 = momP^2
225  // dist = r^2 - fR^2
226  // -> rho1 = (-RxP+sqrt(RxP^2-normP2*dist))/normP2 extrapolation factor for:
227  // intersection point crossP = startP + rho1 * momP
228  Double_t RxP=(momP.x()*(startP.x()-fXm)+momP.y()*(startP.y()-fYm)+momP.z()*(startP.z()-fZm));
229  Double_t normP2=(momP.x()*momP.x()+momP.y()*momP.y()+momP.z()*momP.z());
230  Double_t dist=(startP.x()*startP.x()+fXm*fXm+startP.y()*startP.y()+fYm*fYm+startP.z()*startP.z()+fZm*fZm-2*startP.x()*fXm-2*startP.y()*fYm-2*startP.z()*fZm-fR*fR);
231 
232  if ((RxP*RxP-normP2*dist) > 0.) {
233  if (normP2!=0.) rho1=(-RxP+TMath::Sqrt(RxP*RxP-normP2*dist))/normP2;
234  if (normP2 == 0) cout << " Error in track extrapolation: momentum = 0 " << endl;
235  } else {
236  cout << " -E- RichProjectionProducer: RxP*RxP-normP2*dist = " << RxP*RxP-normP2*dist << endl;
237  }
238 
239  Double_t crossPx = startP.x() + rho1*momP.x();
240  Double_t crossPy = startP.y() + rho1*momP.y();
241  Double_t crossPz = startP.z() + rho1*momP.z();
242  crossP.SetXYZ(crossPx, crossPy, crossPz);
243 
244  // check if crosspoint with mirror and chosen mirrorcenter (y) are in same hemisphere
245  // if not recalculate crossing point
246  if ((fYm*crossP.y())<0) {
247  fYm = -fYm;
248  RxP=(momP.x()*(startP.x()-fXm)+momP.y()*(startP.y()-fYm)+momP.z()*(startP.z()-fZm));
249  normP2=(momP.x()*momP.x()+momP.y()*momP.y()+momP.z()*momP.z());
250  dist=(startP.x()*startP.x()+fXm*fXm+startP.y()*startP.y()+fYm*fYm+startP.z()*startP.z()+fZm*fZm-2*startP.x()*fXm-2*startP.y()*fYm-2*startP.z()*fZm-fR*fR);
251 
252  if ((RxP*RxP-normP2*dist) > 0.) {
253  if (normP2!=0.) rho1=(-RxP+TMath::Sqrt(RxP*RxP-normP2*dist))/normP2;
254  if (normP2 == 0) cout << " Error in track extrapolation: momentum = 0 " << endl;
255  } else{
256  cout << " -E- RichProjectionProducer: RxP*RxP-normP2*dist = " << RxP*RxP-normP2*dist << endl;
257  }
258 
259  crossPx=startP.x()+rho1*momP.x();
260  crossPy=startP.y()+rho1*momP.y();
261  crossPz=startP.z()+rho1*momP.z();
262  crossP.SetXYZ(crossPx,crossPy,crossPz);
263  }
264 
265  centerP.SetXYZ(fXm,fYm,fZm); // mirror center
266  }// if (fZflag ==1)
267 
268  // operate on Rich Mirror point
269  if (fZflag ==2) {
270  Double_t p = 1./TMath::Abs(point->GetQp());
271  Double_t pz;
272  if ((1+point->GetTx()*point->GetTx()+point->GetTy()*point->GetTy()) > 0. ){
273  pz = p/TMath::Sqrt(1+point->GetTx()*point->GetTx()+point->GetTy()*point->GetTy());
274  } else {
275  cout << " -E- RichProjectionProducer: strange value for calculating pz: " <<
276  (1+point->GetTx()*point->GetTx()+point->GetTy()*point->GetTy()) << endl;
277  pz = 0.;
278  }
279  Double_t px = pz*point->GetTx();
280  Double_t py = pz*point->GetTy();
281  momP.SetXYZ(px,py,pz);
282  point->Position(crossP);
283  if ((fYm*crossP.y())<0) fYm = -fYm; // check that mirror center and crossP are in same hemisphere
284 
285  centerP.SetXYZ(fXm,fYm,fZm); // mirror center
286  } // if (fZflag ==2)
287 
288  // calculate normal on crosspoint with mirror
289  TVector3 normP(crossP.x()-centerP.x(),crossP.y()-centerP.y(),crossP.z()-centerP.z());
290  normP=normP.Unit();
291  // check that normal has same z-direction as momentum
292  if ((normP.z()*momP.z())<0.) normP = TVector3(-1.*normP.x(),-1.*normP.y(),-1.*normP.z());
293 
294  // reflect track
295  Double_t np=normP.x()*momP.x()+normP.y()*momP.y()+normP.z()*momP.z();
296 
297  Double_t refX = 2*np*normP.x()-momP.x();
298  Double_t refY = 2*np*normP.y()-momP.y();
299  Double_t refZ = 2*np*normP.z()-momP.z();
300 
301  // crosspoint whith photodetector plane:
302  // calculate intersection between straight line and (tilted) plane:
303  // normal on plane tilted by theta around x-axis: (0,-sin(theta),cos(theta)) = n
304  // normal on plane tilted by phi around y-axis: (-sin(phi),0,cos(phi)) = n
305  // normal on plane tilted by theta around x-axis and phi around y-axis: (-sin(phi),-sin(theta)cos(phi),cos(theta)cos(phi)) = n
306  // point on plane is (fDetX,fDetY,fDetZ) = p as photodetector is tiled around its center
307  // equation of plane for r being point in plane: n(r-p) = 0
308  // calculate intersection point of reflected track with plane: r=intersection point
309  // intersection point = crossP + rho2 * refl_track
310  // take care for all 4 cases:
311  // -> first calculate for case x>0, then check
312  if (refZ!=0.) {
313  if (centerP.y() > 0){
314  rho2 = (-TMath::Sin(fPhiDet)*(fDetX-crossP.x())
315  -TMath::Sin(fThetaDet)*TMath::Cos(fPhiDet)*(fDetY-crossP.y())
316  + TMath::Cos(fThetaDet)*TMath::Cos(fPhiDet)*(fDetZ-crossP.z()))/
317  (-TMath::Sin(fPhiDet)*refX-TMath::Sin(fThetaDet)*TMath::Cos(fPhiDet)*refY + TMath::Cos(fThetaDet)*TMath::Cos(fPhiDet)*refZ);
318  }
319  if (centerP.y() < 0){
320  rho2 = (-TMath::Sin(fPhiDet)*(fDetX-crossP.x())
321  -TMath::Sin(-fThetaDet)*TMath::Cos(fPhiDet)*(-fDetY-crossP.y())
322  + TMath::Cos(-fThetaDet)*TMath::Cos(fPhiDet)*(fDetZ-crossP.z()))/
323  (-TMath::Sin(fPhiDet)*refX-TMath::Sin(-fThetaDet)*TMath::Cos(fPhiDet)*refY + TMath::Cos(-fThetaDet)*TMath::Cos(fPhiDet)*refZ);
324  }
325 
326  //rho2 = -1*(crossP.z() - fDetZ)/refZ; // only for theta = 0, phi=0
327  Double_t xX = crossP.x() + refX * rho2;
328  Double_t yY = crossP.y() + refY * rho2;
329  Double_t zZ = crossP.z() + refZ * rho2;
330 
331  if (xX < 0) {
332  if (centerP.y() > 0){
333  rho2 = (-TMath::Sin(-fPhiDet)*(-fDetX-crossP.x())
334  -TMath::Sin(fThetaDet)*TMath::Cos(-fPhiDet)*(fDetY-crossP.y())
335  + TMath::Cos(fThetaDet)*TMath::Cos(-fPhiDet)*(fDetZ-crossP.z()))/
336  (-TMath::Sin(-fPhiDet)*refX-TMath::Sin(fThetaDet)*TMath::Cos(-fPhiDet)*refY + TMath::Cos(fThetaDet)*TMath::Cos(-fPhiDet)*refZ);
337  }
338  if (centerP.y() < 0){
339  rho2 = (-TMath::Sin(-fPhiDet)*(-fDetX-crossP.x())
340  -TMath::Sin(-fThetaDet)*TMath::Cos(-fPhiDet)*(-fDetY-crossP.y())
341  + TMath::Cos(-fThetaDet)*TMath::Cos(-fPhiDet)*(fDetZ-crossP.z()))/
342  (-TMath::Sin(-fPhiDet)*refX-TMath::Sin(-fThetaDet)*TMath::Cos(-fPhiDet)*refY + TMath::Cos(-fThetaDet)*TMath::Cos(-fPhiDet)*refZ);
343  }
344 
345  xX = crossP.x() + refX * rho2;
346  yY = crossP.y() + refY * rho2;
347  zZ = crossP.z() + refZ * rho2;
348  }
349 
350  // Transform intersection point in same way as MCPoints were
351  // transformed in HitProducer before stored as Hit:
352  TVector3 inPos(xX, yY, zZ);
353  TVector3 outPos;
355  Double_t xDet = outPos.X();
356  Double_t yDet = outPos.Y();
357  Double_t zDet = outPos.Z();
358  //printf("%f %f %f\n", xDet, yDet, zDet);
359 
360  //check that crosspoint inside the plane
361  if( xDet > (-fDetX-fDetWidthX) && xDet < (fDetX+fDetWidthX)){
362  if(TMath::Abs(yDet) > (fDetYTransf-fDetWidthY) && TMath::Abs(yDet) < (fDetYTransf+fDetWidthY)){
363  FairTrackParam richtrack(xDet,yDet,zDet,0.,0.,0.,covMat);
364  * (FairTrackParam*)(richProj->At(j)) = richtrack;
365  }
366  }
367  }// if (refZ!=0.)
368  }// j
369 }