SND@LHC Software
Loading...
Searching...
No Matches
GenieGenerator.cxx
Go to the documentation of this file.
1#include <math.h>
2#include "TSystem.h"
3#include "TROOT.h"
4#include "TMath.h"
5#include "TFile.h"
6#include "TRandom.h"
7#include "TCut.h"
8#include "TEventList.h"
9#include "FairPrimaryGenerator.h"
10#include "GenieGenerator.h"
11#include "TGeoVolume.h"
12#include "TGeoNode.h"
13#include "TGeoManager.h"
14#include "TGeoEltu.h"
15#include "TGeoCompositeShape.h"
16#include "TParticle.h"
17#include "TClonesArray.h"
18
19using std::cout;
20using std::endl;
21// read events from ntuples produced with GENIE
22// http://genie.hepforge.org/manuals/GENIE_PhysicsAndUserManual_20130615.pdf
23// Genie momentum GeV
24// Vertex in SI units, assume this means m
25// important to read back number of events to give to FairRoot
26
27// ----- Default constructor -------------------------------------------
29 fGenOption = 0; //default value, standard Genie Generator used in SHiP
30}
31// -------------------------------------------------------------------------
32// ----- Default constructor -------------------------------------------
33Bool_t GenieGenerator::Init(const char* fileName) {
34 return Init(fileName, 0);
35}
36// ----- Default constructor -------------------------------------------
37Bool_t GenieGenerator::Init(const char* fileName, const int firstEvent) {
38 fNuOnly = false;
39 fcrossingangle = 0.; //default value, no xsec angle
40 if (0 == strncmp("/eos",fileName,4) ) {
41 TString tmp = gSystem->Getenv("EOSSHIP");
42 tmp+=fileName;
43 fInputFile = TFile::Open(tmp);
44 LOGF(info, "Opening input file on eos %s", tmp.Data());
45 }else{
46 fInputFile = new TFile(fileName);
47 LOGF(info, "Opening input file %s", fileName);
48 }
49 if (fInputFile->IsZombie() or !fInputFile) {
50 LOG(FATAL) << "Error opening input file";
51 return kFALSE; }
52 fTree = (TTree *)fInputFile->Get("gst");
53 fNevents = fTree->GetEntries();
54 fn = firstEvent;
55 fTree->SetBranchAddress("Ev",&Ev); // incoming neutrino energy
56 fTree->SetBranchAddress("pxv",&pxv);
57 fTree->SetBranchAddress("pyv",&pyv);
58 fTree->SetBranchAddress("pzv",&pzv);
59 fTree->SetBranchAddress("neu",&neu); // incoming neutrino PDG code
60 fTree->SetBranchAddress("cc",&cc); // Is it a CC event?
61 fTree->SetBranchAddress("nuel",&nuel); // Is it a NUEEL event?
62 fTree->SetBranchAddress("vtxx",&vtxx); // vertex in SI units
63 fTree->SetBranchAddress("vtxy",&vtxy);
64 fTree->SetBranchAddress("vtxz",&vtxz);
65 fTree->SetBranchAddress("vtxt",&vtxt);
66 fTree->SetBranchAddress("El",&El); // outgoing lepton momentum
67 fTree->SetBranchAddress("pxl",&pxl);
68 fTree->SetBranchAddress("pyl",&pyl);
69 fTree->SetBranchAddress("pzl",&pzl);
70 fTree->SetBranchAddress("Ef",&Ef); // outgoing hadronic momenta
71 fTree->SetBranchAddress("pxf",&pxf);
72 fTree->SetBranchAddress("pyf",&pyf);
73 fTree->SetBranchAddress("pzf",&pzf);
74 fTree->SetBranchAddress("nf",&nf); // nr of outgoing hadrons
75 fTree->SetBranchAddress("pdgf",&pdgf); // pdg code of hadron
76
77 if (fGenOption < 0 || fGenOption > 3){
78 LOG(FATAL) <<"Invalid GenieGen Option: "<<fGenOption<<" Please check the option provided with --Genie "<<endl;
79 return kFALSE;
80 }
81 if (fGenOption == 1){
82 fFLUKANuTree = (TTree *)fInputFile->Get("fluka_neutrinos_selected");
83 //check if TTree is actually present
84 if (!fFLUKANuTree){
85 LOG(FATAL) <<"No TTree of interacting neutrinos present. Did you run extract_interacting_neutrinos.py from the macro folder?";
86 return kFALSE;
87 }
88 fNevents = fFLUKANuTree->GetEntries();
89 //setting branches for input neutrino tree.
90 //Nota Bene: I get only the angles and positions,
91 //for the energy I keep the GENIE one, otherwise I lose energy conservation
92 fFLUKANuTree->SetBranchAddress("x_cos",&FLUKA_x_cos);
93 fFLUKANuTree->SetBranchAddress("y_cos",&FLUKA_y_cos);
94 fFLUKANuTree->SetBranchAddress("x",&FLUKA_x);
95 fFLUKANuTree->SetBranchAddress("y",&FLUKA_y);
96 fDeltaE_GenieFLUKA_nu = 10.; //Default value for energy range of FLUKA->Genie matching, 10 GeV.
97 }
98 if (fGenOption == 2){ //pythia tree
99 ancstr = new TClonesArray("TParticle");
100 fFLUKANuTree = (TTree*)fInputFile->Get("NuTauTree");
101 fNevents = fFLUKANuTree->GetEntries();
102 fFLUKANuTree->SetBranchAddress("Ancstr",&ancstr);
103 fDeltaE_GenieFLUKA_nu = 10.; //Default value for energy range of FLUKA->Genie matching, 10 GeV.
104 }
105 fFirst=kTRUE;
106 return kTRUE;
107}
108// -------------------------------------------------------------------------
109Double_t GenieGenerator::MeanMaterialBudget(const Double_t *start, const Double_t *end, Double_t *mparam)
110{
111 //
112 // Calculate mean material budget and material properties between
113 // the points "start" and "end".
114 //
115 // "mparam" - parameters used for the energy and multiple scattering
116 // corrections:
117 //
118 // mparam[0] - mean density: sum(x_i*rho_i)/sum(x_i) [g/cm3]
119 // mparam[1] - equivalent rad length fraction: sum(x_i/X0_i) [adimensional]
120 // mparam[2] - mean A: sum(x_i*A_i)/sum(x_i) [adimensional]
121 // mparam[3] - mean Z: sum(x_i*Z_i)/sum(x_i) [adimensional]
122 // mparam[4] - length: sum(x_i) [cm]
123 // mparam[5] - Z/A mean: sum(x_i*Z_i/A_i)/sum(x_i) [adimensional]
124 // mparam[6] - number of boundary crosses
125 // mparam[7] - maximum density encountered (g/cm^3)
126 // mparam[8] - equivalent interaction length fraction: sum(x_i/I0_i) [adimensional]
127 // mparam[9] - maximum cross section encountered (mbarn)
128 //
129 // Origin: Marian Ivanov, Marian.Ivanov@cern.ch
130 //
131 // Corrections and improvements by
132 // Andrea Dainese, Andrea.Dainese@lnl.infn.it,
133 // Andrei Gheata, Andrei.Gheata@cern.ch
134 // Thomas Ruf, Thomas.Ruf@cern.ch
135 //
136
137 mparam[0]=0; mparam[1]=1; mparam[2]=0; mparam[3]=0; mparam[4]=0;
138 mparam[5]=0; mparam[6]=0; mparam[7]=0; mparam[8]=0; mparam[9]=0;
139 //
140 Double_t bparam[7]; // total parameters
141 Double_t lparam[7]; // local parameters
142 Double_t mbarn = 1E-3*1E-24*TMath::Na(); // cm^2 * Avogadro
143
144 for (Int_t i=0;i<7;i++) bparam[i]=0;
145
146 if (!gGeoManager) {
147 //AliFatalClass("No TGeo\n");
148 return 0.;
149 }
150 //
151 Double_t length;
152 Double_t dir[3];
153 length = TMath::Sqrt((end[0]-start[0])*(end[0]-start[0])+
154 (end[1]-start[1])*(end[1]-start[1])+
155 (end[2]-start[2])*(end[2]-start[2]));
156 mparam[4]=length;
157 if (length<TGeoShape::Tolerance()) return 0.0;
158 Double_t invlen = 1./length;
159 dir[0] = (end[0]-start[0])*invlen;
160 dir[1] = (end[1]-start[1])*invlen;
161 dir[2] = (end[2]-start[2])*invlen;
162
163 // Initialize start point and direction
164 TGeoNode *currentnode = 0;
165 TGeoNode *startnode = gGeoManager->InitTrack(start, dir);
166 if (!startnode) {
167 //AliErrorClass(Form("start point out of geometry: x %f, y %f, z %f",
168 // start[0],start[1],start[2]));
169 return 0.0;
170 }
171 TGeoMaterial *material = startnode->GetVolume()->GetMedium()->GetMaterial();
172 lparam[0] = material->GetDensity();
173 if (lparam[0] > mparam[7]) mparam[7]=lparam[0];
174 lparam[1] = material->GetRadLen();
175 lparam[2] = material->GetA();
176 lparam[3] = material->GetZ();
177 lparam[4] = length;
178 lparam[5] = lparam[3]/lparam[2];
179 lparam[6] = material->GetIntLen();
180 Double_t n = lparam[0]/lparam[2];
181 Double_t sigma = 1./(n*lparam[6])/mbarn;
182 if (sigma > mparam[9]) mparam[9]=sigma;
183 if (material->IsMixture()) {
184 TGeoMixture * mixture = (TGeoMixture*)material;
185 lparam[5] =0;
186 Double_t sum =0;
187 for (Int_t iel=0;iel<mixture->GetNelements();iel++){
188 sum += mixture->GetWmixt()[iel];
189 lparam[5]+= mixture->GetZmixt()[iel]*mixture->GetWmixt()[iel]/mixture->GetAmixt()[iel];
190 }
191 lparam[5]/=sum;
192 }
193
194 // Locate next boundary within length without computing safety.
195 // Propagate either with length (if no boundary found) or just cross boundary
196 gGeoManager->FindNextBoundaryAndStep(length, kFALSE);
197 Double_t step = 0.0; // Step made
198 Double_t snext = gGeoManager->GetStep();
199 // If no boundary within proposed length, return current density
200 if (!gGeoManager->IsOnBoundary()) {
201 mparam[0] = lparam[0];
202 mparam[1] = lparam[4]/lparam[1];
203 mparam[2] = lparam[2];
204 mparam[3] = lparam[3];
205 mparam[4] = lparam[4];
206 mparam[8] = lparam[4]/lparam[6];
207 return lparam[0];
208 }
209 // Try to cross the boundary and see what is next
210 Int_t nzero = 0;
211 while (length>TGeoShape::Tolerance()) {
212 currentnode = gGeoManager->GetCurrentNode();
213 if (snext<2.*TGeoShape::Tolerance()) nzero++;
214 else nzero = 0;
215 if (nzero>3) {
216 // This means navigation has problems on one boundary
217 // Try to cross by making a small step
218 //AliErrorClass("Cannot cross boundary\n");
219 mparam[0] = bparam[0]/step;
220 mparam[1] = bparam[1];
221 mparam[2] = bparam[2]/step;
222 mparam[3] = bparam[3]/step;
223 mparam[5] = bparam[5]/step;
224 mparam[8] = bparam[6];
225 mparam[4] = step;
226 mparam[0] = 0.; // if crash of navigation take mean density 0
227 mparam[1] = 1000000; // and infinite rad length
228 return bparam[0]/step;
229 }
230 mparam[6]+=1.;
231 step += snext;
232 bparam[1] += snext/lparam[1];
233 bparam[2] += snext*lparam[2];
234 bparam[3] += snext*lparam[3];
235 bparam[5] += snext*lparam[5];
236 bparam[6] += snext/lparam[6];
237 bparam[0] += snext*lparam[0];
238
239 if (snext>=length) break;
240 if (!currentnode) break;
241 length -= snext;
242 material = currentnode->GetVolume()->GetMedium()->GetMaterial();
243 lparam[0] = material->GetDensity();
244 if (lparam[0] > mparam[7]) mparam[7]=lparam[0];
245 lparam[1] = material->GetRadLen();
246 lparam[2] = material->GetA();
247 lparam[3] = material->GetZ();
248 lparam[5] = lparam[3]/lparam[2];
249 lparam[6] = material->GetIntLen();
250 n = lparam[0]/lparam[2];
251 sigma = 1./(n*lparam[6])/mbarn;
252 if (sigma > mparam[9]) mparam[9]=sigma;
253 if (material->IsMixture()) {
254 TGeoMixture * mixture = (TGeoMixture*)material;
255 lparam[5]=0;
256 Double_t sum =0;
257 for (Int_t iel=0;iel<mixture->GetNelements();iel++){
258 sum+= mixture->GetWmixt()[iel];
259 lparam[5]+= mixture->GetZmixt()[iel]*mixture->GetWmixt()[iel]/mixture->GetAmixt()[iel];
260 }
261 lparam[5]/=sum;
262 }
263 gGeoManager->FindNextBoundaryAndStep(length, kFALSE);
264 snext = gGeoManager->GetStep();
265 }
266 mparam[0] = bparam[0]/step;
267 mparam[1] = bparam[1];
268 mparam[2] = bparam[2]/step;
269 mparam[3] = bparam[3]/step;
270 mparam[5] = bparam[5]/step;
271 mparam[8] = bparam[6];
272 return bparam[0]/step;
273}
274
275std::vector<double> GenieGenerator::Rotate(Double_t x, Double_t y, Double_t z, Double_t px, Double_t py, Double_t pz)
276{
277 //rotate vector px,py,pz to point at x,y,z at origin.
278 Double_t theta=atan(sqrt(x*x+y*y)/z);
279 Double_t c=cos(theta);
280 Double_t s=sin(theta);
281 //rotate around y-axis
282 Double_t px1=c*px+s*pz;
283 Double_t pzr=-s*px+c*pz;
284 Double_t phi=atan2(y,x);
285 c=cos(phi);
286 s=sin(phi);
287 //rotate around z-axis
288 Double_t pxr=c*px1-s*py;
289 Double_t pyr=s*px1+c*py;
290 std::vector<double> pout;
291 pout.push_back(pxr);
292 pout.push_back(pyr);
293 pout.push_back(pzr);
294 //cout << "Info GenieGenerator: rotated" << pout[0] << " " << pout[1] << " " << pout[2] << " " << x << " " << y << " " << z <<endl;
295 return pout;
296}
297
298Int_t GenieGenerator::ExtractEvent_Ekin(Double_t Ekin, Double_t DeltaE){
299
300 //TCut nucut(Form("Ekin >= (%f - %f) && Ekin < (%f + %f)",Ekin,DeltaE/2., Ekin, DeltaE/2.));
301 //fFLUKANuTree->Draw(">>nulist", nucut );
302 TCut nucut(Form("pzv >= (%f - %f) && pzv < (%f + %f)",Ekin,DeltaE/2., Ekin, DeltaE/2.));
303 fTree->Draw(">>nulist",nucut);
304
305 TEventList *nulist = (TEventList*)gDirectory->GetList()->FindObject("nulist");
306 Int_t nselectedevents = nulist->GetN();
307 //integer function returns randomly integer between 0 e nselectedevents;
308 Int_t myevent = nulist->GetEntry(gRandom->Integer(nselectedevents+1));
309
310 nulist->Reset(); //resetting list for next call of the function.
311 return myevent;
312}
313
314
315// ----- Destructor ----------------------------------------------------
317{
318 fInputFile->Close();
319 fInputFile->Delete();
320 delete fInputFile;
321}
322// -------------------------------------------------------------------------
323void GenieGenerator::AddBox(TVector3 dVec, TVector3 box)
324{
325 dVecs.push_back(dVec);
326 boxs.push_back(box);
327 cout << "Debug GenieGenerator: "<< dVec.X() << " "<<box.Z() << endl;
328}
329// ----- Passing the event ---------------------------------------------
330Bool_t GenieGenerator::OldReadEvent(FairPrimaryGenerator* cpg)
331{
332 if (fFirst){
333 TGeoVolume *top=gGeoManager->GetTopVolume();
334 TGeoNode *linner = top->FindNode("lidT1I_1");
335 TGeoNode *scint = top->FindNode("lidT1lisci_1");
336 TGeoNode *louter = top->FindNode("lidT1O_1");
337 TGeoNode *lidT6I = top->FindNode("lidT6I_1");
338 TGeoNode *t2I = top->FindNode("T2I_1");
339 TGeoNode *t1I = top->FindNode("T1I_1");
340 TGeoEltu *temp = dynamic_cast<TGeoEltu*>(linner->GetVolume()->GetShape());
341 fEntrDz_inner = temp->GetDZ();
342 temp = dynamic_cast<TGeoEltu*>(louter->GetVolume()->GetShape());
343 fEntrDz_outer = temp->GetDZ();
344 fEntrA = temp->GetA();
345 fEntrB = temp->GetB();
346 fEntrZ_inner = linner->GetMatrix()->GetTranslation()[2];
347 fEntrZ_outer = louter->GetMatrix()->GetTranslation()[2];
348 Lvessel = lidT6I->GetMatrix()->GetTranslation()[2] - fEntrZ_inner - fEntrDz_inner;
349 TGeoCompositeShape *tempC = dynamic_cast<TGeoCompositeShape*>(t2I->GetVolume()->GetShape());
350 Xvessel = tempC->GetDX() - 2*fEntrDz_inner ;
351 Yvessel = tempC->GetDY() - 2*fEntrDz_inner ;
352 tempC = dynamic_cast<TGeoCompositeShape*>(t1I->GetVolume()->GetShape());
353 fL1z = tempC->GetDZ()*2;
354 temp = dynamic_cast<TGeoEltu*>(scint->GetVolume()->GetShape());
355 fScintDz = temp->GetDZ()*2;
356 cout << "Info GenieGenerator: geo inner " << fEntrDz_inner << " "<< fEntrZ_inner << endl;
357 cout << "Info GenieGenerator: geo outer " << fEntrDz_outer << " "<< fEntrZ_outer << endl;
358 cout << "Info GenieGenerator: A and B " << fEntrA << " "<< fEntrB << endl;
359 cout << "Info GenieGenerator: vessel length heig<ht width " << Lvessel << " "<<Yvessel<< " "<< Xvessel << endl;
360 cout << "Info GenieGenerator: scint thickness " << fScintDz << endl;
361 cout << "Info GenieGenerator: rextra " << fScintDz/2.+2*fEntrDz_inner << " "<< 2*fEntrDz_outer << " "<<2*fEntrDz_inner << endl;
362 for(int j=0; j<boxs.size(); j++) {
363 cout << "Info GenieGenerator: nuMu X" << j << " - "<< -boxs[j].X()+dVecs[j].X() << " "<< boxs[j].X()+dVecs[j].X() << endl;
364 cout << "Info GenieGenerator: nuMu Y" << j << " - "<< -boxs[j].Y()+dVecs[j].Y() << " "<< boxs[j].Y()+dVecs[j].Y() << endl;
365 cout << "Info GenieGenerator: nuMu Z" << j << " - "<< -boxs[j].Z()+dVecs[j].Z() << " "<< boxs[j].Z()+dVecs[j].Z() << endl;
366 }
367 fFirst = kFALSE;
368 }
369 if (fn==fNevents) {LOG(WARNING) << "End of input file. Rewind.";}
370 fTree->GetEntry(fn%fNevents);
371 fn++;
372 if (fn%1000==0) {
373 cout << "Info GenieGenerator: neutrino event-nr "<< fn << endl;
374 }
375// generate a random point on the vessel, take veto z, and calculate outer lid position
376 //Double_t Yvessel=500.;
377 //Double_t Lvessel=5.*Yvessel+3500.;
378 //Double_t ztarget=zentrance-6350.;
379 //Double_t ea=250.; //inside inner wall vessel
380 Double_t eb=Yvessel; //inside inner wall vessel
381 Double_t x;
382 Double_t y;
383 Double_t z;
384 Double_t where=gRandom->Uniform(0.,1.);
385 if (where<0.03) {
386 // point on entrance lid
387 Double_t ellip=2.;
388 while (ellip>1.){
389 x = gRandom->Uniform(-fEntrA,fEntrA);
390 y = gRandom->Uniform(-fEntrB,fEntrB);
391 ellip=(x*x)/(fEntrA*fEntrA)+(y*y)/(fEntrB*fEntrB);
392 }
393 if (gRandom->Uniform(0.,1.)<0.5){
394 z=fEntrZ_inner + gRandom->Uniform(-fEntrDz_inner,fEntrDz_inner);
395 }else{
396 z=fEntrZ_outer + gRandom->Uniform(-fEntrDz_outer,fEntrDz_outer);
397 }
398 }else if (where<0.64){
399 //point on tube, place over 1 cm radius at 2 radii, separated by 10. cm
400 // first vessel has smaller size
401 Double_t ea = Xvessel;
402 Double_t zrand = gRandom->Uniform(0,Lvessel);
403 if (zrand < fL1z){
404 ea = fEntrA;
405 }
406 z = fEntrZ_outer-fEntrDz_outer + zrand;
407 Double_t theta = gRandom->Uniform(0.,TMath::Pi());
408 Double_t rextra;
409 if ( gRandom->Uniform(0.,1.)>0.5) {
410 // outer vessel
411 rextra = fScintDz/2.+2*fEntrDz_inner + gRandom->Uniform(0,2*fEntrDz_outer);
412 }else{
413 // inner vessel
414 rextra = gRandom->Uniform(0.,2*fEntrDz_inner);
415 }
416 x = (ea+rextra)*cos(theta);
417 y = sqrt(1.-(x*x)/((ea+rextra)*(ea+rextra)))*(eb+rextra);
418 if (gRandom->Uniform(-1.,1.)>0.) y=-y;
419 }else{
420 //point in nu-tau detector area
421 int j = int(gRandom->Uniform(0.,boxs.size()+0.5));
422 x=gRandom->Uniform(-boxs[j].X()+dVecs[j].X(),boxs[j].X()+dVecs[j].X());
423 y=gRandom->Uniform(-boxs[j].Y()+dVecs[j].Y(),boxs[j].Y()+dVecs[j].Y());
424 z=gRandom->Uniform(-boxs[j].Z()+dVecs[j].Z(),boxs[j].Z()+dVecs[j].Z());
425 }
426// first, incoming neutrino
427 //rotate the particle
428 Double_t zrelative=z-ztarget;
429 //cout << "Info GenieGenerator: x,y,z " << x <<" " << y << " " << zrelative << endl;
430 //cout << "Info GenieGenerator: neutrino " << neu << "p-in "<< pxv << pyv << pzv << " nf "<< nf << endl;
431 std::vector<double> pout = Rotate(x,y,zrelative,pxv,pyv,pzv);
432 //cpg->AddTrack(neu,pxv,pyv,pzv,vtxx,vtxy,vtxz,-1,false);
433 //cout << "Info GenieGenerator: neutrino " << neu << "p-rot "<< pout[0] << " fn "<< fn << endl;
434 cpg->AddTrack(neu,pout[0],pout[1],pout[2],x,y,z,-1,false);
435
436
437// second, outgoing lepton
438 pout = Rotate(x,y,zrelative,pxl,pyl,pzl);
439 Int_t oLPdgCode = neu;
440 if (cc){oLPdgCode = copysign(TMath::Abs(neu)-1,neu);}
441 if (nuel){oLPdgCode = 11;}
442 cpg->AddTrack(oLPdgCode,pout[0],pout[1],pout[2],x,y,z,0,true);
443// last, all others
444 for(int i=0; i<nf; i++)
445 {
446 pout = Rotate(x,y,zrelative,pxf[i],pyf[i],pzf[i]);
447 cpg->AddTrack(pdgf[i],pout[0],pout[1],pout[2],x,y,z,0,true);
448 // cout << "f " << pdgf[i] << " pz "<< pzf[i] << endl;
449 }
450
451 return kTRUE;
452}
453
454// ----- Passing the event ---------------------------------------------
455Bool_t GenieGenerator::ReadEvent(FairPrimaryGenerator* cpg)
456{
457
458 if (fGenOption == 3){
459 // Use GENIE geometry driver.
460
461 // Get event from GENIE TTree. If we reach the end of the file, return false.
462 if (fTree->GetEntry(fn) == 0 ) return kFALSE;
463
464 if (fn%100==0) {
465 cout << "Info GenieGenerator: neutrino event-nr "<< fn << endl;
466 }
467
468 fn++;
469
470 // Add the neutrino to the MCTrack stack:
471 cpg->AddTrack(neu, // Neutrino PDG
472 pxv, pyv, pzv, // Neutrino momentum
473 vtxx*100, vtxy*100, vtxz*100, // Event vertex [in cm!]
474 -1, // Parent
475 false); // Don't track this particle
476
477 // Add final state lepton to the MCTrack stack:
478 int outgoing_lepton_pdg = neu;
479 if (cc) outgoing_lepton_pdg = copysign(TMath::Abs(neu)-1,neu);
480 if (nuel) outgoing_lepton_pdg = 11;
481
482 bool track_outgoing_lepton = (cc || nuel) ? true : false;
483 cpg->AddTrack(outgoing_lepton_pdg,
484 pxl, pyl, pzl,
485 vtxx*100, vtxy*100, vtxz*100,
486 0,
487 track_outgoing_lepton);
488
489 // Add final state hadrons to the MCTrack stack
490 for (int i_hadron = 0; i_hadron < nf; i_hadron++){
491 cpg->AddTrack(pdgf[i_hadron],
492 pxf[i_hadron], pyf[i_hadron], pzf[i_hadron],
493 vtxx*100, vtxy*100, vtxz*100,
494 0,
495 true);
496 }
497
498 return kTRUE;
499
500 } else {
501 //some start/end positions in z (emulsion to Tracker 1)
502 Double_t start[3]={0.,0.,startZ};
503 Double_t end[3]={0.,0.,endZ};
504 char ts[20];
505 //cout << "Enter GenieGenerator " << endl;
506 //pick histogram: 1100=100 momentum bins, 1200=25 momentum bins.
507 Int_t idbase=1200;
508 if (fFirst){
509 Double_t bparam=0.;
510 Double_t mparam[10];
511 bparam=MeanMaterialBudget(start, end, mparam);
512 cout << "Info GenieGenerator: MaterialBudget " << start[2] << " - "<< end[2] << endl;
513 cout << "Info GenieGenerator: MaterialBudget " << bparam << endl;
514 cout << "Info GenieGenerator: MaterialBudget 0 " << mparam[0] << endl;
515 cout << "Info GenieGenerator: MaterialBudget 1 " << mparam[1] << endl;
516 cout << "Info GenieGenerator: MaterialBudget 2 " << mparam[2] << endl;
517 cout << "Info GenieGenerator: MaterialBudget 3 " << mparam[3] << endl;
518 cout << "Info GenieGenerator: MaterialBudget 4 " << mparam[4] << endl;
519 cout << "Info GenieGenerator: MaterialBudget 5 " << mparam[5] << endl;
520 cout << "Info GenieGenerator: MaterialBudget 6 " << mparam[6] << endl;
521 cout << "Info GenieGenerator: MaterialBudget " << mparam[0]*mparam[4] << endl;
522 //read the (log10(p),log10(pt)) hists to be able to draw a pt for every neutrino momentum
523 //loop over neutrino types
524 printf("Reading (log10(p),log10(pt)) Hists from file: %s\n",fInputFile->GetName());
525 for (Int_t idnu=12;idnu<17;idnu+=2){
526 for (Int_t idadd=-1;idadd<2;idadd+=2){
527 Int_t idhnu=idbase+idnu;
528 if (idadd<0) idhnu+=1000;
529 sprintf(ts,"%d",idhnu);
530 //pickup corresponding (log10(p),log10(pt)) histogram
531 if (fInputFile->FindObjectAny(ts)){
532 TH2F* h2tmp = (TH2F*) fInputFile->Get(ts);
533 printf("HISTID=%d, Title:%s\n",idhnu,h2tmp->GetTitle());
534 sprintf(ts,"px_%d",idhnu);
535 //make its x-projection, to later be able to convert log10(p) to its bin-number
536 pxhist[idhnu]=h2tmp->ProjectionX(ts,1,-1);
537 Int_t nbinx=h2tmp->GetNbinsX();
538 //printf("idhnu=%d ts=%s nbinx=%d\n",idhnu,ts,nbinx);
539 //project all slices on the y-axis
540 for (Int_t k=1;k<nbinx+1;k+=1){
541 sprintf(ts,"h%d%d",idhnu,k);
542 //printf("idnu %d idhnu %d bin%d ts=%s\n",idnu,idhnu,k,ts);
543 pyslice[idhnu][k]=h2tmp->ProjectionY(ts,k,k);
544 }
545 }
546 }
547 }
548 fFirst = kFALSE;
549 }
550 if (fn==fNevents) {LOG(WARNING) << "End of input file. Rewind.";}
551 if (fGenOption>0) fFLUKANuTree->GetEntry(fn%fNevents); //1 or 2
552 else if (fGenOption == 0) fTree->GetEntry(fn%fNevents);
553 fn++;
554 if (fn%100==0) {
555 cout << "Info GenieGenerator: neutrino event-nr "<< fn << endl;
556 }
557
558// Incoming neutrino, get a random px,py
559 //cout << "Info GenieGenerator: neutrino " << neu << "p-in "<< pzv << " nf "<< nf << endl;
560 //cout << "Info GenieGenerator: ztarget " << ztarget << endl;
561 Double_t bparam=0.;
562 Double_t mparam[10];
563 Double_t pout[3];
564 pout[2]=-1.;
565 Double_t txnu=0;
566 Double_t tynu=0;
567 //Does this neutrino fly through material? Otherwise draw another pt..
568 //cout << "Info GenieGenerator Start bparam while loop" << endl;
569 while (pout[2]<0.) {
570 //***OLD**** Keep for comparison maybe??
571 //generate pt of ~0.3 GeV
572 //pout[0] = gRandom->Exp(0.2);
573 //pout[1] = gRandom->Exp(0.2);
574 //pout[2] = pzv*pzv-pout[0]*pout[0]-pout[1]*pout[1];
575 if(fGenOption == 1){
576 int nuevent = ExtractEvent_Ekin(pzv, 10.);
577 fTree->GetEntry(nuevent);
578 //getting tri-momentum
579 pout[0] = FLUKA_x_cos * pzv;
580 pout[1] = FLUKA_y_cos * pzv;
581 double pt = sqrt(pout[0]*pout[0] + pout[1]*pout[1]);
582 pout[2] = pzv*pzv-pt*pt;
583 }
584 else if(fGenOption == 2){
585 //getting neutrino information
586 TParticle *nuparticle = (TParticle*) ancstr->At(0);
587 TVector3 * nup = new TVector3(nuparticle->Px(), nuparticle->Py(), nuparticle->Pz());
588 //rotating of xsec
589 nup->RotateX(fcrossingangle);
590 //getting associated interaction
591 int nuevent = ExtractEvent_Ekin(nup->Mag(), 10.);
592 fTree->GetEntry(nuevent);
593 //getting GENIE tri-momentum in this reference
594 pout[0] = nup->Px()/nup->Pz() * pzv;
595 pout[1] = nup->Py()/nup->Pz() * pzv;
596 double pt = sqrt(pout[0]*pout[0] + pout[1]*pout[1]);
597 pout[2] = pzv*pzv-pt*pt;
598 }
599
600 else if( fGenOption == 0 ){
601 //**NEW** get pt of this neutrino from 2D hists.
602 Int_t idhnu=TMath::Abs(neu)+idbase;
603 if (neu<0) idhnu+=1000;
604 Int_t nbinmx=pxhist[idhnu]->GetNbinsX();
605 Double_t pl10=log10(pzv);
606 Int_t nbx=pxhist[idhnu]->FindBin(pl10);
607 //printf("idhnu %d, p %f log10(p) %f bin,binmx %d %d \n",idhnu,pzv,pl10,nbx,nbinmx);
608 if (nbx<1) nbx=1;
609 if (nbx>nbinmx) nbx=nbinmx;
610 Double_t ptlog10=pyslice[idhnu][nbx]->GetRandom();
611 //hist was filled with: log10(pt+0.01)
612 Double_t pt=pow(10.,ptlog10)-0.01;
613 //rotate pt in phi:
614 Double_t phi=gRandom->Uniform(0.,2*TMath::Pi());
615 pout[0] = cos(phi)*pt;
616 pout[1] = sin(phi)*pt;
617 pout[2] = pzv*pzv-pt*pt;
618 //printf("p= %f pt=%f px,py,pz**2=%f,%f,%f\n",pzv,pt,pout[0],pout[1],pout[2]);
619 }
620 if (pout[2]>=0.) {
621 pout[2]=TMath::Sqrt(pout[2]);
622 //random chance of reflecting momentum, we do not modify momentum components in SND@LHC physics case
623 if (fGenOption == 0){
624 if (gRandom->Uniform(-1.,1.)<0.) pout[0]=-pout[0];
625 if (gRandom->Uniform(-1.,1.)<0.) pout[1]=-pout[1];
626 }
627 //cout << "Info GenieGenerator: neutrino pxyz " << pout[0] << ", " << pout[1] << ", " << pout[2] << endl;
628 // xyz at start and end
629 start[0]=(pout[0]/pout[2])*(start[2]-ztarget);
630 start[1]=(pout[1]/pout[2])*(start[2]-ztarget);
631 //adding offset from neutrino production point (only for SND@LHC physics case)
632 if (fGenOption == 1){
633 start[0] += FLUKA_x;
634 start[1] += FLUKA_y;
635 }
636 //cout << "Info GenieGenerator: neutrino xyz-start " << start[0] << "-" << start[1] << "-" << start[2] << endl;
637 txnu=pout[0]/pout[2];
638 tynu=pout[1]/pout[2];
639 end[0]=txnu*(end[2]-ztarget);
640 end[1]=tynu*(end[2]-ztarget);
641 //adding offset from neutrino production point (only for SND@LHC phys
642 if (fGenOption == 1){
643 end[0] += FLUKA_x;
644 end[1] += FLUKA_y;
645 }
646 //cout << "Info GenieGenerator: neutrino xyz-end " << end[0] << "-" << end[1] << "-" << end[2] << endl;
647 //get material density between these two points
648 bparam=MeanMaterialBudget(start, end, mparam);
649 //printf("param %e %e %e \n",bparam,mparam[6],mparam[7]);
650 }
651 }
652 //loop over trajectory between start and end to pick an interaction point
653 Double_t prob2int = -1.;
654 Double_t x;
655 Double_t y;
656 Double_t z;
657 Int_t count=0;
658 while (prob2int<gRandom->Uniform(0.,1.)) {
659 //place x,y,z uniform along path
660 z=gRandom->Uniform(start[2],end[2]);
661 x=txnu*(z-ztarget);
662 y=tynu*(z-ztarget);
663 //adding offset from neutrino production point (only for SND@LHC physics case)
664 if (fGenOption == 1){
665 x += FLUKA_x;
666 y += FLUKA_y;
667 }
668 if (mparam[6]<0.5){
669 //mparam is number of boundaries along path. mparam[6]=0.: uniform material budget along path, use present x,y,z
670 prob2int=2.;
671 }else{
672 //get local material at this point, to calculate probability that interaction is at this point.
673 TGeoNode *node = gGeoManager->FindNode(x,y,z);
674 TGeoMaterial *mat = 0;
675 if (node && !gGeoManager->IsOutside()) {
676 mat = node->GetVolume()->GetMaterial();
677 //cout << "Info GenieGenerator: mat " << count << ", " << mat->GetName() << ", " << mat->GetDensity() << endl;
678 //density relative to Prob largest density along this trajectory, i.e. use rho(Pt)
679 prob2int= mat->GetDensity()/mparam[7];
680 if (prob2int>1.) cout << "***WARNING*** GenieGenerator: prob2int > Maximum density????" << prob2int << " maxrho:" << mparam[7] << " material: " << mat->GetName() << endl;
681 count+=1;
682 }else{
683 prob2int=0.;
684 }
685 }
686 }
687 //cout << "Info GenieGenerator: prob2int " << prob2int << ", " << count << endl;
688
689 Double_t zrelative=z-ztarget;
690 Double_t tof=TMath::Sqrt(x*x+y*y+zrelative*zrelative)/2.99792458e+6;
691 cpg->AddTrack(neu,pout[0],pout[1],pout[2],x,y,z,-1,false,TMath::Sqrt(pout[0]*pout[0]+pout[1]*pout[1]+pout[2]*pout[2]),tof,mparam[0]*mparam[4]);
692 if (not fNuOnly){
693 // second, outgoing lepton
694 std::vector<double> pp = Rotate(x,y,zrelative,pxl,pyl,pzl);
695 Int_t oLPdgCode = neu;
696 if (cc){oLPdgCode = copysign(TMath::Abs(neu)-1,neu);}
697 if (nuel){oLPdgCode = 11;}
698 cpg->AddTrack(oLPdgCode,pp[0],pp[1],pp[2],x,y,z,0,true,El,tof,mparam[0]*mparam[4]);
699// last, all others
700 for(int i=0; i<nf; i++)
701 {
702 pp = Rotate(x,y,zrelative,pxf[i],pyf[i],pzf[i]);
703 cpg->AddTrack(pdgf[i],pp[0],pp[1],pp[2],x,y,z,0,true,Ef[i],tof,mparam[0]*mparam[4]);
704 //cout << "f " << pdgf[i] << " pz "<< pzf[i] << endl;
705 }
706 //cout << "Info GenieGenerator Return from GenieGenerator" << endl;
707 }
708 }
709 return kTRUE;
710}
711// -------------------------------------------------------------------------
713{
714 return fNevents;
715}
716
717
const Double_t mbarn
Double_t pzf[500]
TClonesArray * ancstr
Double_t fDeltaE_GenieFLUKA_nu
Double_t fEntrZ_inner
std::vector< double > Rotate(Double_t x, Double_t y, Double_t z, Double_t px, Double_t py, Double_t pz)
virtual ~GenieGenerator()
Double_t fEntrDz_outer
Double_t fEntrZ_outer
Int_t pdgf[500]
virtual Bool_t Init(const char *, int)
TH1D * pyslice[3000][500]
std::vector< TVector3 > dVecs
TH1D * pxhist[3000]
Double_t pxf[500]
Bool_t ReadEvent(FairPrimaryGenerator *)
Double_t MeanMaterialBudget(const Double_t *start, const Double_t *end, Double_t *mparam)
Bool_t OldReadEvent(FairPrimaryGenerator *)
Double_t FLUKA_y_cos
Double_t pyf[500]
void AddBox(TVector3 dVec, TVector3 box)
std::vector< TVector3 > boxs
Double_t fcrossingangle
Int_t ExtractEvent_Ekin(Double_t Ekin, Double_t DeltaE)
Double_t fEntrDz_inner
Double_t Ef[500]
don't make it persistent, magic ROOT command
Double_t FLUKA_x_cos
ClassImp(ecalContFact) ecalContFact