Fork me on GitHub

source: svn/trunk/src/VeryForward.cc@ 439

Last change on this file since 439 was 411, checked in by severine ovyn, 15 years ago

remove segmentation

File size: 14.2 KB
Line 
1/***********************************************************************
2** **
3** /----------------------------------------------\ **
4** | Delphes, a framework for the fast simulation | **
5** | of a generic collider experiment | **
6** \------------- arXiv:0903.2225v1 ------------/ **
7** **
8** **
9** This package uses: **
10** ------------------ **
11** FastJet algorithm: Phys. Lett. B641 (2006) [hep-ph/0512210] **
12** Hector: JINST 2:P09005 (2007) [physics.acc-ph:0707.1198v2] **
13** FROG: [hep-ex/0901.2718v1] **
14** **
15** ------------------------------------------------------------------ **
16** **
17** Main authors: **
18** ------------- **
19** **
20** Severine Ovyn Xavier Rouby **
21** severine.ovyn@uclouvain.be xavier.rouby@cern **
22** **
23** Center for Particle Physics and Phenomenology (CP3) **
24** Universite catholique de Louvain (UCL) **
25** Louvain-la-Neuve, Belgium **
26** **
27** Copyright (C) 2008-2009, **
28** All rights reserved. **
29** **
30***********************************************************************/
31
32#include "VeryForward.h"
33#include "PdgParticle.h"
34#include "H_RomanPot.h"
35
36#include <iostream>
37#include <fstream>
38#include<cmath>
39
40using namespace std;
41
42/* Notes on the correct initialisation for Hector
43 * -- these notes apply to the LHC beamlines
44 *
45 * beam1 : forward direction is for increasing 's' values
46 * beamline1 = new H_BeamLine(1,...);
47 * beamline1->fill(DET->RP_beam1Card,1,DET->RP_IP_name);
48 * beam2 : forward direction is for decreasing 's' values
49 * beamline2 = new H_BeamLine(-1,...);
50 * beamline2->fill(DET->RP_beam2Card,-1,DET->RP_IP_name);
51 *
52 * relative_energy should be false -- kickers_on should be 1
53 *
54 */
55
56//------------------------------------------------------------------------------
57VeryForward::VeryForward() :
58 DET(new RESOLution()), d_max(1.+std::max(DET->RP_420_s,DET->RP_220_s)),
59 beamline1(new H_BeamLine(1,d_max)), beamline2(new H_BeamLine(-1,d_max)),
60 rel_energy(true), // should always be true
61 kickers(1) // should always be 1
62 {
63 init(); //Initialisation of Hector
64}
65
66VeryForward::VeryForward(const string& DetDatacard) :
67 DET(new RESOLution())
68 {
69 DET->ReadDataCard(DetDatacard);
70 const float d_max = 1.+std::max(DET->RP_420_s,DET->RP_220_s);
71 beamline1 = new H_BeamLine(1,d_max);
72 beamline2 = new H_BeamLine(-1,d_max);
73 init(); //Initialisation of Hector
74 rel_energy = true; // should always be true
75 kickers = 1; // should always be 1
76}
77
78VeryForward::VeryForward(const RESOLution * DetDatacard) :
79 DET(new RESOLution(*DetDatacard)), d_max(1.+std::max(DET->RP_420_s,DET->RP_220_s)),
80 beamline1(new H_BeamLine(1,d_max)), beamline2(new H_BeamLine(-1,d_max)),
81 rel_energy(true), // should always be true
82 kickers(1) // should always be 1
83 {
84 init(); //Initialisation of Hector
85}
86
87VeryForward::VeryForward(const VeryForward& vf) :
88 DET(new RESOLution(*(vf.DET))), d_max(vf.d_max),
89 beamline1(new H_BeamLine(*(vf.beamline1))), beamline2(new H_BeamLine(*(vf.beamline2))),
90 rel_energy(vf.rel_energy),
91 kickers(vf.kickers) {
92}
93
94VeryForward& VeryForward::operator=(const VeryForward& vf){
95 if (this==&vf) return *this;
96 DET = new RESOLution(*(vf.DET));
97 d_max = vf.d_max;
98 beamline1 = new H_BeamLine(*(vf.beamline1));
99 beamline2 = new H_BeamLine(*(vf.beamline2));
100 rel_energy =vf.rel_energy;
101 kickers = vf.kickers;
102 return *this;
103}
104
105
106void VeryForward::init() {
107 //Initialisation of Hector
108 static unsigned int counter;
109 counter =0;
110 extern bool relative_energy;
111 extern int kickers_on;
112 relative_energy = rel_energy; // should always be true
113 kickers_on = kickers; // should always be 1
114 beamline1->fill(DET->RP_beam1Card,1,DET->RP_IP_name);
115 beamline1->offsetElements(DET->RP_offsetEl_s,-DET->RP_offsetEl_x); // relative energy: does not change anything
116 H_RomanPot * rp220_1 = new H_RomanPot("rp220_1",DET->RP_220_s,DET->RP_220_x*1E6);
117 // RP 220m, 2mm, beam 1
118 H_RomanPot * rp420_1 = new H_RomanPot("rp420_1",DET->RP_420_s,DET->RP_420_x*1E6);
119 // RP 420m, 4mm, beam 1
120 //rp220_1->printProperties();
121 //rp420_1->printProperties();
122 beamline1->add(rp220_1);
123 beamline1->add(rp420_1);
124
125 beamline2->fill(DET->RP_beam2Card,-1,DET->RP_IP_name);
126 beamline2->offsetElements(DET->RP_offsetEl_s,-DET->RP_offsetEl_x); // relative energy: does not change anything
127 H_RomanPot * rp220_2 = new H_RomanPot("rp220_2",DET->RP_220_s,DET->RP_220_x*1E6);
128 // RP 220m, 2mm, beam 2
129 H_RomanPot * rp420_2 = new H_RomanPot("rp420_2",DET->RP_420_s,DET->RP_420_x*1E6);
130 // RP 420m, 4mm, beam 2
131 //rp220_2->printProperties();
132 //rp420_2->printProperties();
133 beamline2->add(rp220_2);
134 beamline2->add(rp420_2);
135 // rp220_1, rp220_2, rp420_1 and rp420_2 will be deallocated in ~H_AbstractBeamLine
136 // do not put explicit delete
137}
138
139
140float VeryForward::time_of_flight(TRootGenParticle *particle, const float detector_s, const float detector_etamin, const float detector_t_resolution) {
141 // time of flight t is t = T + d/[ cos(theta) v ]
142 float cos_theta = 1; //very good approximation, if detector_etamin >3
143 if (detector_etamin<3) { // if smaller eta -> make the complete calculation
144 double tx = atan(particle->Px/particle->Pz);
145 double ty = atan(particle->Py/particle->Pz);
146 double theta = sqrt( pow(tx,2) + pow(ty,2) );
147 // cout << "tx = " << tx << " ty = " << ty << " theta = " << theta << " cos(theta) = " << cos(theta) << endl;
148 // NB: in practice, eta= 8 <-> theta 0.038° <-> 7x10^-4 rad <-> cos(theta) ~1
149 // eta = 2.6 <-> cos(theta) = 0.99
150 // eta = 3.0 <-> cos(theta) = 0.995
151 cos_theta = cos(theta);
152 }
153 // units from StdHEP : Z [mm] T[mm/c]
154 // units from Delphes : VFD_s_zdc [m] speed_of_light [m/s]
155 double flight_distance = (detector_s - particle->Z*(1E-3))/cos_theta ;
156 // assumes highly relativistic particles
157 double flight_time = (flight_distance + 1E-3 * particle->T )/speed_of_light;
158 double timeS = gRandom->Gaus(flight_time,detector_t_resolution);
159 return timeS;
160}
161
162
163
164void VeryForward::ZDC(ExRootTreeWriter *treeWriter, ExRootTreeBranch *branchZDC, TRootGenParticle *particle)
165{
166 TRootZdcHits *elementZdc;
167 float energy = particle->E;
168
169 // Zero degree calorimeter, for forward neutrons and photons
170 if (particle->Status ==1 && ( (particle->PID==pN && energy>DET->ZDC_n_E) ||
171 (particle->PID==pGAMMA && energy>DET->ZDC_gamma_E) )
172 && fabs(particle->Eta) > DET->VFD_min_zdc ) {
173 elementZdc = (TRootZdcHits*) branchZDC->NewEntry();
174
175 TLorentzVector genMomentum;
176 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
177 elementZdc->Set(genMomentum); // initialises the gen-level data
178 elementZdc->pid = particle->PID;
179
180 // 1) energy smearing
181 float energyS = -1.;
182 if (particle->PID == pGAMMA)
183 energyS = gRandom->Gaus(particle->E, sqrt( pow(DET->ELG_Nzdc,2) +
184 pow(DET->ELG_Czdc*particle->E,2) +
185 pow(DET->ELG_Szdc*sqrt(particle->E),2) ));
186 else // smearing with hadronic resolution
187 energyS = gRandom->Gaus(particle->E, sqrt( pow(DET->HAD_Nzdc,2) +
188 pow(DET->HAD_Czdc*particle->E,2) +
189 pow(DET->HAD_Szdc*sqrt(particle->E),2) ));
190 elementZdc->E = energyS;
191
192 // 2) time of flight t is t = T + d/[ cos(theta) v ] + detector smearing on time
193 elementZdc->T = time_of_flight(particle, DET->VFD_s_zdc, DET->VFD_min_zdc, DET->ZDC_T_resolution);
194
195 // 3) side: which ZDC has been hit?
196 elementZdc->side = sign(particle->Eta);
197
198 // 4) object nature : e.m. (photon) or had (neutron) ?
199 elementZdc->hadronic_hit = (bool) (particle->PID!=pGAMMA);
200 } // if neutrons or photons over E_threshold
201
202}
203
204
205void VeryForward::RomanPots(ExRootTreeWriter *treeWriter, ExRootTreeBranch *branchRP220,ExRootTreeBranch *branchFP420,TRootGenParticle *particle)
206{
207 if(particle->Status != 1) return; // reject particles that are not final ones
208 extern bool relative_energy;
209 relative_energy = rel_energy;
210 extern int kickers_on;
211 kickers_on = kickers;
212
213 float charge = particle->Charge, mass = particle->M;
214 //float charge, mass, ctau;
215 //charge = mass = ctau = UNDEFINED;
216 if (mass<-999) { // unitialised!
217 PdgParticle pdg_part = DET->PDGtable[particle->PID];
218 charge = pdg_part.charge(); // e+
219 mass = pdg_part.mass(); // GeV
220 // ctau = pdg_part.ctau(); // m
221 // cout << "ctau = " << ctau << endl;
222 }
223
224
225 if(particle->Charge==0) return; // only particles with Q=+1 can hope to reach RP200/FP420
226 //cout << "particle ("<< particle->PID << "): m = " << mass << " \t Q= " << charge << endl;
227 TRootRomanPotHits* elementRP220;
228 //TRootForwardTaggerHits* elementFP420;
229 TRootRomanPotHits* elementFP420;
230
231 TLorentzVector genMomentum;
232 genMomentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
233
234 // to go faster, why not rejecting particles already going into the ZDC?
235
236 // K_L^0 has a ctau of ~16m ; only pi+ and p+ can beat it ;
237 // so if RP/FP too far away, the particle must be a proton or a mu+
238 if( std::min(DET->RP_420_s,DET->RP_220_s) > 17 && (particle->PID != pP && particle->PID != 13)) return;
239
240 if( fabs(genMomentum.Eta()) > DET->CEN_max_calo_fwd )
241 {
242 H_BeamParticle p1(mass,charge);
243 double tx = 1E6*atan(particle->Px/particle->Pz); // in microrad
244 double ty = 1E6*atan(particle->Py/particle->Pz); // in microrad
245 //p1.smearAng(); p1.smearPos(); // vertex smearing do not put it here !!!
246
247 /* cout << "x = " << particle->X << " + " << p1.getX() << " + " << DET->RP_cross_x
248 << " y= " << particle->Y << " + " << p1.getY() << " + " << DET->RP_cross_y
249 << " tx= " << tx << " + " << p1.getTX() << " - " << kickers_on*DET->RP_cross_ang_x
250 << " ty= " << ty << " + " << p1.getTY() << " - " << kickers_on*DET->RP_cross_ang_y
251 << " z= " << particle->Z << endl;*/
252
253 // here below, p1.getX(), p1.getY(), p1.getTX() and p1.getTY() =0 unless some smearing is done
254 // all in micrometers or microradians
255 p1.setPosition((1E3)*particle->X + p1.getX() + DET->RP_cross_x,
256 (1E3)*particle->Y + p1.getY() + DET->RP_cross_y,
257 tx + p1.getTX()- kickers_on*DET->RP_cross_ang_x,
258 ty + p1.getTY()- kickers_on*DET->RP_cross_ang_y,
259 0*(1E3)*particle->Z);
260 p1.setE(particle->E);
261
262 H_BeamLine *beamline;
263 if(genMomentum.Eta() >0) beamline = beamline1;
264 else beamline = beamline2;
265
266 p1.computePath(beamline,1);
267
268 if(p1.stopped(beamline)) {
269 // roman pots at 220 m
270 if (p1.getStoppingElement()->getName()=="rp220_1" || p1.getStoppingElement()->getName()=="rp220_2") {
271
272 /*static unsigned int counter;
273 counter++;
274 if (counter==1) {
275 p1.getPath(0,"p1path.txt");
276 cout << "RP : " << particle->PID << "\t" << charge << "=" << particle->Charge
277 << "\t" << mass << "=" << particle->M << "\t E=" << particle->E << endl;
278 }*/
279
280 p1.propagate(DET->RP_220_s);
281 elementRP220 = (TRootRomanPotHits*) branchRP220->NewEntry();
282
283 // detector measurements
284 elementRP220->X = (1E-6)*p1.getX(); // [m]
285 elementRP220->Y = (1E-6)*p1.getY(); // [m]
286 elementRP220->Tx = (1E-6)*p1.getTX(); // [rad]
287 elementRP220->Ty = (1E-6)*p1.getTY(); // [rad]
288 elementRP220->S = p1.getS(); // [m]
289 elementRP220->T = time_of_flight(particle, DET->RP_220_s, DET->CEN_max_calo_fwd, DET->RP220_T_resolution);
290 elementRP220->side = sign(particle->Eta);
291
292 // reconstructed data
293 float sE = p1.getE(); // apply the smearing here!!!
294 elementRP220->E = sE; // not yet implemented
295 elementRP220->q2 = UNDEFINED; // not yet implemented
296
297 // generator level data
298 elementRP220->pid = particle->PID;
299 elementRP220->Set(genMomentum);
300 } // if RP220
301
302 // proton taggers at 420 m
303 else if (p1.getStoppingElement()->getName()=="rp420_1" || p1.getStoppingElement()->getName()=="rp420_2") {
304
305 p1.propagate(DET->RP_420_s);
306 elementFP420 = (TRootRomanPotHits*) branchFP420->NewEntry();
307 //elementFP420 = (TRootForwardTaggerHits*) branchFP420->NewEntry();
308
309 // detector measurements
310 elementFP420->X = (1E-6)*p1.getX(); // [m]
311 elementFP420->Y = (1E-6)*p1.getY(); // [m]
312 elementFP420->Tx = (1E-6)*p1.getTX(); // [rad]
313 elementFP420->Ty = (1E-6)*p1.getTY(); // [rad]
314 elementFP420->S = p1.getS(); // [m]
315 elementFP420->T = time_of_flight(particle, DET->RP_420_s, DET->CEN_max_calo_fwd, DET->RP420_T_resolution);
316 elementFP420->side = sign(particle->Eta);
317
318 // reconstructed data
319 elementFP420->E = p1.getE(); // not yet implemented
320 elementFP420->q2 = UNDEFINED; // not yet implemented
321
322 // generator level data
323 elementFP420->pid = particle->PID;
324 elementFP420->Set(genMomentum);
325 } // if FP420
326
327 } // if stopped
328 } // if forward proton
329
330}
Note: See TracBrowser for help on using the repository browser.