1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 #include <Riostream.h>
20 #include "TClonesArray.h"
22 #include "TCallf77.h" //For the fortran calls
23 #include "Fdblprc.h" //(DBLPRC) fluka common
24 #include "Fepisor.h" //(EPISOR) fluka common
25 #include "Ffinuc.h" //(FINUC) fluka common
26 #include "Fiounit.h" //(IOUNIT) fluka common
27 #include "Fpaprop.h" //(PAPROP) fluka common
28 #include "Fpart.h" //(PART) fluka common
29 #include "Ftrackr.h" //(TRACKR) fluka common
30 #include "Fpaprop.h" //(PAPROP) fluka common
31 #include "Ffheavy.h" //(FHEAVY) fluka common
33 #include "TVirtualMC.h"
34 #include "TG4GeometryManager.h" //For the geometry management
35 #include "TG4DetConstruction.h" //For the detector construction
37 #include "FGeometryInit.hh"
38 #include "TLorentzVector.h"
39 #include "FlukaVolume.h"
41 // Fluka methods that may be needed.
43 # define flukam flukam_
44 # define fluka_openinp fluka_openinp_
45 # define fluka_closeinp fluka_closeinp_
46 # define mcihad mcihad_
47 # define mpdgha mpdgha_
49 # define flukam FLUKAM
50 # define fluka_openinp FLUKA_OPENINP
51 # define fluka_closeinp FLUKA_CLOSEINP
52 # define mcihad MCIHAD
53 # define mpdgha MPDGHA
59 // Prototypes for FLUKA functions
61 void type_of_call flukam(const int&);
62 void type_of_call fluka_openinp(const int&, DEFCHARA);
63 void type_of_call fluka_closeinp(const int&);
64 int type_of_call mcihad(const int&);
65 int type_of_call mpdgha(const int&);
69 // Class implementation for ROOT
74 //----------------------------------------------------------------------------
75 // TFluka constructors and destructors.
76 //____________________________________________________________________________
82 fCurrentFlukaRegion(-1)
85 // Default constructor
89 TFluka::TFluka(const char *title, Int_t verbosity)
90 :TVirtualMC("TFluka",title),
91 fVerbosityLevel(verbosity),
96 fCurrentFlukaRegion(-1)
98 if (fVerbosityLevel >=3)
99 cout << "==> TFluka::TFluka(" << title << ") constructor called." << endl;
102 // create geometry manager
103 if (fVerbosityLevel >=2)
104 cout << "\t* Creating G4 Geometry manager..." << endl;
105 fGeometryManager = new TG4GeometryManager();
106 if (fVerbosityLevel >=2)
107 cout << "\t* Creating G4 Detector..." << endl;
108 fDetector = new TG4DetConstruction();
109 FGeometryInit* geominit = FGeometryInit::GetInstance();
111 geominit->setDetConstruction(fDetector);
113 cerr << "ERROR: Could not create FGeometryInit!" << endl;
114 cerr << " Exiting!!!" << endl;
118 if (fVerbosityLevel >=3)
119 cout << "<== TFluka::TFluka(" << title << ") constructor called." << endl;
121 fVolumeMediaMap = new TClonesArray("FlukaVolume",1000);
127 if (fVerbosityLevel >=3)
128 cout << "==> TFluka::~TFluka() destructor called." << endl;
130 delete fGeometryManager;
131 fVolumeMediaMap->Delete();
132 delete fVolumeMediaMap;
135 if (fVerbosityLevel >=3)
136 cout << "<== TFluka::~TFluka() destructor called." << endl;
140 //_____________________________________________________________________________
141 // TFluka control methods
142 //____________________________________________________________________________
143 void TFluka::Init() {
145 if (fVerbosityLevel >=3)
146 cout << "==> TFluka::Init() called." << endl;
148 cout << "\t* InitPhysics() - Prepare input file to be called" << endl;
149 InitPhysics(); // prepare input file
150 cout << "\t* InitPhysics() - Prepare input file called" << endl;
152 if (fVerbosityLevel >=2)
153 cout << "\t* Changing lfdrtr = (" << (GLOBAL.lfdrtr?'T':'F')
154 << ") in fluka..." << endl;
155 GLOBAL.lfdrtr = true;
157 if (fVerbosityLevel >=2)
158 cout << "\t* Opening file " << sInputFileName << endl;
159 const char* fname = sInputFileName;
160 fluka_openinp(lunin, PASSCHARA(fname));
162 if (fVerbosityLevel >=2)
163 cout << "\t* Calling flukam..." << endl;
166 if (fVerbosityLevel >=2)
167 cout << "\t* Closing file " << sInputFileName << endl;
168 fluka_closeinp(lunin);
172 if (fVerbosityLevel >=3)
173 cout << "<== TFluka::Init() called." << endl;
177 void TFluka::FinishGeometry() {
179 // Build-up table with region to medium correspondance
183 if (fVerbosityLevel >=3)
184 cout << "==> TFluka::FinishGeometry() called." << endl;
186 // fGeometryManager->Ggclos();
188 FGeometryInit* flugg = FGeometryInit::GetInstance();
190 fMediaByRegion = new Int_t[fNVolumes+2];
191 for (Int_t i = 0; i < fNVolumes; i++)
193 FlukaVolume* vol = dynamic_cast<FlukaVolume*>((*fVolumeMediaMap)[i]);
194 TString volName = vol->GetName();
195 Int_t media = vol->GetMedium();
196 if (fVerbosityLevel >= 3)
197 printf("Finish Geometry: volName, media %d %s %d \n", i, volName.Data(), media);
198 strcpy(tmp, volName.Data());
200 flugg->SetMediumFromName(tmp, media, i+1);
201 fMediaByRegion[i] = media;
204 flugg->BuildMediaMap();
206 if (fVerbosityLevel >=3)
207 cout << "<== TFluka::FinishGeometry() called." << endl;
210 void TFluka::BuildPhysics() {
211 if (fVerbosityLevel >=3)
212 cout << "==> TFluka::BuildPhysics() called." << endl;
215 if (fVerbosityLevel >=3)
216 cout << "<== TFluka::BuildPhysics() called." << endl;
219 void TFluka::ProcessEvent() {
220 if (fVerbosityLevel >=3)
221 cout << "==> TFluka::ProcessEvent() called." << endl;
222 fApplication->GeneratePrimaries();
223 EPISOR.lsouit = true;
225 if (fVerbosityLevel >=3)
226 cout << "<== TFluka::ProcessEvent() called." << endl;
230 void TFluka::ProcessRun(Int_t nevent) {
231 if (fVerbosityLevel >=3)
232 cout << "==> TFluka::ProcessRun(" << nevent << ") called."
235 if (fVerbosityLevel >=2) {
236 cout << "\t* GLOBAL.fdrtr = " << (GLOBAL.lfdrtr?'T':'F') << endl;
237 cout << "\t* Calling flukam again..." << endl;
239 fApplication->InitGeometry();
240 fApplication->BeginEvent();
242 fApplication->FinishEvent();
243 if (fVerbosityLevel >=3)
244 cout << "<== TFluka::ProcessRun(" << nevent << ") called."
249 //_____________________________________________________________________________
250 // methods for building/management of geometry
251 //____________________________________________________________________________
252 // functions from GCONS
253 void TFluka::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
254 Float_t &dens, Float_t &radl, Float_t &absl,
255 Float_t* ubuf, Int_t& nbuf) {
257 fGeometryManager->Gfmate(imat, name, a, z, dens, radl, absl, ubuf, nbuf);
260 void TFluka::Gfmate(Int_t imat, char *name, Double_t &a, Double_t &z,
261 Double_t &dens, Double_t &radl, Double_t &absl,
262 Double_t* ubuf, Int_t& nbuf) {
264 fGeometryManager->Gfmate(imat, name, a, z, dens, radl, absl, ubuf, nbuf);
267 // detector composition
268 void TFluka::Material(Int_t& kmat, const char* name, Double_t a,
269 Double_t z, Double_t dens, Double_t radl, Double_t absl,
270 Float_t* buf, Int_t nwbuf) {
273 ->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
275 void TFluka::Material(Int_t& kmat, const char* name, Double_t a,
276 Double_t z, Double_t dens, Double_t radl, Double_t absl,
277 Double_t* buf, Int_t nwbuf) {
280 ->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
283 void TFluka::Mixture(Int_t& kmat, const char *name, Float_t *a,
284 Float_t *z, Double_t dens, Int_t nlmat, Float_t *wmat) {
287 ->Mixture(kmat, name, a, z, dens, nlmat, wmat);
289 void TFluka::Mixture(Int_t& kmat, const char *name, Double_t *a,
290 Double_t *z, Double_t dens, Int_t nlmat, Double_t *wmat) {
293 ->Mixture(kmat, name, a, z, dens, nlmat, wmat);
296 void TFluka::Medium(Int_t& kmed, const char *name, Int_t nmat,
297 Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd,
298 Double_t stemax, Double_t deemax, Double_t epsil,
299 Double_t stmin, Float_t* ubuf, Int_t nbuf) {
302 ->Medium(kmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax,
303 epsil, stmin, ubuf, nbuf);
305 void TFluka::Medium(Int_t& kmed, const char *name, Int_t nmat,
306 Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd,
307 Double_t stemax, Double_t deemax, Double_t epsil,
308 Double_t stmin, Double_t* ubuf, Int_t nbuf) {
311 ->Medium(kmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax,
312 epsil, stmin, ubuf, nbuf);
315 void TFluka::Matrix(Int_t& krot, Double_t thetaX, Double_t phiX,
316 Double_t thetaY, Double_t phiY, Double_t thetaZ,
320 ->Matrix(krot, thetaX, phiX, thetaY, phiY, thetaZ, phiZ);
323 void TFluka::Gstpar(Int_t itmed, const char *param, Double_t parval) {
325 fGeometryManager->Gstpar(itmed, param, parval);
328 // functions from GGEOM
329 Int_t TFluka::Gsvolu(const char *name, const char *shape, Int_t nmed,
330 Float_t *upar, Int_t np) {
332 // fVolumeMediaMap[TString(name)] = nmed;
333 if (fVerbosityLevel >= 3)
334 printf("TFluka::Gsvolu() name = %s, nmed = %d\n", name, nmed);
336 TClonesArray &lvols = *fVolumeMediaMap;
337 new(lvols[fNVolumes++])
338 FlukaVolume(name, nmed);
339 return fGeometryManager->Gsvolu(name, shape, nmed, upar, np);
341 Int_t TFluka::Gsvolu(const char *name, const char *shape, Int_t nmed,
342 Double_t *upar, Int_t np) {
344 TClonesArray &lvols = *fVolumeMediaMap;
345 new(lvols[fNVolumes++])
346 FlukaVolume(name, nmed);
348 return fGeometryManager->Gsvolu(name, shape, nmed, upar, np);
351 void TFluka::Gsdvn(const char *name, const char *mother, Int_t ndiv,
354 // The medium of the daughter is the one of the mother
355 Int_t volid = TFluka::VolId(mother);
356 Int_t med = TFluka::VolId2Mate(volid);
357 TClonesArray &lvols = *fVolumeMediaMap;
358 new(lvols[fNVolumes++])
359 FlukaVolume(name, med);
360 fGeometryManager->Gsdvn(name, mother, ndiv, iaxis);
363 void TFluka::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
364 Int_t iaxis, Double_t c0i, Int_t numed) {
366 TClonesArray &lvols = *fVolumeMediaMap;
367 new(lvols[fNVolumes++])
368 FlukaVolume(name, numed);
369 fGeometryManager->Gsdvn2(name, mother, ndiv, iaxis, c0i, numed);
372 void TFluka::Gsdvt(const char *name, const char *mother, Double_t step,
373 Int_t iaxis, Int_t numed, Int_t ndvmx) {
375 TClonesArray &lvols = *fVolumeMediaMap;
376 new(lvols[fNVolumes++])
377 FlukaVolume(name, numed);
378 fGeometryManager->Gsdvt(name, mother, step, iaxis, numed, ndvmx);
381 void TFluka::Gsdvt2(const char *name, const char *mother, Double_t step,
382 Int_t iaxis, Double_t c0, Int_t numed, Int_t ndvmx) {
384 TClonesArray &lvols = *fVolumeMediaMap;
385 new(lvols[fNVolumes++])
386 FlukaVolume(name, numed);
387 fGeometryManager->Gsdvt2(name, mother, step, iaxis, c0, numed, ndvmx);
390 void TFluka::Gsord(const char *name, Int_t iax) {
392 fGeometryManager->Gsord(name, iax);
395 void TFluka::Gspos(const char *name, Int_t nr, const char *mother,
396 Double_t x, Double_t y, Double_t z, Int_t irot,
399 fGeometryManager->Gspos(name, nr, mother, x, y, z, irot, konly);
402 void TFluka::Gsposp(const char *name, Int_t nr, const char *mother,
403 Double_t x, Double_t y, Double_t z, Int_t irot,
404 const char *konly, Float_t *upar, Int_t np) {
406 fGeometryManager->Gsposp(name, nr, mother, x, y, z, irot, konly, upar, np);
408 void TFluka::Gsposp(const char *name, Int_t nr, const char *mother,
409 Double_t x, Double_t y, Double_t z, Int_t irot,
410 const char *konly, Double_t *upar, Int_t np) {
412 fGeometryManager->Gsposp(name, nr, mother, x, y, z, irot, konly, upar, np);
415 void TFluka::Gsbool(const char* onlyVolName, const char* manyVolName) {
417 fGeometryManager->Gsbool(onlyVolName, manyVolName);
420 void TFluka::SetCerenkov(Int_t itmed, Int_t npckov, Float_t *ppckov,
421 Float_t *absco, Float_t *effic, Float_t *rindex) {
423 fGeometryManager->SetCerenkov(itmed, npckov, ppckov, absco, effic, rindex);
425 void TFluka::SetCerenkov(Int_t itmed, Int_t npckov, Double_t *ppckov,
426 Double_t *absco, Double_t *effic, Double_t *rindex) {
428 fGeometryManager->SetCerenkov(itmed, npckov, ppckov, absco, effic, rindex);
432 void TFluka::WriteEuclid(const char* fileName, const char* topVol,
433 Int_t number, Int_t nlevel) {
435 fGeometryManager->WriteEuclid(fileName, topVol, number, nlevel);
440 //_____________________________________________________________________________
441 // methods needed by the stepping
442 //____________________________________________________________________________
444 Int_t TFluka::GetMedium() const {
446 // Get the medium number for the current fluka region
448 FGeometryInit* flugg = FGeometryInit::GetInstance();
449 return flugg->GetMedium(fCurrentFlukaRegion);
454 //____________________________________________________________________________
455 // particle table usage
456 // ID <--> PDG transformations
457 //_____________________________________________________________________________
458 Int_t TFluka::IdFromPDG(Int_t pdg) const
461 // Return Fluka code from PDG and pseudo ENDF code
463 // MCIHAD() goes from pdg to fluka internal.
464 Int_t intfluka = mcihad(pdg);
465 // KPTOIP array goes from internal to official
466 return GetFlukaKPTOIP(intfluka);
469 Int_t TFluka::PDGFromId(Int_t id) const
473 // Return PDG code and pseudo ENDF code from Fluka code
475 //IPTOKP array goes from official to internal
479 if (fVerbosityLevel >= 1)
480 printf("\n PDGFromId: Cerenkov Photon \n");
485 if (fVerbosityLevel >= 1)
486 printf("PDGFromId: Error id = 0\n");
490 Int_t intfluka = GetFlukaIPTOKP(id);
492 if (fVerbosityLevel >= 1)
493 printf("PDGFromId: Error intfluka = 0: %d\n", id);
495 } else if (intfluka < 0) {
496 if (fVerbosityLevel >= 1)
497 printf("PDGFromId: Error intfluka < 0: %d\n", id);
500 if (fVerbosityLevel >= 3)
501 printf("mpdgha called with %d %d \n", id, intfluka);
502 return mpdgha(intfluka);
505 //_____________________________________________________________________________
506 // methods for physics management
507 //____________________________________________________________________________
512 void TFluka::SetProcess(const char* flagName, Int_t flagValue)
515 if (iNbOfProc < 100) {
516 for (i=0; i<iNbOfProc; i++) {
517 if (strcmp(&sProcessFlag[i][0],flagName) == 0) {
518 iProcessValue[iNbOfProc] = flagValue;
522 strcpy(&sProcessFlag[iNbOfProc][0],flagName);
523 iProcessValue[iNbOfProc++] = flagValue;
526 cout << "Nb of SetProcess calls exceeds 100 - ignored" << endl;
528 iNbOfProc = iNbOfProc;
531 void TFluka::SetCut(const char* cutName, Double_t cutValue)
534 if (iNbOfCut < 100) {
535 for (i=0; i<iNbOfCut; i++) {
536 if (strcmp(&sCutFlag[i][0],cutName) == 0) {
537 fCutValue[iNbOfCut] = cutValue;
541 strcpy(&sCutFlag[iNbOfCut][0],cutName);
542 fCutValue[iNbOfCut++] = cutValue;
545 cout << "Nb of SetCut calls exceeds 100 - ignored" << endl;
550 Double_t TFluka::Xsec(char*, Double_t, Int_t, Int_t)
552 printf("WARNING: Xsec not yet implemented !\n"); return -1.;
556 void TFluka::InitPhysics()
558 // Last material number taken from the "corealice.inp" file, presently 31
559 // !!! it should be available from Flugg !!!
562 Float_t fLastMaterial = 31.0;
564 // construct file names
565 TString sAliceInp = getenv("ALICE_ROOT");
566 sAliceInp +="/TFluka/input/";
567 TString sAliceCoreInp = sAliceInp;
568 sAliceInp += GetInputFileName();
569 sAliceCoreInp += GetCoreInputFileName();
570 ifstream AliceCoreInp(sAliceCoreInp.Data());
571 ofstream AliceInp(sAliceInp.Data());
573 // copy core input file until (not included) START card
575 Float_t fEventsPerRun;
576 while (AliceCoreInp.getline(sLine,255)) {
577 if (strncmp(sLine,"START",5) != 0)
578 AliceInp << sLine << endl;
580 sscanf(sLine+10,"%10f",&fEventsPerRun);
586 // in G3 the process control values meaning can be different for
587 // different processes, but for most of them is:
588 // 0 process is not activated
589 // 1 process is activated WITH generation of secondaries
590 // 2 process is activated WITHOUT generation of secondaries
591 // if process does not generate secondaries => 1 same as 2
600 // Loop over number of SetProcess calls
601 AliceInp << "*----------------------------------------------------------------------------- ";
603 AliceInp << "*----- The following data are generated from SetProcess and SetCut calls ----- ";
605 AliceInp << "*----------------------------------------------------------------------------- ";
607 for (i=0; i<iNbOfProc; i++) {
610 // G3 default value: 1
611 // G4 processes: G4eplusAnnihilation/G4IeplusAnnihilation
614 // flag = 0 no annihilation
615 // flag = 1 annihilation, decays processed
616 // flag = 2 annihilation, no decay product stored
617 // gMC ->SetProcess("ANNI",1); // EMFCUT -1. 0. 0. 3. lastmat 0. ANNH-THR
618 if (strncmp(&sProcessFlag[i][0],"ANNI",4) == 0) {
619 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
622 AliceInp << "*Kinetic energy threshold (GeV) for e+ annihilation - resets to default=0.";
624 AliceInp << "*Generated from call: SetProcess('ANNI',1) or SetProcess('ANNI',2)";
626 AliceInp << setw(10) << "EMFCUT ";
627 AliceInp << setiosflags(ios::scientific) << setprecision(5);
628 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
629 AliceInp << setw(10) << -1.0; // kinetic energy threshold (GeV) for e+ annihilation (resets to default=0)
630 AliceInp << setw(10) << 0.0; // not used
631 AliceInp << setw(10) << 0.0; // not used
632 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
633 AliceInp << setw(10) << setprecision(2);
634 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
635 AliceInp << setprecision(1);
636 AliceInp << setw(10) << 1.0; // step length in assigning indices
637 AliceInp << setw(8) << "ANNH-THR";
640 else if (iProcessValue[i] == 0) {
643 AliceInp << "*No annihilation - no FLUKA card generated";
645 AliceInp << "*Generated from call: SetProcess('ANNI',0)";
651 AliceInp << "*Illegal flag value in SetProcess('ANNI',?) call.";
653 AliceInp << "*No FLUKA card generated";
658 // bremsstrahlung and pair production are both activated
659 // G3 default value: 1
660 // G4 processes: G4eBremsstrahlung/G4IeBremsstrahlung,
661 // G4MuBremsstrahlung/G4IMuBremsstrahlung,
662 // G4LowEnergyBremstrahlung
663 // Particles: e-/e+; mu+/mu-
665 // flag = 0 no bremsstrahlung
666 // flag = 1 bremsstrahlung, photon processed
667 // flag = 2 bremsstrahlung, no photon stored
668 // gMC ->SetProcess("BREM",1); // PAIRBREM 2. 0. 0. 3. lastmat
669 // EMFCUT -1. 0. 0. 3. lastmat 0. ELPO-THR
670 // G3 default value: 1
671 // G4 processes: G4GammaConversion,
672 // G4MuPairProduction/G4IMuPairProduction
673 // G4LowEnergyGammaConversion
674 // Particles: gamma, mu
676 // flag = 0 no delta rays
677 // flag = 1 delta rays, secondaries processed
678 // flag = 2 delta rays, no secondaries stored
679 // gMC ->SetProcess("PAIR",1); // PAIRBREM 1. 0. 0. 3. lastmat
680 // EMFCUT 0. 0. -1. 3. lastmat 0. PHOT-THR
681 else if ((strncmp(&sProcessFlag[i][0],"PAIR",4) == 0) && (iProcessValue[i] == 1 || iProcessValue[i] == 2)) {
682 for (j=0; j<iNbOfProc; j++) {
683 if ((strncmp(&sProcessFlag[j][0],"BREM",4) == 0) && (iProcessValue[j] == 1 || iProcessValue[j] == 2)) {
686 AliceInp << "*Bremsstrahlung and pair production by muons and charged hadrons both activated";
688 AliceInp << "*Generated from call: SetProcess('BREM',1) and SetProcess('PAIR',1)";
690 AliceInp << "*Energy threshold set by call SetCut('BCUTM',cut) or set to 0.";
692 AliceInp << "*Energy threshold set by call SetCut('PPCUTM',cut) or set to 0.";
694 AliceInp << setw(10) << "PAIRBREM ";
695 AliceInp << setiosflags(ios::scientific) << setprecision(5);
696 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
697 AliceInp << setw(10) << 3.0; // bremsstrahlung and pair production by muons and charged hadrons both are activated
698 // direct pair production by muons
699 // G4 particles: "e-", "e+"
700 // G3 default value: 0.01 GeV
701 //gMC ->SetCut("PPCUTM",cut); // total energy cut for direct pair prod. by muons
703 for (k=0; k<iNbOfCut; k++) {
704 if (strncmp(&sCutFlag[k][0],"PPCUTM",6) == 0) fCut = fCutValue[k];
706 AliceInp << setiosflags(ios::scientific) << setprecision(5);
707 AliceInp << setw(10) << fCut; // e+, e- kinetic energy threshold (in GeV) for explicit pair production.
708 // muon and hadron bremsstrahlung
709 // G4 particles: "gamma"
710 // G3 default value: CUTGAM=0.001 GeV
711 //gMC ->SetCut("BCUTM",cut); // cut for muon and hadron bremsstrahlung
713 for (k=0; k<iNbOfCut; k++) {
714 if (strncmp(&sCutFlag[k][0],"BCUTM",5) == 0) fCut = fCutValue[k];
716 AliceInp << setiosflags(ios::scientific) << setprecision(5);
717 AliceInp << setw(10) << fCut; // photon energy threshold (GeV) for explicit bremsstrahlung production
718 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
719 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
720 AliceInp << setw(10) << setprecision(2);
721 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
727 AliceInp << "*Kinetic energy threshold (GeV) for e+/e- bremsstrahlung - resets to default=0.";
729 AliceInp << "*Generated from call: SetProcess('BREM',1);";
731 AliceInp << setw(10) << "EMFCUT ";
733 for (k=0; k<iNbOfCut; k++) {
734 if (strncmp(&sCutFlag[k][0],"BCUTE",5) == 0) fCut = fCutValue[k];
736 AliceInp << setiosflags(ios::scientific) << setprecision(5);
737 AliceInp << setw(10) << fCut; // kinetic energy threshold (GeV) for e+/e- bremsstrahlung (resets to default=0)
738 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
739 AliceInp << setw(10) << 0.0; // not used
740 AliceInp << setw(10) << 0.0; // not used
741 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
742 AliceInp << setw(10) << setprecision(2);
743 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
744 AliceInp << setprecision(1);
745 AliceInp << setw(10) << 1.0; // step length in assigning indices
746 AliceInp << setw(8) << "ELPO-THR";
752 AliceInp << "*Pair production by electrons is activated";
754 AliceInp << "*Generated from call: SetProcess('PAIR',1);";
756 AliceInp << setw(10) << "EMFCUT ";
757 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
758 AliceInp << setw(10) << 0.0; // energy threshold (GeV) for Compton scattering (= 0.0 : ignored)
759 AliceInp << setw(10) << 0.0; // energy threshold (GeV) for Photoelectric (= 0.0 : ignored)
761 for (j=0; j<iNbOfCut; j++) {
762 if (strncmp(&sCutFlag[j][0],"CUTGAM",6) == 0) fCut = fCutValue[j];
764 AliceInp << setiosflags(ios::scientific) << setprecision(5);
765 AliceInp << setw(10) << fCut; // energy threshold (GeV) for gamma pair production (< 0.0 : resets to default, = 0.0 : ignored)
766 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
767 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
768 AliceInp << setprecision(2);
769 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
770 AliceInp << setprecision(1);
771 AliceInp << setw(10) << 1.0; // step length in assigning indices
772 AliceInp << setw(8) << "PHOT-THR";
775 } // end of if for BREM
776 } // end of loop for BREM
778 // only pair production by muons and charged hadrons is activated
781 AliceInp << "*Pair production by muons and charged hadrons is activated";
783 AliceInp << "*Generated from call: SetProcess('PAIR',1) or SetProcess('PAIR',2)";
785 AliceInp << "*Energy threshold set by call SetCut('PPCUTM',cut) or set to 0.";
787 AliceInp << setw(10) << "PAIRBREM ";
788 AliceInp << setiosflags(ios::scientific) << setprecision(5);
789 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
790 AliceInp << setw(10) << 1.0; // pair production by muons and charged hadrons is activated
791 // direct pair production by muons
792 // G4 particles: "e-", "e+"
793 // G3 default value: 0.01 GeV
794 //gMC ->SetCut("PPCUTM",cut); // total energy cut for direct pair prod. by muons
795 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
796 AliceInp << setw(10) << 0.0; // e+, e- kinetic energy threshold (in GeV) for explicit pair production.
797 AliceInp << setw(10) << 0.0; // no explicit bremsstrahlung production is simulated
798 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
799 AliceInp << setprecision(2);
800 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
806 AliceInp << "*Pair production by electrons is activated";
808 AliceInp << "*Generated from call: SetProcess('PAIR',1) or SetProcess('PAIR',2)";
810 AliceInp << setw(10) << "EMFCUT ";
811 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
812 AliceInp << setw(10) << 0.0; // energy threshold (GeV) for Compton scattering (= 0.0 : ignored)
813 AliceInp << setw(10) << 0.0; // energy threshold (GeV) for Photoelectric (= 0.0 : ignored)
816 for (j=0; j<iNbOfCut; j++) {
817 if (strncmp(&sCutFlag[j][0],"CUTGAM",6) == 0) fCut = fCutValue[j];
819 AliceInp << setiosflags(ios::scientific) << setprecision(5);
820 AliceInp << setw(10) << fCut; // energy threshold (GeV) for gamma pair production (< 0.0 : resets to default, = 0.0 : ignored)
821 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
822 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
823 AliceInp << setprecision(2);
824 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
825 AliceInp << setprecision(1);
826 AliceInp << setw(10) << 1.0; // step length in assigning indices
827 AliceInp << setw(8) << "PHOT-THR";
832 } // end of if for PAIR
837 // G3 default value: 1
838 // G4 processes: G4eBremsstrahlung/G4IeBremsstrahlung,
839 // G4MuBremsstrahlung/G4IMuBremsstrahlung,
840 // G4LowEnergyBremstrahlung
841 // Particles: e-/e+; mu+/mu-
843 // flag = 0 no bremsstrahlung
844 // flag = 1 bremsstrahlung, photon processed
845 // flag = 2 bremsstrahlung, no photon stored
846 // gMC ->SetProcess("BREM",1); // PAIRBREM 2. 0. 0. 3. lastmat
847 // EMFCUT -1. 0. 0. 3. lastmat 0. ELPO-THR
848 else if (strncmp(&sProcessFlag[i][0],"BREM",4) == 0) {
849 for (j=0; j<iNbOfProc; j++) {
850 if ((strncmp(&sProcessFlag[j][0],"PAIR",4) == 0) && iProcessValue[j] == 1) goto NOBREM;
852 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
855 AliceInp << "*Bremsstrahlung by muons and charged hadrons is activated";
857 AliceInp << "*Generated from call: SetProcess('BREM',1) or SetProcess('BREM',2)";
859 AliceInp << "*Energy threshold set by call SetCut('BCUTM',cut) or set to 0.";
861 AliceInp << setw(10) << "PAIRBREM ";
862 AliceInp << setiosflags(ios::scientific) << setprecision(5);
863 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
864 AliceInp << setw(10) << 2.0; // bremsstrahlung by muons and charged hadrons is activated
865 AliceInp << setw(10) << 0.0; // no meaning
866 // muon and hadron bremsstrahlung
867 // G4 particles: "gamma"
868 // G3 default value: CUTGAM=0.001 GeV
869 //gMC ->SetCut("BCUTM",cut); // cut for muon and hadron bremsstrahlung
871 for (j=0; j<iNbOfCut; j++) {
872 if (strncmp(&sCutFlag[j][0],"BCUTM",5) == 0) fCut = fCutValue[j];
874 AliceInp << setw(10) << fCut; // photon energy threshold (GeV) for explicit bremsstrahlung production
875 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
876 AliceInp << setw(10) << setprecision(2);
877 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
883 AliceInp << "*Kinetic energy threshold (GeV) for e+/e- bremsstrahlung - resets to default=0.";
885 AliceInp << "*Generated from call: SetProcess('BREM',1);";
887 AliceInp << setw(10) << "EMFCUT ";
888 AliceInp << setiosflags(ios::scientific) << setprecision(5);
889 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
890 AliceInp << setw(10) << -1.0; // kinetic energy threshold (GeV) for e+/e- bremsstrahlung (resets to default=0)
891 AliceInp << setw(10) << 0.0; // not used
892 AliceInp << setw(10) << 0.0; // not used
893 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
894 AliceInp << setw(10) << setprecision(2);
895 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
896 AliceInp << setprecision(1);
897 AliceInp << setw(10) << 1.0; // step length in assigning indices
898 AliceInp << setw(8) << "ELPO-THR";
901 else if (iProcessValue[i] == 0) {
904 AliceInp << "*No bremsstrahlung - no FLUKA card generated";
906 AliceInp << "*Generated from call: SetProcess('BREM',0)";
912 AliceInp << "*Illegal flag value in SetProcess('BREM',?) call.";
914 AliceInp << "*No FLUKA card generated";
919 } // end of else if (strncmp(&sProcessFlag[i][0],"BREM",4) == 0)
922 // Cerenkov photon generation
923 // G3 default value: 0
924 // G4 process: G4Cerenkov
926 // Particles: charged
928 // flag = 0 no Cerenkov photon generation
929 // flag = 1 Cerenkov photon generation
930 // flag = 2 Cerenkov photon generation with primary stopped at each step
931 //xx gMC ->SetProcess("CKOV",1); // ??? Cerenkov photon generation
932 else if (strncmp(&sProcessFlag[i][0],"CKOV",4) == 0) {
933 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
936 AliceInp << "*Cerenkov photon generation";
938 AliceInp << "*Generated from call: SetProcess('CKOV',1) or SetProcess('CKOV',2)";
940 AliceInp << setw(10) << "OPT-PROD ";
941 AliceInp << setiosflags(ios::scientific) << setprecision(5);
942 AliceInp << setw(10) << 2.07e-9 ; // minimum Cerenkov photon emission energy (in GeV!). Default: 2.07E-9 GeV (corresponding to 600 nm)
943 AliceInp << setw(10) << 4.96e-9; // maximum Cerenkov photon emission energy (in GeV!). Default: 4.96E-9 GeV (corresponding to 250 nm)
944 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
945 AliceInp << setw(10) << 0.0; // not used
946 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
947 AliceInp << setprecision(2);
948 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
949 AliceInp << setprecision(1);
950 AliceInp << setw(10) << 1.0; // step length in assigning indices
951 AliceInp << setw(8) << "CERENKOV";
954 else if (iProcessValue[i] == 0) {
957 AliceInp << "*No Cerenkov photon generation";
959 AliceInp << "*Generated from call: SetProcess('CKOV',0)";
961 AliceInp << setw(10) << "OPT-PROD ";
962 AliceInp << setiosflags(ios::scientific) << setprecision(5);
963 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
964 AliceInp << setw(10) << 0.0; // not used
965 AliceInp << setw(10) << 0.0; // not used
966 AliceInp << setw(10) << 0.0; // not used
967 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
968 AliceInp << setprecision(2);
969 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
970 AliceInp << setprecision(1);
971 AliceInp << setw(10) << 1.0; // step length in assigning indices
972 AliceInp << setw(8) << "CERE-OFF";
978 AliceInp << "*Illegal flag value in SetProcess('CKOV',?) call.";
980 AliceInp << "*No FLUKA card generated";
983 } // end of else if (strncmp(&sProcessFlag[i][0],"CKOV",4) == 0)
986 // Compton scattering
987 // G3 default value: 1
988 // G4 processes: G4ComptonScattering,
989 // G4LowEnergyCompton,
990 // G4PolarizedComptonScattering
993 // flag = 0 no Compton scattering
994 // flag = 1 Compton scattering, electron processed
995 // flag = 2 Compton scattering, no electron stored
996 // gMC ->SetProcess("COMP",1); // EMFCUT -1. 0. 0. 3. lastmat 0. PHOT-THR
997 else if (strncmp(&sProcessFlag[i][0],"COMP",4) == 0) {
998 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
1001 AliceInp << "*Energy threshold (GeV) for Compton scattering - resets to default=0.";
1003 AliceInp << "*Generated from call: SetProcess('COMP',1);";
1005 AliceInp << setw(10) << "EMFCUT ";
1006 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1007 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1008 AliceInp << setw(10) << -1.0; // energy threshold (GeV) for Compton scattering - resets to default=0.
1009 AliceInp << setw(10) << 0.0; // not used
1010 AliceInp << setw(10) << 0.0; // not used
1011 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1012 AliceInp << setprecision(2);
1013 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1014 AliceInp << setprecision(1);
1015 AliceInp << setw(10) << 1.0; // step length in assigning indices
1016 AliceInp << setw(8) << "PHOT-THR";
1019 else if (iProcessValue[i] == 0) {
1022 AliceInp << "*No Compton scattering - no FLUKA card generated";
1024 AliceInp << "*Generated from call: SetProcess('COMP',0)";
1030 AliceInp << "*Illegal flag value in SetProcess('COMP',?) call.";
1032 AliceInp << "*No FLUKA card generated";
1035 } // end of else if (strncmp(&sProcessFlag[i][0],"COMP",4) == 0)
1038 // G3 default value: 1
1039 // G4 process: G4Decay
1041 // Particles: all which decay is applicable for
1043 // flag = 0 no decays
1044 // flag = 1 decays, secondaries processed
1045 // flag = 2 decays, no secondaries stored
1046 //gMC ->SetProcess("DCAY",1); // not available
1047 else if ((strncmp(&sProcessFlag[i][0],"DCAY",4) == 0) && iProcessValue[i] == 1)
1048 cout << "SetProcess for flag=" << &sProcessFlag[i][0] << " value=" << iProcessValue[i] << " not avaliable!" << endl;
1051 // G3 default value: 2
1052 // !! G4 treats delta rays in different way
1053 // G4 processes: G4eIonisation/G4IeIonization,
1054 // G4MuIonisation/G4IMuIonization,
1055 // G4hIonisation/G4IhIonisation
1056 // Particles: charged
1058 // flag = 0 no energy loss
1059 // flag = 1 restricted energy loss fluctuations
1060 // flag = 2 complete energy loss fluctuations
1061 // flag = 3 same as 1
1062 // flag = 4 no energy loss fluctuations
1063 // gMC ->SetProcess("DRAY",0); // DELTARAY 1.E+6 0. 0. 3. lastmat 0.
1064 else if (strncmp(&sProcessFlag[i][0],"DRAY",4) == 0) {
1065 if (iProcessValue[i] == 0 || iProcessValue[i] == 4) {
1068 AliceInp << "*Kinetic energy threshold (GeV) for delta ray production";
1070 AliceInp << "*Generated from call: SetProcess('DRAY',0) or SetProcess('DRAY',4)";
1072 AliceInp << "*No delta ray production by muons - threshold set artificially high";
1074 AliceInp << setw(10) << "DELTARAY ";
1075 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1076 AliceInp << setw(10) << 1.0e+6; // kinetic energy threshold (GeV) for delta ray production (discrete energy transfer)
1077 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1078 AliceInp << setw(10) << 0.0; // ignored
1079 AliceInp << setw(10) << 0.0; // ignored
1080 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1081 AliceInp << setw(10) << setprecision(2);
1082 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1083 AliceInp << setprecision(1);
1084 AliceInp << setw(10) << 1.0; // step length in assigning indices
1087 else if (iProcessValue[i] == 1 || iProcessValue[i] == 2 || iProcessValue[i] == 3) {
1090 AliceInp << "*Kinetic energy threshold (GeV) for delta ray production";
1092 AliceInp << "*Generated from call: SetProcess('DRAY',flag), flag=1,2,3";
1094 AliceInp << "*Delta ray production by muons switched on";
1096 AliceInp << "*Energy threshold set by call SetCut('DCUTM',cut) or set to 0.";
1098 AliceInp << setw(10) << "DELTARAY ";
1099 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1101 for (j=0; j<iNbOfCut; j++) {
1102 if (strncmp(&sCutFlag[j][0],"DCUTM",5) == 0) fCut = fCutValue[j];
1104 AliceInp << setw(10) << fCut; // kinetic energy threshold (GeV) for delta ray production (discrete energy transfer)
1105 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1106 AliceInp << setw(10) << 0.0; // ignored
1107 AliceInp << setw(10) << 0.0; // ignored
1108 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1109 AliceInp << setw(10) << setprecision(2);
1110 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1111 AliceInp << setprecision(1);
1112 AliceInp << setw(10) << 1.0; // step length in assigning indices
1118 AliceInp << "*Illegal flag value in SetProcess('DRAY',?) call.";
1120 AliceInp << "*No FLUKA card generated";
1123 } // end of else if (strncmp(&sProcessFlag[i][0],"DRAY",4) == 0)
1126 // G3 default value: 1
1127 // G4 processes: all defined by TG4PhysicsConstructorHadron
1129 // Particles: hadrons
1131 // flag = 0 no multiple scattering
1132 // flag = 1 hadronic interactions, secondaries processed
1133 // flag = 2 hadronic interactions, no secondaries stored
1134 // gMC ->SetProcess("HADR",1); // ??? hadronic process
1135 //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3) ?????
1136 else if (strncmp(&sProcessFlag[i][0],"HADR",4) == 0) {
1137 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
1140 AliceInp << "*Hadronic interaction is ON by default in FLUKA";
1142 AliceInp << "*No FLUKA card generated";
1145 else if (iProcessValue[i] == 0) {
1148 AliceInp << "*Hadronic interaction is set OFF";
1150 AliceInp << "*Generated from call: SetProcess('HADR',0);";
1152 AliceInp << setw(10) << "MULSOPT ";
1153 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1154 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1155 AliceInp << setw(10) << 0.0; // ignored
1156 AliceInp << setw(10) << 3.0; // multiple scattering for hadrons and muons is completely suppressed
1157 AliceInp << setw(10) << 0.0; // no spin-relativistic corrections
1158 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1159 AliceInp << setprecision(2);
1160 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1167 AliceInp << "*Illegal flag value in SetProcess('HADR',?) call.";
1169 AliceInp << "*No FLUKA card generated";
1172 } // end of else if (strncmp(&sProcessFlag[i][0],"HADR",4) == 0)
1176 // G3 default value: 2
1177 // G4 processes: G4eIonisation/G4IeIonization,
1178 // G4MuIonisation/G4IMuIonization,
1179 // G4hIonisation/G4IhIonisation
1181 // Particles: charged
1183 // flag=0 no energy loss
1184 // flag=1 restricted energy loss fluctuations
1185 // flag=2 complete energy loss fluctuations
1187 // flag=4 no energy loss fluctuations
1188 // If the value ILOSS is changed, then (in G3) cross-sections and energy
1189 // loss tables must be recomputed via the command 'PHYSI'
1190 // gMC ->SetProcess("LOSS",2); // ??? IONFLUCT ? energy loss
1191 else if (strncmp(&sProcessFlag[i][0],"LOSS",4) == 0) {
1192 if (iProcessValue[i] == 2) { // complete energy loss fluctuations
1195 AliceInp << "*Complete energy loss fluctuations do not exist in FLUKA";
1197 AliceInp << "*Generated from call: SetProcess('LOSS',2);";
1199 AliceInp << "*flag=2=complete energy loss fluctuations";
1201 AliceInp << "*No input card generated";
1204 else if (iProcessValue[i] == 1 || iProcessValue[i] == 3) { // restricted energy loss fluctuations
1207 AliceInp << "*Restricted energy loss fluctuations";
1209 AliceInp << "*Generated from call: SetProcess('LOSS',1) or SetProcess('LOSS',3)";
1211 AliceInp << setw(10) << "IONFLUCT ";
1212 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1213 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1214 AliceInp << setw(10) << 1.0; // restricted energy loss fluctuations (for hadrons and muons) switched on
1215 AliceInp << setw(10) << 1.0; // restricted energy loss fluctuations (for e+ and e-) switched on
1216 AliceInp << setw(10) << 1.0; // minimal accuracy
1217 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1218 AliceInp << setprecision(2);
1219 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1222 else if (iProcessValue[i] == 4) { // no energy loss fluctuations
1225 AliceInp << "*No energy loss fluctuations";
1227 AliceInp << "*Generated from call: SetProcess('LOSS',4)";
1229 AliceInp << setw(10) << -1.0; // restricted energy loss fluctuations (for hadrons and muons) switched off
1230 AliceInp << setw(10) << -1.0; // restricted energy loss fluctuations (for e+ and e-) switched off
1231 AliceInp << setw(10) << 1.0; // minimal accuracy
1232 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1233 AliceInp << setprecision(2);
1234 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1240 AliceInp << "*Illegal flag value in SetProcess('LOSS',?) call.";
1242 AliceInp << "*No FLUKA card generated";
1245 } // end of else if (strncmp(&sProcessFlag[i][0],"LOSS",4) == 0)
1248 // multiple scattering
1249 // G3 default value: 1
1250 // G4 process: G4MultipleScattering/G4IMultipleScattering
1252 // Particles: charged
1254 // flag = 0 no multiple scattering
1255 // flag = 1 Moliere or Coulomb scattering
1256 // flag = 2 Moliere or Coulomb scattering
1257 // flag = 3 Gaussian scattering
1258 // gMC ->SetProcess("MULS",1); // MULSOPT multiple scattering
1259 else if (strncmp(&sProcessFlag[i][0],"MULS",4) == 0) {
1260 if (iProcessValue[i] == 1 || iProcessValue[i] == 2 || iProcessValue[i] == 3) {
1263 AliceInp << "*Multiple scattering is ON by default for e+e- and for hadrons/muons";
1265 AliceInp << "*No FLUKA card generated";
1268 else if (iProcessValue[i] == 0) {
1271 AliceInp << "*Multiple scattering is set OFF";
1273 AliceInp << "*Generated from call: SetProcess('MULS',0);";
1275 AliceInp << setw(10) << "MULSOPT ";
1276 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1277 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1278 AliceInp << setw(10) << 0.0; // ignored
1279 AliceInp << setw(10) << 3.0; // multiple scattering for hadrons and muons is completely suppressed
1280 AliceInp << setw(10) << 3.0; // multiple scattering for e+ and e- is completely suppressed
1281 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1282 AliceInp << setprecision(2);
1283 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1289 AliceInp << "*Illegal flag value in SetProcess('MULS',?) call.";
1291 AliceInp << "*No FLUKA card generated";
1294 } // end of else if (strncmp(&sProcessFlag[i][0],"MULS",4) == 0)
1297 // muon nuclear interaction
1298 // G3 default value: 0
1299 // G4 processes: G4MuNuclearInteraction,
1300 // G4MuonMinusCaptureAtRest
1304 // flag = 0 no muon-nuclear interaction
1305 // flag = 1 nuclear interaction, secondaries processed
1306 // flag = 2 nuclear interaction, secondaries not processed
1307 // gMC ->SetProcess("MUNU",1); // MUPHOTON 1. 0. 0. 3. lastmat
1308 else if (strncmp(&sProcessFlag[i][0],"MUNU",4) == 0) {
1309 if (iProcessValue[i] == 1) {
1312 AliceInp << "*Muon nuclear interactions with production of secondary hadrons";
1314 AliceInp << "*Generated from call: SetProcess('MUNU',1);";
1316 AliceInp << setw(10) << "MUPHOTON ";
1317 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1318 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1319 AliceInp << setw(10) << 1.0; // full simulation of muon nuclear interactions and production of secondary hadrons
1320 AliceInp << setw(10) << 0.0; // ratio of longitudinal to transverse virtual photon cross-section - Default = 0.25.
1321 AliceInp << setw(10) << 0.0; // fraction of rho-like interactions ( must be < 1) - Default = 0.75.
1322 AliceInp << setprecision(1);
1323 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1324 AliceInp << setprecision(2);
1325 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1328 else if (iProcessValue[i] == 2) {
1331 AliceInp << "*Muon nuclear interactions without production of secondary hadrons";
1333 AliceInp << "*Generated from call: SetProcess('MUNU',2);";
1335 AliceInp << setw(10) << "MUPHOTON ";
1336 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1337 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1338 AliceInp << setw(10) << 2.0; // full simulation of muon nuclear interactions and production of secondary hadrons
1339 AliceInp << setw(10) << 0.0; // ratio of longitudinal to transverse virtual photon cross-section - Default = 0.25.
1340 AliceInp << setw(10) << 0.0; // fraction of rho-like interactions ( must be < 1) - Default = 0.75.
1341 AliceInp << setprecision(1);
1342 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1343 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1346 else if (iProcessValue[i] == 0) {
1349 AliceInp << "*No muon nuclear interaction - no FLUKA card generated";
1351 AliceInp << "*Generated from call: SetProcess('MUNU',0)";
1357 AliceInp << "*Illegal flag value in SetProcess('MUNU',?) call.";
1359 AliceInp << "*No FLUKA card generated";
1362 } // end of else if (strncmp(&sProcessFlag[i][0],"MUNU",4) == 0)
1366 // G3 default value: 0
1371 // gMC ->SetProcess("PFIS",0); // PHOTONUC -1. 0. 0. 3. lastmat 0.
1372 // flag = 0 no photon fission
1373 // flag = 1 photon fission, secondaries processed
1374 // flag = 2 photon fission, no secondaries stored
1375 else if (strncmp(&sProcessFlag[i][0],"PFIS",4) == 0) {
1376 if (iProcessValue[i] == 0) {
1379 AliceInp << "*No photonuclear interactions";
1381 AliceInp << "*Generated from call: SetProcess('PFIS',0);";
1383 AliceInp << setw(10) << "PHOTONUC ";
1384 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1385 AliceInp << setw(10) << -1.0; // no photonuclear interactions
1386 AliceInp << setw(10) << 0.0; // not used
1387 AliceInp << setw(10) << 0.0; // not used
1388 AliceInp << setw(10) << 3.0; // upper bound of the material indices in which the respective thresholds apply
1389 AliceInp << setprecision(2);
1390 AliceInp << setw(10) << fLastMaterial;
1391 AliceInp << setprecision(1); // upper bound of the material indices in which the respective thresholds apply
1392 AliceInp << setprecision(1);
1393 AliceInp << setw(10) << 1.0; // step length in assigning indices
1396 else if (iProcessValue[i] == 1) {
1399 AliceInp << "*Photon nuclear interactions are activated at all energies";
1401 AliceInp << "*Generated from call: SetProcess('PFIS',1);";
1403 AliceInp << setw(10) << "PHOTONUC ";
1404 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1405 AliceInp << setw(10) << 1.0; // photonuclear interactions are activated at all energies
1406 AliceInp << setw(10) << 0.0; // not used
1407 AliceInp << setw(10) << 0.0; // not used
1408 AliceInp << setprecision(2);
1409 AliceInp << setw(10) << 3.0; // upper bound of the material indices in which the respective thresholds apply
1410 AliceInp << setw(10) << fLastMaterial;
1411 AliceInp << setprecision(1); // upper bound of the material indices in which the respective thresholds apply
1412 AliceInp << setprecision(1);
1413 AliceInp << setw(10) << 1.0; // step length in assigning indices
1416 else if (iProcessValue[i] == 0) {
1419 AliceInp << "*No photofission - no FLUKA card generated";
1421 AliceInp << "*Generated from call: SetProcess('PFIS',0)";
1427 AliceInp << "*Illegal flag value in SetProcess('PFIS',?) call.";
1429 AliceInp << "*No FLUKA card generated";
1435 // photo electric effect
1436 // G3 default value: 1
1437 // G4 processes: G4PhotoElectricEffect
1438 // G4LowEnergyPhotoElectric
1441 // flag = 0 no photo electric effect
1442 // flag = 1 photo electric effect, electron processed
1443 // flag = 2 photo electric effect, no electron stored
1444 // gMC ->SetProcess("PHOT",1); // EMFCUT 0. -1. 0. 3. lastmat 0. PHOT-THR
1445 else if (strncmp(&sProcessFlag[i][0],"PHOT",4) == 0) {
1446 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
1449 AliceInp << "*Photo electric effect is activated";
1451 AliceInp << "*Generated from call: SetProcess('PHOT',1);";
1453 AliceInp << setw(10) << "EMFCUT ";
1454 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1455 AliceInp << setw(10) << 0.0; // ignored
1456 AliceInp << setw(10) << -1.0; // resets to default=0.
1457 AliceInp << setw(10) << 0.0; // ignored
1458 AliceInp << setw(10) << 3.0; // upper bound of the material indices in which the respective thresholds apply
1459 AliceInp << setprecision(2);
1460 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1461 AliceInp << setprecision(1);
1462 AliceInp << setw(10) << 1.0; // step length in assigning indices
1463 AliceInp << setw(8) << "PHOT-THR";
1466 else if (iProcessValue[i] == 0) {
1469 AliceInp << "*No photo electric effect - no FLUKA card generated";
1471 AliceInp << "*Generated from call: SetProcess('PHOT',0)";
1477 AliceInp << "*Illegal flag value in SetProcess('PHOT',?) call.";
1479 AliceInp << "*No FLUKA card generated";
1482 } // else if (strncmp(&sProcessFlag[i][0],"PHOT",4) == 0)
1484 // Rayleigh scattering
1485 // G3 default value: 0
1486 // G4 process: G4OpRayleigh
1488 // Particles: optical photon
1490 // flag = 0 Rayleigh scattering off
1491 // flag = 1 Rayleigh scattering on
1492 //xx gMC ->SetProcess("RAYL",1);
1493 else if (strncmp(&sProcessFlag[i][0],"RAYL",4) == 0) {
1494 if (iProcessValue[i] == 1) {
1497 AliceInp << "*Rayleigh scattering is ON by default in FLUKA";
1499 AliceInp << "*No FLUKA card generated";
1502 else if (iProcessValue[i] == 0) {
1505 AliceInp << "*Rayleigh scattering is set OFF";
1507 AliceInp << "*Generated from call: SetProcess('RAYL',0);";
1509 AliceInp << setw(10) << "EMFRAY ";
1510 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1511 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1512 AliceInp << setw(10) << -1.0; // no Rayleigh scattering and no binding corrections for Compton
1513 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1514 AliceInp << setprecision(2);
1515 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1521 AliceInp << "*Illegal flag value in SetProcess('RAYL',?) call.";
1523 AliceInp << "*No FLUKA card generated";
1526 } // end of else if (strncmp(&sProcessFlag[i][0],"RAYL",4) == 0)
1529 else { // processes not yet treated
1531 // Automatic calculation of tracking medium parameters
1532 // flag = 0 no automatic calculation
1533 // flag = 1 automatic calculation
1534 //xx gMC ->SetProcess("AUTO",1); // ??? automatic computation of the tracking medium parameters
1537 // light photon absorption (Cerenkov photons)
1538 // it is turned on when Cerenkov process is turned on
1539 // G3 default value: 0
1540 // G4 process: G4OpAbsorption, G4OpBoundaryProcess
1542 // Particles: optical photon
1544 // flag = 0 no absorption of Cerenkov photons
1545 // flag = 1 absorption of Cerenkov photons
1546 // gMC ->SetProcess("LABS",2); // ??? Cerenkov light absorption
1549 // To control energy loss fluctuation model
1550 // flag = 0 Urban model
1551 // flag = 1 PAI model
1552 // flag = 2 PAI+ASHO model (not active at the moment)
1553 //xx gMC ->SetProcess("STRA",1); // ??? energy fluctuation model
1555 // synchrotron radiation in magnetic field
1556 // G3 default value: 0
1557 // G4 process: G4SynchrotronRadiation
1561 // flag = 0 no synchrotron radiation
1562 // flag = 1 synchrotron radiation
1563 //xx gMC ->SetProcess("SYNC",1); // ??? synchrotron radiation generation
1565 cout << "SetProcess for flag=" << &sProcessFlag[i][0] << " value=" << iProcessValue[i] << " not yet implemented!" << endl;
1567 } //end of loop number of SetProcess calls
1570 // Loop over number of SetCut calls
1571 for (Int_t i=0; i<iNbOfCut; i++) {
1573 // cuts used in SetProcess calls
1574 if (strncmp(&sCutFlag[i][0],"BCUTM",5) == 0) continue;
1575 else if (strncmp(&sCutFlag[i][0],"BCUTE",5) == 0) continue;
1576 else if (strncmp(&sCutFlag[i][0],"DCUTM",5) == 0) continue;
1577 else if (strncmp(&sCutFlag[i][0],"PPCUTM",6) == 0) continue;
1580 // G4 particles: "gamma"
1581 // G3 default value: 0.001 GeV
1582 //gMC ->SetCut("CUTGAM",cut); // cut for gammas
1583 else if (strncmp(&sCutFlag[i][0],"CUTGAM",6) == 0) {
1586 AliceInp << "*Cut for gamma";
1588 AliceInp << "*Generated from call: SetCut('CUTGAM',cut);";
1590 AliceInp << setw(10) << "PART-THR ";
1591 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1592 AliceInp << setw(10) << -fCutValue[i];
1593 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1594 AliceInp << setw(10) << 7.0;
1599 // G4 particles: "e-"
1601 // G3 default value: 0.001 GeV
1602 //gMC ->SetCut("CUTELE",cut); // cut for e+,e-
1603 else if (strncmp(&sCutFlag[i][0],"CUTELE",6) == 0) {
1606 AliceInp << "*Cut for electrons";
1608 AliceInp << "*Generated from call: SetCut('CUTELE',cut);";
1610 AliceInp << setw(10) << "PART-THR ";
1611 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1612 AliceInp << setw(10) << -fCutValue[i];
1613 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1614 AliceInp << setw(10) << 3.0;
1615 AliceInp << setw(10) << 4.0;
1616 AliceInp << setw(10) << 1.0;
1621 // G4 particles: of type "baryon", "meson", "nucleus" with zero charge
1622 // G3 default value: 0.01 GeV
1623 //gMC ->SetCut("CUTNEU",cut); // cut for neutral hadrons
1624 else if (strncmp(&sCutFlag[i][0],"CUTNEU",6) == 0) {
1627 AliceInp << "*Cut for neutral hadrons";
1629 AliceInp << "*Generated from call: SetCut('CUTNEU',cut);";
1631 AliceInp << setw(10) << "PART-THR ";
1632 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1633 AliceInp << setw(10) << -fCutValue[i];
1634 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1635 AliceInp << setw(10) << 8.0; // Neutron
1636 AliceInp << setw(10) << 9.0; // Antineutron
1638 AliceInp << setw(10) << "PART-THR ";
1639 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1640 AliceInp << setw(10) << -fCutValue[i];
1641 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1642 AliceInp << setw(10) << 12.0; // Kaon zero long
1643 AliceInp << setw(10) << 12.0; // Kaon zero long
1645 AliceInp << setw(10) << "PART-THR ";
1646 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1647 AliceInp << setw(10) << -fCutValue[i];
1648 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1649 AliceInp << setw(10) << 17.0; // Lambda, 18=Antilambda
1650 AliceInp << setw(10) << 19.0; // Kaon zero short
1652 AliceInp << setw(10) << "PART-THR ";
1653 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1654 AliceInp << setw(10) << -fCutValue[i];
1655 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1656 AliceInp << setw(10) << 22.0; // Sigma zero, Pion zero, Kaon zero
1657 AliceInp << setw(10) << 25.0; // Antikaon zero
1659 AliceInp << setw(10) << "PART-THR ";
1660 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1661 AliceInp << setw(10) << -fCutValue[i];
1662 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1663 AliceInp << setw(10) << 32.0; // Antisigma zero
1664 AliceInp << setw(10) << 32.0; // Antisigma zero
1666 AliceInp << setw(10) << "PART-THR ";
1667 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1668 AliceInp << setw(10) << -fCutValue[i];
1669 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1670 AliceInp << setw(10) << 34.0; // Xi zero
1671 AliceInp << setw(10) << 35.0; // AntiXi zero
1673 AliceInp << setw(10) << "PART-THR ";
1674 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1675 AliceInp << setw(10) << -fCutValue[i];
1676 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1677 AliceInp << setw(10) << 47.0; // D zero
1678 AliceInp << setw(10) << 48.0; // AntiD zero
1680 AliceInp << setw(10) << "PART-THR ";
1681 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1682 AliceInp << setw(10) << -fCutValue[i];
1683 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1684 AliceInp << setw(10) << 53.0; // Xi_c zero
1685 AliceInp << setw(10) << 53.0; // Xi_c zero
1687 AliceInp << setw(10) << "PART-THR ";
1688 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1689 AliceInp << setw(10) << -fCutValue[i];
1690 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1691 AliceInp << setw(10) << 55.0; // Xi'_c zero
1692 AliceInp << setw(10) << 56.0; // Omega_c zero
1694 AliceInp << setw(10) << "PART-THR ";
1695 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1696 AliceInp << setw(10) << -fCutValue[i];
1697 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1698 AliceInp << setw(10) << 59.0; // AntiXi_c zero
1699 AliceInp << setw(10) << 59.0; // AntiXi_c zero
1701 AliceInp << setw(10) << "PART-THR ";
1702 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1703 AliceInp << setw(10) << -fCutValue[i];
1704 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1705 AliceInp << setw(10) << 61.0; // AntiXi'_c zero
1706 AliceInp << setw(10) << 62.0; // AntiOmega_c zero
1711 // G4 particles: of type "baryon", "meson", "nucleus" with non-zero charge
1712 // G3 default value: 0.01 GeV
1713 //gMC ->SetCut("CUTHAD",cut); // cut for charged hadrons
1714 else if (strncmp(&sCutFlag[i][0],"CUTHAD",6) == 0) {
1717 AliceInp << "*Cut for charged hadrons";
1719 AliceInp << "*Generated from call: SetCut('CUTHAD',cut);";
1721 AliceInp << setw(10) << "PART-THR ";
1722 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1723 AliceInp << setw(10) << -fCutValue[i];
1724 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1725 AliceInp << setw(10) << 1.0; // Proton
1726 AliceInp << setw(10) << 2.0; // Antiproton
1728 AliceInp << setw(10) << "PART-THR ";
1729 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1730 AliceInp << setw(10) << -fCutValue[i];
1731 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1732 AliceInp << setw(10) << 13.0; // Positive Pion, Negative Pion, Positive Kaon
1733 AliceInp << setw(10) << 16.0; // Negative Kaon
1735 AliceInp << setw(10) << "PART-THR ";
1736 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1737 AliceInp << setw(10) << -fCutValue[i];
1738 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1739 AliceInp << setw(10) << 20.0; // Negative Sigma
1740 AliceInp << setw(10) << 16.0; // Positive Sigma
1742 AliceInp << setw(10) << "PART-THR ";
1743 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1744 AliceInp << setw(10) << -fCutValue[i];
1745 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1746 AliceInp << setw(10) << 31.0; // Antisigma minus
1747 AliceInp << setw(10) << 33.0; // Antisigma plus
1748 AliceInp << setprecision(1);
1749 AliceInp << setw(10) << 2.0; // step length
1751 AliceInp << setw(10) << "PART-THR ";
1752 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1753 AliceInp << setw(10) << -fCutValue[i];
1754 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1755 AliceInp << setw(10) << 36.0; // Negative Xi, Positive Xi, Omega minus
1756 AliceInp << setw(10) << 39.0; // Antiomega
1758 AliceInp << setw(10) << "PART-THR ";
1759 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1760 AliceInp << setw(10) << -fCutValue[i];
1761 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1762 AliceInp << setw(10) << 45.0; // D plus
1763 AliceInp << setw(10) << 46.0; // D minus
1765 AliceInp << setw(10) << "PART-THR ";
1766 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1767 AliceInp << setw(10) << -fCutValue[i];
1768 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1769 AliceInp << setw(10) << 49.0; // D_s plus, D_s minus, Lambda_c plus
1770 AliceInp << setw(10) << 52.0; // Xi_c plus
1772 AliceInp << setw(10) << "PART-THR ";
1773 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1774 AliceInp << setw(10) << -fCutValue[i];
1775 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1776 AliceInp << setw(10) << 54.0; // Xi'_c plus
1777 AliceInp << setw(10) << 60.0; // AntiXi'_c minus
1778 AliceInp << setprecision(1);
1779 AliceInp << setw(10) << 6.0; // step length
1781 AliceInp << setw(10) << "PART-THR ";
1782 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1783 AliceInp << setw(10) << -fCutValue[i];
1784 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1785 AliceInp << setw(10) << 57.0; // Antilambda_c minus
1786 AliceInp << setw(10) << 58.0; // AntiXi_c minus
1791 // G4 particles: "mu+", "mu-"
1792 // G3 default value: 0.01 GeV
1793 //gMC ->SetCut("CUTMUO",cut); // cut for mu+, mu-
1794 else if (strncmp(&sCutFlag[i][0],"CUTMUO",6) == 0) {
1797 AliceInp << "*Cut for muons";
1799 AliceInp << "*Generated from call: SetCut('CUTMUO',cut);";
1801 AliceInp << setw(10) << "PART-THR ";
1802 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1803 AliceInp << setw(10) << -fCutValue[i];
1804 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1805 AliceInp << setprecision(2);
1806 AliceInp << setw(10) << 10.0;
1807 AliceInp << setw(10) << 11.0;
1810 // delta-rays by electrons
1811 // G4 particles: "e-"
1812 // G3 default value: 10**4 GeV
1813 // gMC ->SetCut("DCUTE",cut); // cut for deltarays by electrons ???????????????
1814 else if (strncmp(&sCutFlag[i][0],"DCUTE",5) == 0) {
1817 AliceInp << "*Cut for delta rays by electrons ????????????";
1819 AliceInp << "*Generated from call: SetCut('DCUTE',cut);";
1821 AliceInp << setw(10) << "EMFCUT ";
1822 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1823 AliceInp << setw(10) << -fCutValue[i];
1824 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1825 AliceInp << setw(10) << 0.0;
1826 AliceInp << setw(10) << 0.0;
1827 AliceInp << setw(10) << 3.0;
1828 AliceInp << setprecision(2);
1829 AliceInp << setw(10) << fLastMaterial;
1830 AliceInp << setprecision(1);
1831 AliceInp << setw(10) << 1.0;
1836 // time of flight cut in seconds
1837 // G4 particles: all
1838 // G3 default value: 0.01 GeV
1839 //gMC ->SetCut("TOFMAX",tofmax); // time of flight cuts in seconds
1840 else if (strncmp(&sCutFlag[i][0],"TOFMAX",6) == 0) {
1843 AliceInp << "*Time of flight cuts in seconds";
1845 AliceInp << "*Generated from call: SetCut('TOFMAX',tofmax);";
1847 AliceInp << setw(10) << "TIME-CUT ";
1848 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1849 AliceInp << setw(10) << fCutValue[i]*1.e9;
1850 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1851 AliceInp << setw(10) << 0.0;
1852 AliceInp << setw(10) << 0.0;
1853 AliceInp << setw(10) << -6.0; // lower bound of the particle numbers for which the transport time cut-off and/or the start signal is to be applied
1854 AliceInp << setprecision(2);
1855 AliceInp << setw(10) << 64.0; // upper bound of the particle numbers for which the transport time cut-off and/or the start signal is to be applied
1856 AliceInp << setprecision(1);
1857 AliceInp << setw(10) << 1.0; // step length in assigning numbers
1862 cout << "SetCut for flag=" << &sCutFlag[i][0] << " value=" << fCutValue[i] << " not yet implemented!" << endl;
1864 } //end of loop over SeCut calls
1866 // Add START and STOP card
1867 AliceInp << setw(10) << "START ";
1868 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint);
1869 AliceInp << setw(10) << fEventsPerRun;
1871 AliceInp << setw(10) << "STOP ";
1877 void TFluka::SetMaxStep(Double_t)
1879 // SetMaxStep is dummy procedure in TFluka !
1880 if (fVerbosityLevel >=3)
1881 cout << "SetMaxStep is dummy procedure in TFluka !" << endl;
1884 void TFluka::SetMaxNStep(Int_t)
1886 // SetMaxNStep is dummy procedure in TFluka !
1887 if (fVerbosityLevel >=3)
1888 cout << "SetMaxNStep is dummy procedure in TFluka !" << endl;
1891 void TFluka::SetUserDecay(Int_t)
1893 // SetUserDecay is dummy procedure in TFluka !
1894 if (fVerbosityLevel >=3)
1895 cout << "SetUserDecay is dummy procedure in TFluka !" << endl;
1899 // dynamic properties
1901 void TFluka::TrackPosition(TLorentzVector& position) const
1903 // Return the current position in the master reference frame of the
1904 // track being transported
1905 // TRACKR.atrack = age of the particle
1906 // TRACKR.xtrack = x-position of the last point
1907 // TRACKR.ytrack = y-position of the last point
1908 // TRACKR.ztrack = z-position of the last point
1909 Int_t caller = GetCaller();
1910 if (caller == 1 || caller == 3 || caller == 6 || caller == 11 || caller == 12) { //bxdraw,endraw,usdraw
1911 position.SetX(GetXsco());
1912 position.SetY(GetYsco());
1913 position.SetZ(GetZsco());
1914 position.SetT(TRACKR.atrack);
1916 else if (caller == 4) { // mgdraw
1917 position.SetX(TRACKR.xtrack[TRACKR.ntrack]);
1918 position.SetY(TRACKR.ytrack[TRACKR.ntrack]);
1919 position.SetZ(TRACKR.ztrack[TRACKR.ntrack]);
1920 position.SetT(TRACKR.atrack);
1922 else if (caller == 5) { // sodraw
1923 position.SetX(TRACKR.xtrack[TRACKR.ntrack]);
1924 position.SetY(TRACKR.ytrack[TRACKR.ntrack]);
1925 position.SetZ(TRACKR.ztrack[TRACKR.ntrack]);
1929 Warning("TrackPosition","position not available");
1933 void TFluka::TrackPosition(Double_t& x, Double_t& y, Double_t& z) const
1935 // Return the current position in the master reference frame of the
1936 // track being transported
1937 // TRACKR.atrack = age of the particle
1938 // TRACKR.xtrack = x-position of the last point
1939 // TRACKR.ytrack = y-position of the last point
1940 // TRACKR.ztrack = z-position of the last point
1941 Int_t caller = GetCaller();
1942 if (caller == 1 || caller == 3 || caller == 6 || caller == 11 || caller == 12) { //bxdraw,endraw,usdraw
1947 else if (caller == 4) { // mgdraw
1948 x = TRACKR.xtrack[TRACKR.ntrack];
1949 y = TRACKR.ytrack[TRACKR.ntrack];
1950 z = TRACKR.ztrack[TRACKR.ntrack];
1952 else if (caller == 5) { // sodraw
1953 x = TRACKR.xtrack[TRACKR.ntrack];
1954 y = TRACKR.ytrack[TRACKR.ntrack];
1955 z = TRACKR.ztrack[TRACKR.ntrack];
1958 Warning("TrackPosition","position not available");
1961 void TFluka::TrackMomentum(TLorentzVector& momentum) const
1963 // Return the direction and the momentum (GeV/c) of the track
1964 // currently being transported
1965 // TRACKR.ptrack = momentum of the particle (not always defined, if
1966 // < 0 must be obtained from etrack)
1967 // TRACKR.cx,y,ztrck = direction cosines of the current particle
1968 // TRACKR.etrack = total energy of the particle
1969 // TRACKR.jtrack = identity number of the particle
1970 // PAPROP.am[TRACKR.jtrack] = particle mass in gev
1971 Int_t caller = GetCaller();
1972 if (caller != 2) { // not eedraw
1973 if (TRACKR.ptrack >= 0) {
1974 momentum.SetPx(TRACKR.ptrack*TRACKR.cxtrck);
1975 momentum.SetPy(TRACKR.ptrack*TRACKR.cytrck);
1976 momentum.SetPz(TRACKR.ptrack*TRACKR.cztrck);
1977 momentum.SetE(TRACKR.etrack);
1981 Double_t p = sqrt(TRACKR.etrack*TRACKR.etrack - PAPROP.am[TRACKR.jtrack+6]*PAPROP.am[TRACKR.jtrack+6]);
1982 momentum.SetPx(p*TRACKR.cxtrck);
1983 momentum.SetPy(p*TRACKR.cytrck);
1984 momentum.SetPz(p*TRACKR.cztrck);
1985 momentum.SetE(TRACKR.etrack);
1990 Warning("TrackMomentum","momentum not available");
1993 void TFluka::TrackMomentum(Double_t& px, Double_t& py, Double_t& pz, Double_t& e) const
1995 // Return the direction and the momentum (GeV/c) of the track
1996 // currently being transported
1997 // TRACKR.ptrack = momentum of the particle (not always defined, if
1998 // < 0 must be obtained from etrack)
1999 // TRACKR.cx,y,ztrck = direction cosines of the current particle
2000 // TRACKR.etrack = total energy of the particle
2001 // TRACKR.jtrack = identity number of the particle
2002 // PAPROP.am[TRACKR.jtrack] = particle mass in gev
2003 Int_t caller = GetCaller();
2004 if (caller != 2) { // not eedraw
2005 if (TRACKR.ptrack >= 0) {
2006 px = TRACKR.ptrack*TRACKR.cxtrck;
2007 py = TRACKR.ptrack*TRACKR.cytrck;
2008 pz = TRACKR.ptrack*TRACKR.cztrck;
2013 Double_t p = sqrt(TRACKR.etrack*TRACKR.etrack - PAPROP.am[TRACKR.jtrack+6]*PAPROP.am[TRACKR.jtrack+6]);
2014 px = p*TRACKR.cxtrck;
2015 py = p*TRACKR.cytrck;
2016 pz = p*TRACKR.cztrck;
2022 Warning("TrackMomentum","momentum not available");
2025 Double_t TFluka::TrackStep() const
2027 // Return the length in centimeters of the current step
2028 // TRACKR.ctrack = total curved path
2029 Int_t caller = GetCaller();
2030 if (caller == 1 || caller == 3 || caller == 6) //bxdraw,endraw,usdraw
2032 else if (caller == 4) //mgdraw
2033 return TRACKR.ctrack;
2038 Double_t TFluka::TrackLength() const
2041 // This is the sum of substeps !!!
2042 // TRACKR.ctrack = total curved path of the current step
2043 // Sum of the substeps is identical to TRACKR.ctrack if the is no mag. field
2044 // The sum of all step length starting from the beginning of the track
2045 // for the time being returns only the length in centimeters of the current step
2047 Int_t caller = GetCaller();
2048 if (caller == 1 || caller == 3 || caller == 4 || caller == 6) { //bxdraw,endraw,mgdraw,usdraw
2049 for ( Int_t j=0;j<TRACKR.ntrack;j++) {
2050 sum +=TRACKR.ttrack[j];
2058 Double_t TFluka::TrackTime() const
2060 // Return the current time of flight of the track being transported
2061 // TRACKR.atrack = age of the particle
2062 Int_t caller = GetCaller();
2063 if (caller == 1 || caller == 3 || caller == 4 || caller == 6) //bxdraw,endraw,mgdraw,usdraw
2064 return TRACKR.atrack;
2069 Double_t TFluka::Edep() const
2071 // Energy deposition
2072 // if TRACKR.ntrack = 0, TRACKR.mtrack = 0:
2073 // -->local energy deposition (the value and the point are not recorded in TRACKR)
2074 // but in the variable "rull" of the procedure "endraw.cxx"
2075 // if TRACKR.ntrack > 0, TRACKR.mtrack = 0:
2076 // -->no energy loss along the track
2077 // if TRACKR.ntrack > 0, TRACKR.mtrack > 0:
2078 // -->energy loss distributed along the track
2079 // TRACKR.dtrack = energy deposition of the jth deposition even
2081 for ( Int_t j=0;j<TRACKR.mtrack;j++) {
2082 sum +=TRACKR.dtrack[j];
2084 if (TRACKR.ntrack == 0 && TRACKR.mtrack == 0)
2091 Int_t TFluka::TrackPid() const
2093 // Return the id of the particle transported
2094 // TRACKR.jtrack = identity number of the particle
2095 Int_t caller = GetCaller();
2096 if (caller != 2) // not eedraw
2097 return PDGFromId(TRACKR.jtrack);
2102 Double_t TFluka::TrackCharge() const
2104 // Return charge of the track currently transported
2105 // PAPROP.ichrge = electric charge of the particle
2106 // TRACKR.jtrack = identity number of the particle
2107 Int_t caller = GetCaller();
2108 if (caller != 2) // not eedraw
2109 return PAPROP.ichrge[TRACKR.jtrack+6];
2114 Double_t TFluka::TrackMass() const
2116 // PAPROP.am = particle mass in GeV
2117 // TRACKR.jtrack = identity number of the particle
2118 Int_t caller = GetCaller();
2119 if (caller != 2) // not eedraw
2120 return PAPROP.am[TRACKR.jtrack+6];
2125 Double_t TFluka::Etot() const
2127 // TRACKR.etrack = total energy of the particle
2128 Int_t caller = GetCaller();
2129 if (caller != 2) // not eedraw
2130 return TRACKR.etrack;
2138 Bool_t TFluka::IsNewTrack() const
2141 // True if the track is not at the boundary of the current volume
2142 // Not true in some cases in bxdraw - to be solved
2143 Int_t caller = GetCaller();
2145 return 1; // how to handle double step ?????????????
2147 return 0; // ??????????????
2150 Bool_t TFluka::IsTrackInside() const
2152 // True if the track is not at the boundary of the current volume
2153 // In Fluka a step is always inside one kind of material
2154 // If the step would go behind the region of one material,
2155 // it will be shortened to reach only the boundary.
2156 // Therefore IsTrackInside() is always true.
2157 Int_t caller = GetCaller();
2158 if (caller == 1) // bxdraw
2164 Bool_t TFluka::IsTrackEntering() const
2166 // True if this is the first step of the track in the current volume
2168 Int_t caller = GetCaller();
2169 if (caller == 11) // bxdraw entering
2174 Bool_t TFluka::IsTrackExiting() const
2176 Int_t caller = GetCaller();
2177 if (caller == 12) // bxdraw exiting
2182 Bool_t TFluka::IsTrackOut() const
2184 // True if the track is out of the setup
2186 // Icode = 14: escape - call from Kaskad
2187 // Icode = 23: escape - call from Emfsco
2188 // Icode = 32: escape - call from Kasneu
2189 // Icode = 40: escape - call from Kashea
2190 // Icode = 51: escape - call from Kasoph
2195 iIcode == 51) return 1;
2199 Bool_t TFluka::IsTrackDisappeared() const
2201 // means all inelastic interactions and decays
2202 // iIcode from usdraw
2203 if (iIcode == 101 || // inelastic interaction
2204 iIcode == 102 || // particle decay
2205 iIcode == 214 || // in-flight annihilation
2206 iIcode == 215 || // annihilation at rest
2207 iIcode == 217 || // pair production
2208 iIcode == 221) return 1;
2212 Bool_t TFluka::IsTrackStop() const
2214 // True if the track energy has fallen below the threshold
2215 // means stopped by signal or below energy threshold
2216 // Icode = 12: stopping particle - call from Kaskad
2217 // Icode = 15: time kill - call from Kaskad
2218 // Icode = 21: below threshold, iarg=1 - call from Emfsco
2219 // Icode = 22: below threshold, iarg=2 - call from Emfsco
2220 // Icode = 24: time kill - call from Emfsco
2221 // Icode = 31: below threshold - call from Kasneu
2222 // Icode = 33: time kill - call from Kasneu
2223 // Icode = 41: time kill - call from Kashea
2224 // Icode = 52: time kill - call from Kasoph
2233 iIcode == 52) return 1;
2237 Bool_t TFluka::IsTrackAlive() const
2239 // means not disappeared or not out
2240 if (IsTrackDisappeared() || IsTrackOut() ) return 0;
2248 Int_t TFluka::NSecondaries() const
2249 // Number of secondary particles generated in the current step
2250 // FINUC.np = number of secondaries except light and heavy ions
2251 // FHEAVY.npheav = number of secondaries for light and heavy secondary ions
2253 Int_t caller = GetCaller();
2254 if (caller == 6) // valid only after usdraw
2255 return FINUC.np + FHEAVY.npheav;
2258 } // end of NSecondaries
2260 void TFluka::GetSecondary(Int_t isec, Int_t& particleId,
2261 TLorentzVector& position, TLorentzVector& momentum)
2263 Int_t caller = GetCaller();
2264 if (caller == 6) { // valid only after usdraw
2265 if (isec >= 0 && isec < FINUC.np) {
2266 particleId = PDGFromId(FINUC.kpart[isec]);
2267 position.SetX(fXsco);
2268 position.SetY(fYsco);
2269 position.SetZ(fZsco);
2270 position.SetT(TRACKR.atrack);
2271 // position.SetT(TRACKR.atrack+FINUC.agesec[isec]); //not yet implem.
2272 momentum.SetPx(FINUC.plr[isec]*FINUC.cxr[isec]);
2273 momentum.SetPy(FINUC.plr[isec]*FINUC.cyr[isec]);
2274 momentum.SetPz(FINUC.plr[isec]*FINUC.czr[isec]);
2275 momentum.SetE(FINUC.tki[isec] + PAPROP.am[FINUC.kpart[isec]+6]);
2277 else if (isec >= FINUC.np && isec < FINUC.np + FHEAVY.npheav) {
2278 Int_t jsec = isec - FINUC.np;
2279 particleId = FHEAVY.kheavy[jsec]; // this is Fluka id !!!
2280 position.SetX(fXsco);
2281 position.SetY(fYsco);
2282 position.SetZ(fZsco);
2283 position.SetT(TRACKR.atrack);
2284 // position.SetT(TRACKR.atrack+FHEAVY.agheav[jsec]); //not yet implem.
2285 momentum.SetPx(FHEAVY.pheavy[jsec]*FHEAVY.cxheav[jsec]);
2286 momentum.SetPy(FHEAVY.pheavy[jsec]*FHEAVY.cyheav[jsec]);
2287 momentum.SetPz(FHEAVY.pheavy[jsec]*FHEAVY.czheav[jsec]);
2288 if (FHEAVY.tkheav[jsec] >= 3 && FHEAVY.tkheav[jsec] <= 6)
2289 momentum.SetE(FHEAVY.tkheav[jsec] + PAPROP.am[jsec+6]);
2290 else if (FHEAVY.tkheav[jsec] > 6)
2291 momentum.SetE(FHEAVY.tkheav[jsec] + FHEAVY.amnhea[jsec]); // to be checked !!!
2294 Warning("GetSecondary","isec out of range");
2297 Warning("GetSecondary","no secondaries available");
2298 } // end of GetSecondary
2300 TMCProcess TFluka::ProdProcess(Int_t isec) const
2301 // Name of the process that has produced the secondary particles
2302 // in the current step
2304 const TMCProcess kIpNoProc = kPNoProcess;
2305 const TMCProcess kIpPDecay = kPDecay;
2306 const TMCProcess kIpPPair = kPPair;
2307 // const TMCProcess kIpPPairFromPhoton = kPPairFromPhoton;
2308 // const TMCProcess kIpPPairFromVirtualPhoton = kPPairFromVirtualPhoton;
2309 const TMCProcess kIpPCompton = kPCompton;
2310 const TMCProcess kIpPPhotoelectric = kPPhotoelectric;
2311 const TMCProcess kIpPBrem = kPBrem;
2312 // const TMCProcess kIpPBremFromHeavy = kPBremFromHeavy;
2313 // const TMCProcess kIpPBremFromElectronOrPositron = kPBremFromElectronOrPositron;
2314 const TMCProcess kIpPDeltaRay = kPDeltaRay;
2315 // const TMCProcess kIpPMoller = kPMoller;
2316 // const TMCProcess kIpPBhabha = kPBhabha;
2317 const TMCProcess kIpPAnnihilation = kPAnnihilation;
2318 // const TMCProcess kIpPAnnihilInFlight = kPAnnihilInFlight;
2319 // const TMCProcess kIpPAnnihilAtRest = kPAnnihilAtRest;
2320 const TMCProcess kIpPHadronic = kPHadronic;
2321 const TMCProcess kIpPMuonNuclear = kPMuonNuclear;
2322 const TMCProcess kIpPPhotoFission = kPPhotoFission;
2323 const TMCProcess kIpPRayleigh = kPRayleigh;
2324 // const TMCProcess kIpPCerenkov = kPCerenkov;
2325 // const TMCProcess kIpPSynchrotron = kPSynchrotron;
2327 Int_t mugamma = TRACKR.jtrack == 7 || TRACKR.jtrack == 10 || TRACKR.jtrack == 11;
2328 if (iIcode == 102) return kIpPDecay;
2329 else if (iIcode == 104 || iIcode == 217) return kIpPPair;
2330 // else if (iIcode == 104) return kIpPairFromPhoton;
2331 // else if (iIcode == 217) return kIpPPairFromVirtualPhoton;
2332 else if (iIcode == 219) return kIpPCompton;
2333 else if (iIcode == 221) return kIpPPhotoelectric;
2334 else if (iIcode == 105 || iIcode == 208) return kIpPBrem;
2335 // else if (iIcode == 105) return kIpPBremFromHeavy;
2336 // else if (iIcode == 208) return kPBremFromElectronOrPositron;
2337 else if (iIcode == 103 || iIcode == 400) return kIpPDeltaRay;
2338 else if (iIcode == 210 || iIcode == 212) return kIpPDeltaRay;
2339 // else if (iIcode == 210) return kIpPMoller;
2340 // else if (iIcode == 212) return kIpPBhabha;
2341 else if (iIcode == 214 || iIcode == 215) return kIpPAnnihilation;
2342 // else if (iIcode == 214) return kIpPAnnihilInFlight;
2343 // else if (iIcode == 215) return kIpPAnnihilAtRest;
2344 else if (iIcode == 101) return kIpPHadronic;
2345 else if (iIcode == 101) {
2346 if (!mugamma) return kIpPHadronic;
2347 else if (TRACKR.jtrack == 7) return kIpPPhotoFission;
2348 else return kIpPMuonNuclear;
2350 else if (iIcode == 225) return kIpPRayleigh;
2351 // Fluka codes 100, 300 and 400 still to be investigasted
2352 else return kIpNoProc;
2355 //Int_t StepProcesses(TArrayI &proc) const
2356 // Return processes active in the current step
2358 //ck = total energy of the particl ????????????????
2362 Int_t TFluka::VolId2Mate(Int_t id) const
2365 // Returns the material number for a given volume ID
2367 if (fVerbosityLevel >= 3)
2368 printf("VolId2Mate %d %d\n", id, fMediaByRegion[id-1]);
2369 return fMediaByRegion[id-1];
2372 const char* TFluka::VolName(Int_t id) const
2375 // Returns the volume name for a given volume ID
2377 FlukaVolume* vol = dynamic_cast<FlukaVolume*>((*fVolumeMediaMap)[id-1]);
2378 const char* name = vol->GetName();
2379 if (fVerbosityLevel >= 3)
2380 printf("VolName %d %s \n", id, name);
2384 Int_t TFluka::VolId(const Text_t* volName) const
2387 // Converts from volume name to volume ID.
2388 // Time consuming. (Only used during set-up)
2389 // Could be replaced by hash-table
2393 for (i = 0; i < fNVolumes; i++)
2395 FlukaVolume* vol = dynamic_cast<FlukaVolume*>((*fVolumeMediaMap)[i]);
2396 TString name = vol->GetName();
2397 strcpy(tmp, name.Data());
2399 if (!strcmp(tmp, volName)) break;
2407 Int_t TFluka::CurrentVolID(Int_t& copyNo) const
2410 // Return the logical id and copy number corresponding to the current fluka region
2412 int ir = fCurrentFlukaRegion;
2413 int id = (FGeometryInit::GetInstance())->CurrentVolID(ir, copyNo);
2415 if (fVerbosityLevel >= 3)
2416 printf("CurrentVolID: %d %d %d \n", ir, id, copyNo);
2420 Int_t TFluka::CurrentVolOffID(Int_t off, Int_t& copyNo) const
2423 // Return the logical id and copy number of off'th mother
2424 // corresponding to the current fluka region
2427 return CurrentVolID(copyNo);
2429 int ir = fCurrentFlukaRegion;
2430 int id = (FGeometryInit::GetInstance())->CurrentVolOffID(ir, off, copyNo);
2432 if (fVerbosityLevel >= 3)
2433 printf("CurrentVolOffID: %d %d %d \n", ir, id, copyNo);
2435 if (fVerbosityLevel >= 0)
2436 printf("CurrentVolOffID: Warning Mother not found !!!\n");
2441 const char* TFluka::CurrentVolName() const
2444 // Return the current volume name
2447 Int_t id = TFluka::CurrentVolID(copy);
2448 const char* name = TFluka::VolName(id);
2449 if (fVerbosityLevel >= 3)
2450 printf("CurrentVolumeName: %d %s \n", fCurrentFlukaRegion, name);
2454 const char* TFluka::CurrentVolOffName(Int_t off) const
2457 // Return the volume name of the off'th mother of the current volume
2460 Int_t id = TFluka::CurrentVolOffID(off, copy);
2461 const char* name = TFluka::VolName(id);
2462 if (fVerbosityLevel >= 3)
2463 printf("CurrentVolumeOffName: %d %s \n", fCurrentFlukaRegion, name);
2467 Int_t TFluka::CurrentMaterial(Float_t &a, Float_t &z,
2468 Float_t &dens, Float_t &radl, Float_t &absl) const
2471 // Return the current medium number
2474 Int_t id = TFluka::CurrentVolID(copy);
2475 Int_t med = TFluka::VolId2Mate(id);
2476 if (fVerbosityLevel >= 3)
2477 printf("CurrentMaterial: %d %d \n", fCurrentFlukaRegion, med);
2481 void TFluka::Gmtod(Float_t* xm, Float_t* xd, Int_t iflag)
2483 // Transforms a position from the world reference frame
2484 // to the current volume reference frame.
2486 // Geant3 desription:
2487 // ==================
2488 // Computes coordinates XD (in DRS)
2489 // from known coordinates XM in MRS
2490 // The local reference system can be initialized by
2491 // - the tracking routines and GMTOD used in GUSTEP
2492 // - a call to GMEDIA(XM,NUMED)
2493 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2494 // (inverse routine is GDTOM)
2496 // If IFLAG=1 convert coordinates
2497 // IFLAG=2 convert direction cosinus
2500 Double_t xmD[3], xdD[3];
2501 xmD[0] = xm[0]; xmD[1] = xm[1]; xmD[2] = xm[2];
2502 (FGeometryInit::GetInstance())->Gmtod(xmD, xdD, iflag);
2503 xd[0] = xdD[0]; xd[1] = xdD[1]; xd[2] = xdD[2];
2507 void TFluka::Gmtod(Double_t* xm, Double_t* xd, Int_t iflag)
2509 // Transforms a position from the world reference frame
2510 // to the current volume reference frame.
2512 // Geant3 desription:
2513 // ==================
2514 // Computes coordinates XD (in DRS)
2515 // from known coordinates XM in MRS
2516 // The local reference system can be initialized by
2517 // - the tracking routines and GMTOD used in GUSTEP
2518 // - a call to GMEDIA(XM,NUMED)
2519 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2520 // (inverse routine is GDTOM)
2522 // If IFLAG=1 convert coordinates
2523 // IFLAG=2 convert direction cosinus
2526 Double_t xmD[3], xdD[3];
2527 xdD[0] = xd[0]; xdD[1] = xd[1]; xdD[2] = xd[2];
2528 (FGeometryInit::GetInstance())->Gdtom(xmD, xdD, iflag);
2529 xm[0] = xmD[0]; xm[1] = xmD[1]; xm[2] = xmD[2];
2532 void TFluka::Gdtom(Float_t* xd, Float_t* xm, Int_t iflag)
2534 // Transforms a position from the current volume reference frame
2535 // to the world reference frame.
2537 // Geant3 desription:
2538 // ==================
2539 // Computes coordinates XM (Master Reference System
2540 // knowing the coordinates XD (Detector Ref System)
2541 // The local reference system can be initialized by
2542 // - the tracking routines and GDTOM used in GUSTEP
2543 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2544 // (inverse routine is GMTOD)
2546 // If IFLAG=1 convert coordinates
2547 // IFLAG=2 convert direction cosinus
2553 void TFluka::Gdtom(Double_t* xd, Double_t* xm, Int_t iflag)
2555 // Transforms a position from the current volume reference frame
2556 // to the world reference frame.
2558 // Geant3 desription:
2559 // ==================
2560 // Computes coordinates XM (Master Reference System
2561 // knowing the coordinates XD (Detector Ref System)
2562 // The local reference system can be initialized by
2563 // - the tracking routines and GDTOM used in GUSTEP
2564 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2565 // (inverse routine is GMTOD)
2567 // If IFLAG=1 convert coordinates
2568 // IFLAG=2 convert direction cosinus
2572 (FGeometryInit::GetInstance())->Gdtom(xm, xd, iflag);
2575 // ===============================================================
2576 void TFluka::FutoTest()
2578 Int_t icode, mreg, newreg, particleId;
2579 Double_t rull, xsco, ysco, zsco;
2580 TLorentzVector position, momentum;
2583 if (fVerbosityLevel >=3)
2584 cout << " icode=" << icode << endl;
2585 } else if (icode > 0 && icode <= 5) {
2588 if (fVerbosityLevel >=3)
2589 cout << " icode=" << icode
2592 TrackPosition(position);
2593 TrackMomentum(momentum);
2594 if (fVerbosityLevel >=3) {
2595 cout << "TLorentzVector positionX=" << position.X()
2596 << "positionY=" << position.Y()
2597 << "positionZ=" << position.Z()
2598 << "timeT=" << position.T() << endl;
2599 cout << "TLorentzVector momentumX=" << momentum.X()
2600 << "momentumY=" << momentum.Y()
2601 << "momentumZ=" << momentum.Z()
2602 << "energyE=" << momentum.E() << endl;
2603 cout << "TrackStep=" << TrackStep() << endl;
2604 cout << "TrackLength=" << TrackLength() << endl;
2605 cout << "TrackTime=" << TrackTime() << endl;
2606 cout << "Edep=" << Edep() << endl;
2607 cout << "TrackPid=" << TrackPid() << endl;
2608 cout << "TrackCharge=" << TrackCharge() << endl;
2609 cout << "TrackMass=" << TrackMass() << endl;
2610 cout << "Etot=" << Etot() << endl;
2611 cout << "IsNewTrack=" << IsNewTrack() << endl;
2612 cout << "IsTrackInside=" << IsTrackInside() << endl;
2613 cout << "IsTrackEntering=" << IsTrackEntering() << endl;
2614 cout << "IsTrackExiting=" << IsTrackExiting() << endl;
2615 cout << "IsTrackOut=" << IsTrackOut() << endl;
2616 cout << "IsTrackDisappeared=" << IsTrackDisappeared() << endl;
2617 cout << "IsTrackAlive=" << IsTrackAlive() << endl;
2620 Float_t x = position.X();
2621 Float_t y = position.Y();
2622 Float_t z = position.Z();
2625 xm[0] = x; xm[1] = y; xm[2] = z;
2626 if (fVerbosityLevel >= 3)
2627 printf("Global trackPosition: %f %f %f \n", x, y, z);
2629 if (fVerbosityLevel >= 3)
2630 printf("Local trackPosition: %f %f %f \n", xd[0], xd[1], xd[2]);
2632 if (fVerbosityLevel >= 3)
2633 printf("New trackPosition: %f %f %f \n", xm[0], xm[1], xm[2]);
2634 } else if((icode >= 10 && icode <= 15) ||
2635 (icode >= 20 && icode <= 24) ||
2636 (icode >= 30 && icode <= 33) ||
2637 (icode >= 40 && icode <= 41) ||
2638 (icode >= 50 && icode <= 52)) {
2646 if (fVerbosityLevel >=3) {
2647 cout << " icode=" << icode
2652 << " zsco=" << zsco << endl;
2654 TrackPosition(position);
2655 TrackMomentum(momentum);
2656 if (fVerbosityLevel >=3) {
2657 cout << "Edep=" << Edep() << endl;
2658 cout << "Etot=" << Etot() << endl;
2659 cout << "TrackPid=" << TrackPid() << endl;
2660 cout << "TrackCharge=" << TrackCharge() << endl;
2661 cout << "TrackMass=" << TrackMass() << endl;
2662 cout << "IsTrackOut=" << IsTrackOut() << endl;
2663 cout << "IsTrackDisappeared=" << IsTrackDisappeared() << endl;
2664 cout << "IsTrackStop=" << IsTrackStop() << endl;
2665 cout << "IsTrackAlive=" << IsTrackAlive() << endl;
2667 } else if((icode >= 100 && icode <= 105) ||
2671 (icode >= 214 && icode <= 215) ||
2684 if (fVerbosityLevel >=3) {
2685 cout << " icode=" << icode
2689 << " zsco=" << zsco << endl;
2690 cout << "TrackPid=" << TrackPid() << endl;
2691 cout << "NSecondaries=" << NSecondaries() << endl;
2694 for (Int_t isec=0; isec< NSecondaries(); isec++) {
2695 TFluka::GetSecondary(isec, particleId, position, momentum);
2696 if (fVerbosityLevel >=3) {
2697 cout << "TLorentzVector positionX=" << position.X()
2698 << "positionY=" << position.Y()
2699 << "positionZ=" << position.Z()
2700 << "timeT=" << position.T() << endl;
2701 cout << "TLorentzVector momentumX=" << momentum.X()
2702 << "momentumY=" << momentum.Y()
2703 << "momentumZ=" << momentum.Z()
2704 << "energyE=" << momentum.E() << endl;
2705 cout << "TrackPid=" << particleId << endl;
2708 } else if((icode == 19) ||
2714 newreg = GetNewreg();
2718 if (fVerbosityLevel >=3) {
2719 cout << " icode=" << icode
2721 << " newreg=" << newreg
2724 << " zsco=" << zsco << endl;
2727 } // end of FutoTest