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"
21 #include "TFlukaGeo.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 "TGeoManager.h"
35 #include "TFlukaMCGeometry.h"
37 #include "TLorentzVector.h"
39 // Fluka methods that may be needed.
41 # define flukam flukam_
42 # define fluka_openinp fluka_openinp_
43 # define fluka_closeinp fluka_closeinp_
44 # define mcihad mcihad_
45 # define mpdgha mpdgha_
47 # define flukam FLUKAM
48 # define fluka_openinp FLUKA_OPENINP
49 # define fluka_closeinp FLUKA_CLOSEINP
50 # define mcihad MCIHAD
51 # define mpdgha MPDGHA
57 // Prototypes for FLUKA functions
59 void type_of_call flukam(const int&);
60 void type_of_call fluka_openinp(const int&, DEFCHARA);
61 void type_of_call fluka_closeinp(const int&);
62 int type_of_call mcihad(const int&);
63 int type_of_call mpdgha(const int&);
67 // Class implementation for ROOT
72 //----------------------------------------------------------------------------
73 // TFluka constructors and destructors.
74 //______________________________________________________________________________
81 // Default constructor
84 fCurrentFlukaRegion = -1;
88 //______________________________________________________________________________
89 TFluka::TFluka(const char *title, Int_t verbosity, Bool_t isRootGeometrySupported)
90 :TVirtualMC("TFluka",title, isRootGeometrySupported),
91 fVerbosityLevel(verbosity),
96 // create geometry interface
97 if (fVerbosityLevel >=3)
98 cout << "<== TFluka::TFluka(" << title << ") constructor called." << endl;
101 fCurrentFlukaRegion = -1;
102 fGeom = new TFlukaMCGeometry("geom", "ALICE geometry");
105 //______________________________________________________________________________
107 if (fVerbosityLevel >=3)
108 cout << "==> TFluka::~TFluka() destructor called." << endl;
112 if (fVerbosityLevel >=3)
113 cout << "<== TFluka::~TFluka() destructor called." << endl;
117 //______________________________________________________________________________
118 // TFluka control methods
119 //______________________________________________________________________________
120 void TFluka::Init() {
122 if (fVerbosityLevel >=3)
123 cout << "==> TFluka::Init() called." << endl;
125 if (!gGeoManager) new TGeoManager("geom", "FLUKA geometry");
126 fApplication->ConstructGeometry();
127 TGeoVolume *top = (TGeoVolume*)gGeoManager->GetListOfVolumes()->First();
128 gGeoManager->SetTopVolume(top);
129 gGeoManager->CloseGeometry("di");
130 gGeoManager->DefaultColors(); // to be removed
131 fNVolumes = fGeom->NofVolumes();
132 printf("== Number of volumes: %i\n ==", fNVolumes);
133 fGeom->CreateFlukaMatFile("flukaMat.inp");
134 cout << "\t* InitPhysics() - Prepare input file to be called" << endl;
135 // now we have TGeo geometry created and we have to patch alice.inp
136 // with the material mapping file FlukaMat.inp
137 InitPhysics(); // prepare input file with the current physics settings
138 cout << "\t* InitPhysics() - Prepare input file was called" << endl;
140 if (fVerbosityLevel >=2)
141 cout << "\t* Changing lfdrtr = (" << (GLOBAL.lfdrtr?'T':'F')
142 << ") in fluka..." << endl;
143 GLOBAL.lfdrtr = true;
145 if (fVerbosityLevel >=2)
146 cout << "\t* Opening file " << sInputFileName << endl;
147 const char* fname = sInputFileName;
148 fluka_openinp(lunin, PASSCHARA(fname));
150 if (fVerbosityLevel >=2)
151 cout << "\t* Calling flukam..." << endl;
154 if (fVerbosityLevel >=2)
155 cout << "\t* Closing file " << sInputFileName << endl;
156 fluka_closeinp(lunin);
160 if (fVerbosityLevel >=3)
161 cout << "<== TFluka::Init() called." << endl;
164 //______________________________________________________________________________
165 void TFluka::FinishGeometry() {
167 // Build-up table with region to medium correspondance
169 if (fVerbosityLevel >=3)
170 cout << "==> TFluka::FinishGeometry() called." << endl;
172 printf("----FinishGeometry - nothing to do with TGeo\n");
174 if (fVerbosityLevel >=3)
175 cout << "<== TFluka::FinishGeometry() called." << endl;
178 //______________________________________________________________________________
179 void TFluka::BuildPhysics() {
180 if (fVerbosityLevel >=3)
181 cout << "==> TFluka::BuildPhysics() called." << endl;
184 if (fVerbosityLevel >=3)
185 cout << "<== TFluka::BuildPhysics() called." << endl;
188 //______________________________________________________________________________
189 void TFluka::ProcessEvent() {
190 if (fVerbosityLevel >=3)
191 cout << "==> TFluka::ProcessEvent() called." << endl;
192 fApplication->GeneratePrimaries();
193 EPISOR.lsouit = true;
195 if (fVerbosityLevel >=3)
196 cout << "<== TFluka::ProcessEvent() called." << endl;
199 //______________________________________________________________________________
200 void TFluka::ProcessRun(Int_t nevent) {
201 if (fVerbosityLevel >=3)
202 cout << "==> TFluka::ProcessRun(" << nevent << ") called."
205 if (fVerbosityLevel >=2) {
206 cout << "\t* GLOBAL.fdrtr = " << (GLOBAL.lfdrtr?'T':'F') << endl;
207 cout << "\t* Calling flukam again..." << endl;
209 fApplication->InitGeometry();
210 fApplication->BeginEvent();
212 fApplication->FinishEvent();
213 if (fVerbosityLevel >=3)
214 cout << "<== TFluka::ProcessRun(" << nevent << ") called."
219 //_____________________________________________________________________________
220 // methods for building/management of geometry
222 // functions from GCONS
223 //____________________________________________________________________________
224 void TFluka::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
225 Float_t &dens, Float_t &radl, Float_t &absl,
226 Float_t* ubuf, Int_t& nbuf) {
228 fGeom->Gfmate(imat, name, a, z, dens, radl, absl, ubuf, nbuf);
231 //______________________________________________________________________________
232 void TFluka::Gfmate(Int_t imat, char *name, Double_t &a, Double_t &z,
233 Double_t &dens, Double_t &radl, Double_t &absl,
234 Double_t* ubuf, Int_t& nbuf) {
236 fGeom->Gfmate(imat, name, a, z, dens, radl, absl, ubuf, nbuf);
239 // detector composition
240 //______________________________________________________________________________
241 void TFluka::Material(Int_t& kmat, const char* name, Double_t a,
242 Double_t z, Double_t dens, Double_t radl, Double_t absl,
243 Float_t* buf, Int_t nwbuf) {
245 fGeom->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
248 //______________________________________________________________________________
249 void TFluka::Material(Int_t& kmat, const char* name, Double_t a,
250 Double_t z, Double_t dens, Double_t radl, Double_t absl,
251 Double_t* buf, Int_t nwbuf) {
253 fGeom->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
256 //______________________________________________________________________________
257 void TFluka::Mixture(Int_t& kmat, const char *name, Float_t *a,
258 Float_t *z, Double_t dens, Int_t nlmat, Float_t *wmat) {
260 fGeom->Mixture(kmat, name, a, z, dens, nlmat, wmat);
263 //______________________________________________________________________________
264 void TFluka::Mixture(Int_t& kmat, const char *name, Double_t *a,
265 Double_t *z, Double_t dens, Int_t nlmat, Double_t *wmat) {
267 fGeom->Mixture(kmat, name, a, z, dens, nlmat, wmat);
270 //______________________________________________________________________________
271 void TFluka::Medium(Int_t& kmed, const char *name, Int_t nmat,
272 Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd,
273 Double_t stemax, Double_t deemax, Double_t epsil,
274 Double_t stmin, Float_t* ubuf, Int_t nbuf) {
276 fGeom->Medium(kmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax,
277 epsil, stmin, ubuf, nbuf);
280 //______________________________________________________________________________
281 void TFluka::Medium(Int_t& kmed, const char *name, Int_t nmat,
282 Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd,
283 Double_t stemax, Double_t deemax, Double_t epsil,
284 Double_t stmin, Double_t* ubuf, Int_t nbuf) {
286 fGeom->Medium(kmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax,
287 epsil, stmin, ubuf, nbuf);
290 //______________________________________________________________________________
291 void TFluka::Matrix(Int_t& krot, Double_t thetaX, Double_t phiX,
292 Double_t thetaY, Double_t phiY, Double_t thetaZ,
295 fGeom->Matrix(krot, thetaX, phiX, thetaY, phiY, thetaZ, phiZ);
298 //______________________________________________________________________________
299 void TFluka::Gstpar(Int_t /*itmed*/, const char */*param*/, Double_t /*parval*/) {
301 // Is it needed with TGeo ??? - to clear-up
302 Warning("Gstpar", "Not implemented with TGeo");
305 // functions from GGEOM
306 //_____________________________________________________________________________
307 void TFluka::Gsatt(const char *name, const char *att, Int_t val)
309 fGeom->Gsatt(name,att, val);
312 //______________________________________________________________________________
313 Int_t TFluka::Gsvolu(const char *name, const char *shape, Int_t nmed,
314 Float_t *upar, Int_t np) {
316 return fGeom->Gsvolu(name, shape, nmed, upar, np);
319 //______________________________________________________________________________
320 Int_t TFluka::Gsvolu(const char *name, const char *shape, Int_t nmed,
321 Double_t *upar, Int_t np) {
323 return fGeom->Gsvolu(name, shape, nmed, upar, np);
326 //______________________________________________________________________________
327 void TFluka::Gsdvn(const char *name, const char *mother, Int_t ndiv,
330 fGeom->Gsdvn(name, mother, ndiv, iaxis);
333 //______________________________________________________________________________
334 void TFluka::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
335 Int_t iaxis, Double_t c0i, Int_t numed) {
337 fGeom->Gsdvn2(name, mother, ndiv, iaxis, c0i, numed);
340 //______________________________________________________________________________
341 void TFluka::Gsdvt(const char *name, const char *mother, Double_t step,
342 Int_t iaxis, Int_t numed, Int_t ndvmx) {
344 fGeom->Gsdvt(name, mother, step, iaxis, numed, ndvmx);
347 //______________________________________________________________________________
348 void TFluka::Gsdvt2(const char *name, const char *mother, Double_t step,
349 Int_t iaxis, Double_t c0, Int_t numed, Int_t ndvmx) {
351 fGeom->Gsdvt2(name, mother, step, iaxis, c0, numed, ndvmx);
354 //______________________________________________________________________________
355 void TFluka::Gsord(const char * /*name*/, Int_t /*iax*/) {
357 // Nothing to do with TGeo
360 //______________________________________________________________________________
361 void TFluka::Gspos(const char *name, Int_t nr, const char *mother,
362 Double_t x, Double_t y, Double_t z, Int_t irot,
365 fGeom->Gspos(name, nr, mother, x, y, z, irot, konly);
368 //______________________________________________________________________________
369 void TFluka::Gsposp(const char *name, Int_t nr, const char *mother,
370 Double_t x, Double_t y, Double_t z, Int_t irot,
371 const char *konly, Float_t *upar, Int_t np) {
373 fGeom->Gsposp(name, nr, mother, x, y, z, irot, konly, upar, np);
376 //______________________________________________________________________________
377 void TFluka::Gsposp(const char *name, Int_t nr, const char *mother,
378 Double_t x, Double_t y, Double_t z, Int_t irot,
379 const char *konly, Double_t *upar, Int_t np) {
381 fGeom->Gsposp(name, nr, mother, x, y, z, irot, konly, upar, np);
384 //______________________________________________________________________________
385 void TFluka::Gsbool(const char* /*onlyVolName*/, const char* /*manyVolName*/) {
387 // Nothing to do with TGeo
388 Warning("Gsbool", "Not implemented with TGeo");
391 //______________________________________________________________________________
392 void TFluka::SetCerenkov(Int_t /*itmed*/, Int_t /*npckov*/, Float_t */*ppckov*/,
393 Float_t * /*absco*/, Float_t * /*effic*/, Float_t * /*rindex*/) {
395 // Not implemented with TGeo - what G4 did ? Any FLUKA card generated?
396 Warning("SetCerenkov", "Not implemented with TGeo");
399 //______________________________________________________________________________
400 void TFluka::SetCerenkov(Int_t /*itmed*/, Int_t /*npckov*/, Double_t * /*ppckov*/,
401 Double_t * /*absco*/, Double_t * /*effic*/, Double_t * /*rindex*/) {
403 // Not implemented with TGeo - what G4 did ? Any FLUKA card generated?
404 Warning("SetCerenkov", "Not implemented with TGeo");
408 //______________________________________________________________________________
409 void TFluka::WriteEuclid(const char* /*fileName*/, const char* /*topVol*/,
410 Int_t /*number*/, Int_t /*nlevel*/) {
413 Warning("WriteEuclid", "Not implemented with TGeo");
418 //_____________________________________________________________________________
419 // methods needed by the stepping
420 //____________________________________________________________________________
422 Int_t TFluka::GetMedium() const {
424 // Get the medium number for the current fluka region
426 return fGeom->GetMedium(); // this I need to check due to remapping !!!
431 //____________________________________________________________________________
432 // particle table usage
433 // ID <--> PDG transformations
434 //_____________________________________________________________________________
435 Int_t TFluka::IdFromPDG(Int_t pdg) const
438 // Return Fluka code from PDG and pseudo ENDF code
440 // Catch the feedback photons
441 if (pdg == 50000051) return (-1);
442 // MCIHAD() goes from pdg to fluka internal.
443 Int_t intfluka = mcihad(pdg);
444 // KPTOIP array goes from internal to official
445 return GetFlukaKPTOIP(intfluka);
448 //______________________________________________________________________________
449 Int_t TFluka::PDGFromId(Int_t id) const
452 // Return PDG code and pseudo ENDF code from Fluka code
454 // IPTOKP array goes from official to internal
458 if (fVerbosityLevel >= 1)
459 printf("\n PDGFromId: Cerenkov Photon \n");
464 if (fVerbosityLevel >= 1)
465 printf("PDGFromId: Error id = 0\n");
469 Int_t intfluka = GetFlukaIPTOKP(id);
471 if (fVerbosityLevel >= 1)
472 printf("PDGFromId: Error intfluka = 0: %d\n", id);
474 } else if (intfluka < 0) {
475 if (fVerbosityLevel >= 1)
476 printf("PDGFromId: Error intfluka < 0: %d\n", id);
479 if (fVerbosityLevel >= 3)
480 printf("mpdgha called with %d %d \n", id, intfluka);
481 // MPDGHA() goes from fluka internal to pdg.
482 return mpdgha(intfluka);
485 //_____________________________________________________________________________
486 // methods for physics management
487 //____________________________________________________________________________
492 //______________________________________________________________________________
493 void TFluka::SetProcess(const char* flagName, Int_t flagValue)
496 if (iNbOfProc < 100) {
497 for (i=0; i<iNbOfProc; i++) {
498 if (strcmp(&sProcessFlag[i][0],flagName) == 0) {
499 iProcessValue[iNbOfProc] = flagValue;
503 strcpy(&sProcessFlag[iNbOfProc][0],flagName);
504 iProcessValue[iNbOfProc++] = flagValue;
507 cout << "Nb of SetProcess calls exceeds 100 - ignored" << endl;
509 iNbOfProc = iNbOfProc;
512 //______________________________________________________________________________
513 void TFluka::SetCut(const char* cutName, Double_t cutValue)
516 if (iNbOfCut < 100) {
517 for (i=0; i<iNbOfCut; i++) {
518 if (strcmp(&sCutFlag[i][0],cutName) == 0) {
519 fCutValue[iNbOfCut] = cutValue;
523 strcpy(&sCutFlag[iNbOfCut][0],cutName);
524 fCutValue[iNbOfCut++] = cutValue;
527 cout << "Nb of SetCut calls exceeds 100 - ignored" << endl;
532 //______________________________________________________________________________
533 Double_t TFluka::Xsec(char*, Double_t, Int_t, Int_t)
535 printf("WARNING: Xsec not yet implemented !\n"); return -1.;
539 //______________________________________________________________________________
540 void TFluka::InitPhysics()
542 // Last material number taken from the "corealice.inp" file, presently 31
543 // !!! it should be available from Flugg !!!
546 Float_t fLastMaterial = fGeom->GetLastMaterialIndex();
547 printf(" last FLUKA material is %g\n", fLastMaterial);
548 // construct file names
549 TString sAliceInp = getenv("ALICE_ROOT");
550 sAliceInp +="/TFluka/input/";
551 TString sAliceTmp ="./flukaMat.inp";
552 TString sAliceCoreInp = sAliceInp;
553 sAliceInp += GetInputFileName();
554 sAliceCoreInp += GetCoreInputFileName();
555 ifstream AliceCoreInp(sAliceCoreInp.Data());
556 ifstream AliceFlukaMat(sAliceTmp.Data());
557 ofstream AliceInp(sAliceInp.Data());
559 // copy core input file
561 Float_t fEventsPerRun;
563 while (AliceCoreInp.getline(sLine,255)) {
564 if (strncmp(sLine,"GEOEND",6) != 0)
565 AliceInp << sLine << endl; // copy until GEOEND card
567 AliceInp << "GEOEND" << endl; // add GEOEND card
570 } // end of while until GEOEND card
573 while (AliceFlukaMat.getline(sLine,255)) { // copy flukaMat.inp file
574 AliceInp << sLine << endl;
577 while (AliceCoreInp.getline(sLine,255)) {
578 if (strncmp(sLine,"START",5) != 0)
579 AliceInp << sLine << endl;
581 sscanf(sLine+10,"%10f",&fEventsPerRun);
584 } //end of while until START card
587 // in G3 the process control values meaning can be different for
588 // different processes, but for most of them is:
589 // 0 process is not activated
590 // 1 process is activated WITH generation of secondaries
591 // 2 process is activated WITHOUT generation of secondaries
592 // if process does not generate secondaries => 1 same as 2
601 // Loop over number of SetProcess calls
602 AliceInp << "*----------------------------------------------------------------------------- ";
604 AliceInp << "*----- The following data are generated from SetProcess and SetCut calls ----- ";
606 AliceInp << "*----------------------------------------------------------------------------- ";
608 for (i=0; i<iNbOfProc; i++) {
611 // G3 default value: 1
612 // G4 processes: G4eplusAnnihilation/G4IeplusAnnihilation
615 // flag = 0 no annihilation
616 // flag = 1 annihilation, decays processed
617 // flag = 2 annihilation, no decay product stored
618 // gMC ->SetProcess("ANNI",1); // EMFCUT -1. 0. 0. 3. lastmat 0. ANNH-THR
619 if (strncmp(&sProcessFlag[i][0],"ANNI",4) == 0) {
620 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
623 AliceInp << "*Kinetic energy threshold (GeV) for e+ annihilation - resets to default=0.";
625 AliceInp << "*Generated from call: SetProcess('ANNI',1) or SetProcess('ANNI',2)";
627 AliceInp << setw(10) << "EMFCUT ";
628 AliceInp << setiosflags(ios::scientific) << setprecision(5);
629 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
630 AliceInp << setw(10) << -1.0; // kinetic energy threshold (GeV) for e+ annihilation (resets to default=0)
631 AliceInp << setw(10) << 0.0; // not used
632 AliceInp << setw(10) << 0.0; // not used
633 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
634 AliceInp << setw(10) << setprecision(2);
635 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
636 AliceInp << setprecision(1);
637 AliceInp << setw(10) << 1.0; // step length in assigning indices
638 AliceInp << setw(8) << "ANNH-THR";
641 else if (iProcessValue[i] == 0) {
644 AliceInp << "*No annihilation - no FLUKA card generated";
646 AliceInp << "*Generated from call: SetProcess('ANNI',0)";
652 AliceInp << "*Illegal flag value in SetProcess('ANNI',?) call.";
654 AliceInp << "*No FLUKA card generated";
659 // bremsstrahlung and pair production are both activated
660 // G3 default value: 1
661 // G4 processes: G4eBremsstrahlung/G4IeBremsstrahlung,
662 // G4MuBremsstrahlung/G4IMuBremsstrahlung,
663 // G4LowEnergyBremstrahlung
664 // Particles: e-/e+; mu+/mu-
666 // flag = 0 no bremsstrahlung
667 // flag = 1 bremsstrahlung, photon processed
668 // flag = 2 bremsstrahlung, no photon stored
669 // gMC ->SetProcess("BREM",1); // PAIRBREM 2. 0. 0. 3. lastmat
670 // EMFCUT -1. 0. 0. 3. lastmat 0. ELPO-THR
671 // G3 default value: 1
672 // G4 processes: G4GammaConversion,
673 // G4MuPairProduction/G4IMuPairProduction
674 // G4LowEnergyGammaConversion
675 // Particles: gamma, mu
677 // flag = 0 no delta rays
678 // flag = 1 delta rays, secondaries processed
679 // flag = 2 delta rays, no secondaries stored
680 // gMC ->SetProcess("PAIR",1); // PAIRBREM 1. 0. 0. 3. lastmat
681 // EMFCUT 0. 0. -1. 3. lastmat 0. PHOT-THR
682 else if ((strncmp(&sProcessFlag[i][0],"PAIR",4) == 0) && (iProcessValue[i] == 1 || iProcessValue[i] == 2)) {
683 for (j=0; j<iNbOfProc; j++) {
684 if ((strncmp(&sProcessFlag[j][0],"BREM",4) == 0) && (iProcessValue[j] == 1 || iProcessValue[j] == 2)) {
687 AliceInp << "*Bremsstrahlung and pair production by muons and charged hadrons both activated";
689 AliceInp << "*Generated from call: SetProcess('BREM',1) and SetProcess('PAIR',1)";
691 AliceInp << "*Energy threshold set by call SetCut('BCUTM',cut) or set to 0.";
693 AliceInp << "*Energy threshold set by call SetCut('PPCUTM',cut) or set to 0.";
695 AliceInp << setw(10) << "PAIRBREM ";
696 AliceInp << setiosflags(ios::scientific) << setprecision(5);
697 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
698 AliceInp << setw(10) << 3.0; // bremsstrahlung and pair production by muons and charged hadrons both are activated
699 // direct pair production by muons
700 // G4 particles: "e-", "e+"
701 // G3 default value: 0.01 GeV
702 //gMC ->SetCut("PPCUTM",cut); // total energy cut for direct pair prod. by muons
704 for (k=0; k<iNbOfCut; k++) {
705 if (strncmp(&sCutFlag[k][0],"PPCUTM",6) == 0) fCut = fCutValue[k];
707 AliceInp << setiosflags(ios::scientific) << setprecision(5);
708 AliceInp << setw(10) << fCut; // e+, e- kinetic energy threshold (in GeV) for explicit pair production.
709 // muon and hadron bremsstrahlung
710 // G4 particles: "gamma"
711 // G3 default value: CUTGAM=0.001 GeV
712 //gMC ->SetCut("BCUTM",cut); // cut for muon and hadron bremsstrahlung
714 for (k=0; k<iNbOfCut; k++) {
715 if (strncmp(&sCutFlag[k][0],"BCUTM",5) == 0) fCut = fCutValue[k];
717 AliceInp << setiosflags(ios::scientific) << setprecision(5);
718 AliceInp << setw(10) << fCut; // photon energy threshold (GeV) for explicit bremsstrahlung production
719 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
720 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
721 AliceInp << setw(10) << setprecision(2);
722 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
728 AliceInp << "*Kinetic energy threshold (GeV) for e+/e- bremsstrahlung - resets to default=0.";
730 AliceInp << "*Generated from call: SetProcess('BREM',1);";
732 AliceInp << setw(10) << "EMFCUT ";
734 for (k=0; k<iNbOfCut; k++) {
735 if (strncmp(&sCutFlag[k][0],"BCUTE",5) == 0) fCut = fCutValue[k];
737 AliceInp << setiosflags(ios::scientific) << setprecision(5);
738 AliceInp << setw(10) << fCut; // kinetic energy threshold (GeV) for e+/e- bremsstrahlung (resets to default=0)
739 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
740 AliceInp << setw(10) << 0.0; // not used
741 AliceInp << setw(10) << 0.0; // not used
742 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
743 AliceInp << setw(10) << setprecision(2);
744 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
745 AliceInp << setprecision(1);
746 AliceInp << setw(10) << 1.0; // step length in assigning indices
747 AliceInp << setw(8) << "ELPO-THR";
753 AliceInp << "*Pair production by electrons is activated";
755 AliceInp << "*Generated from call: SetProcess('PAIR',1);";
757 AliceInp << setw(10) << "EMFCUT ";
758 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
759 AliceInp << setw(10) << 0.0; // energy threshold (GeV) for Compton scattering (= 0.0 : ignored)
760 AliceInp << setw(10) << 0.0; // energy threshold (GeV) for Photoelectric (= 0.0 : ignored)
762 for (j=0; j<iNbOfCut; j++) {
763 if (strncmp(&sCutFlag[j][0],"CUTGAM",6) == 0) fCut = fCutValue[j];
765 AliceInp << setiosflags(ios::scientific) << setprecision(5);
766 AliceInp << setw(10) << fCut; // energy threshold (GeV) for gamma pair production (< 0.0 : resets to default, = 0.0 : ignored)
767 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
768 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
769 AliceInp << setprecision(2);
770 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
771 AliceInp << setprecision(1);
772 AliceInp << setw(10) << 1.0; // step length in assigning indices
773 AliceInp << setw(8) << "PHOT-THR";
776 } // end of if for BREM
777 } // end of loop for BREM
779 // only pair production by muons and charged hadrons is activated
782 AliceInp << "*Pair production by muons and charged hadrons is activated";
784 AliceInp << "*Generated from call: SetProcess('PAIR',1) or SetProcess('PAIR',2)";
786 AliceInp << "*Energy threshold set by call SetCut('PPCUTM',cut) or set to 0.";
788 AliceInp << setw(10) << "PAIRBREM ";
789 AliceInp << setiosflags(ios::scientific) << setprecision(5);
790 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
791 AliceInp << setw(10) << 1.0; // pair production by muons and charged hadrons is activated
792 // direct pair production by muons
793 // G4 particles: "e-", "e+"
794 // G3 default value: 0.01 GeV
795 //gMC ->SetCut("PPCUTM",cut); // total energy cut for direct pair prod. by muons
796 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
797 AliceInp << setw(10) << 0.0; // e+, e- kinetic energy threshold (in GeV) for explicit pair production.
798 AliceInp << setw(10) << 0.0; // no explicit bremsstrahlung production is simulated
799 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
800 AliceInp << setprecision(2);
801 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
807 AliceInp << "*Pair production by electrons is activated";
809 AliceInp << "*Generated from call: SetProcess('PAIR',1) or SetProcess('PAIR',2)";
811 AliceInp << setw(10) << "EMFCUT ";
812 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
813 AliceInp << setw(10) << 0.0; // energy threshold (GeV) for Compton scattering (= 0.0 : ignored)
814 AliceInp << setw(10) << 0.0; // energy threshold (GeV) for Photoelectric (= 0.0 : ignored)
817 for (j=0; j<iNbOfCut; j++) {
818 if (strncmp(&sCutFlag[j][0],"CUTGAM",6) == 0) fCut = fCutValue[j];
820 AliceInp << setiosflags(ios::scientific) << setprecision(5);
821 AliceInp << setw(10) << fCut; // energy threshold (GeV) for gamma pair production (< 0.0 : resets to default, = 0.0 : ignored)
822 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
823 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
824 AliceInp << setprecision(2);
825 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
826 AliceInp << setprecision(1);
827 AliceInp << setw(10) << 1.0; // step length in assigning indices
828 AliceInp << setw(8) << "PHOT-THR";
833 } // end of if for PAIR
838 // G3 default value: 1
839 // G4 processes: G4eBremsstrahlung/G4IeBremsstrahlung,
840 // G4MuBremsstrahlung/G4IMuBremsstrahlung,
841 // G4LowEnergyBremstrahlung
842 // Particles: e-/e+; mu+/mu-
844 // flag = 0 no bremsstrahlung
845 // flag = 1 bremsstrahlung, photon processed
846 // flag = 2 bremsstrahlung, no photon stored
847 // gMC ->SetProcess("BREM",1); // PAIRBREM 2. 0. 0. 3. lastmat
848 // EMFCUT -1. 0. 0. 3. lastmat 0. ELPO-THR
849 else if (strncmp(&sProcessFlag[i][0],"BREM",4) == 0) {
850 for (j=0; j<iNbOfProc; j++) {
851 if ((strncmp(&sProcessFlag[j][0],"PAIR",4) == 0) && iProcessValue[j] == 1) goto NOBREM;
853 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
856 AliceInp << "*Bremsstrahlung by muons and charged hadrons is activated";
858 AliceInp << "*Generated from call: SetProcess('BREM',1) or SetProcess('BREM',2)";
860 AliceInp << "*Energy threshold set by call SetCut('BCUTM',cut) or set to 0.";
862 AliceInp << setw(10) << "PAIRBREM ";
863 AliceInp << setiosflags(ios::scientific) << setprecision(5);
864 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
865 AliceInp << setw(10) << 2.0; // bremsstrahlung by muons and charged hadrons is activated
866 AliceInp << setw(10) << 0.0; // no meaning
867 // muon and hadron bremsstrahlung
868 // G4 particles: "gamma"
869 // G3 default value: CUTGAM=0.001 GeV
870 //gMC ->SetCut("BCUTM",cut); // cut for muon and hadron bremsstrahlung
872 for (j=0; j<iNbOfCut; j++) {
873 if (strncmp(&sCutFlag[j][0],"BCUTM",5) == 0) fCut = fCutValue[j];
875 AliceInp << setw(10) << fCut; // photon energy threshold (GeV) for explicit bremsstrahlung production
876 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
877 AliceInp << setw(10) << setprecision(2);
878 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
884 AliceInp << "*Kinetic energy threshold (GeV) for e+/e- bremsstrahlung - resets to default=0.";
886 AliceInp << "*Generated from call: SetProcess('BREM',1);";
888 AliceInp << setw(10) << "EMFCUT ";
889 AliceInp << setiosflags(ios::scientific) << setprecision(5);
890 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
891 AliceInp << setw(10) << -1.0; // kinetic energy threshold (GeV) for e+/e- bremsstrahlung (resets to default=0)
892 AliceInp << setw(10) << 0.0; // not used
893 AliceInp << setw(10) << 0.0; // not used
894 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
895 AliceInp << setw(10) << setprecision(2);
896 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
897 AliceInp << setprecision(1);
898 AliceInp << setw(10) << 1.0; // step length in assigning indices
899 AliceInp << setw(8) << "ELPO-THR";
902 else if (iProcessValue[i] == 0) {
905 AliceInp << "*No bremsstrahlung - no FLUKA card generated";
907 AliceInp << "*Generated from call: SetProcess('BREM',0)";
913 AliceInp << "*Illegal flag value in SetProcess('BREM',?) call.";
915 AliceInp << "*No FLUKA card generated";
920 } // end of else if (strncmp(&sProcessFlag[i][0],"BREM",4) == 0)
923 // Cerenkov photon generation
924 // G3 default value: 0
925 // G4 process: G4Cerenkov
927 // Particles: charged
929 // flag = 0 no Cerenkov photon generation
930 // flag = 1 Cerenkov photon generation
931 // flag = 2 Cerenkov photon generation with primary stopped at each step
932 //xx gMC ->SetProcess("CKOV",1); // ??? Cerenkov photon generation
933 else if (strncmp(&sProcessFlag[i][0],"CKOV",4) == 0) {
934 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
937 AliceInp << "*Cerenkov photon generation";
939 AliceInp << "*Generated from call: SetProcess('CKOV',1) or SetProcess('CKOV',2)";
941 AliceInp << setw(10) << "OPT-PROD ";
942 AliceInp << setiosflags(ios::scientific) << setprecision(5);
943 AliceInp << setw(10) << 2.07e-9 ; // minimum Cerenkov photon emission energy (in GeV!). Default: 2.07E-9 GeV (corresponding to 600 nm)
944 AliceInp << setw(10) << 4.96e-9; // maximum Cerenkov photon emission energy (in GeV!). Default: 4.96E-9 GeV (corresponding to 250 nm)
945 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
946 AliceInp << setw(10) << 0.0; // not used
947 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
948 AliceInp << setprecision(2);
949 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
950 AliceInp << setprecision(1);
951 AliceInp << setw(10) << 1.0; // step length in assigning indices
952 AliceInp << setw(8) << "CERENKOV";
955 else if (iProcessValue[i] == 0) {
958 AliceInp << "*No Cerenkov photon generation";
960 AliceInp << "*Generated from call: SetProcess('CKOV',0)";
962 AliceInp << setw(10) << "OPT-PROD ";
963 AliceInp << setiosflags(ios::scientific) << setprecision(5);
964 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
965 AliceInp << setw(10) << 0.0; // not used
966 AliceInp << setw(10) << 0.0; // not used
967 AliceInp << setw(10) << 0.0; // not used
968 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
969 AliceInp << setprecision(2);
970 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
971 AliceInp << setprecision(1);
972 AliceInp << setw(10) << 1.0; // step length in assigning indices
973 AliceInp << setw(8) << "CERE-OFF";
979 AliceInp << "*Illegal flag value in SetProcess('CKOV',?) call.";
981 AliceInp << "*No FLUKA card generated";
984 } // end of else if (strncmp(&sProcessFlag[i][0],"CKOV",4) == 0)
987 // Compton scattering
988 // G3 default value: 1
989 // G4 processes: G4ComptonScattering,
990 // G4LowEnergyCompton,
991 // G4PolarizedComptonScattering
994 // flag = 0 no Compton scattering
995 // flag = 1 Compton scattering, electron processed
996 // flag = 2 Compton scattering, no electron stored
997 // gMC ->SetProcess("COMP",1); // EMFCUT -1. 0. 0. 3. lastmat 0. PHOT-THR
998 else if (strncmp(&sProcessFlag[i][0],"COMP",4) == 0) {
999 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
1002 AliceInp << "*Energy threshold (GeV) for Compton scattering - resets to default=0.";
1004 AliceInp << "*Generated from call: SetProcess('COMP',1);";
1006 AliceInp << setw(10) << "EMFCUT ";
1007 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1008 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1009 AliceInp << setw(10) << -1.0; // energy threshold (GeV) for Compton scattering - resets to default=0.
1010 AliceInp << setw(10) << 0.0; // not used
1011 AliceInp << setw(10) << 0.0; // not used
1012 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1013 AliceInp << setprecision(2);
1014 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1015 AliceInp << setprecision(1);
1016 AliceInp << setw(10) << 1.0; // step length in assigning indices
1017 AliceInp << setw(8) << "PHOT-THR";
1020 else if (iProcessValue[i] == 0) {
1023 AliceInp << "*No Compton scattering - no FLUKA card generated";
1025 AliceInp << "*Generated from call: SetProcess('COMP',0)";
1031 AliceInp << "*Illegal flag value in SetProcess('COMP',?) call.";
1033 AliceInp << "*No FLUKA card generated";
1036 } // end of else if (strncmp(&sProcessFlag[i][0],"COMP",4) == 0)
1039 // G3 default value: 1
1040 // G4 process: G4Decay
1042 // Particles: all which decay is applicable for
1044 // flag = 0 no decays
1045 // flag = 1 decays, secondaries processed
1046 // flag = 2 decays, no secondaries stored
1047 //gMC ->SetProcess("DCAY",1); // not available
1048 else if ((strncmp(&sProcessFlag[i][0],"DCAY",4) == 0) && iProcessValue[i] == 1)
1049 cout << "SetProcess for flag=" << &sProcessFlag[i][0] << " value=" << iProcessValue[i] << " not avaliable!" << endl;
1052 // G3 default value: 2
1053 // !! G4 treats delta rays in different way
1054 // G4 processes: G4eIonisation/G4IeIonization,
1055 // G4MuIonisation/G4IMuIonization,
1056 // G4hIonisation/G4IhIonisation
1057 // Particles: charged
1059 // flag = 0 no energy loss
1060 // flag = 1 restricted energy loss fluctuations
1061 // flag = 2 complete energy loss fluctuations
1062 // flag = 3 same as 1
1063 // flag = 4 no energy loss fluctuations
1064 // gMC ->SetProcess("DRAY",0); // DELTARAY 1.E+6 0. 0. 3. lastmat 0.
1065 else if (strncmp(&sProcessFlag[i][0],"DRAY",4) == 0) {
1066 if (iProcessValue[i] == 0 || iProcessValue[i] == 4) {
1069 AliceInp << "*Kinetic energy threshold (GeV) for delta ray production";
1071 AliceInp << "*Generated from call: SetProcess('DRAY',0) or SetProcess('DRAY',4)";
1073 AliceInp << "*No delta ray production by muons - threshold set artificially high";
1075 AliceInp << setw(10) << "DELTARAY ";
1076 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1077 AliceInp << setw(10) << 1.0e+6; // kinetic energy threshold (GeV) for delta ray production (discrete energy transfer)
1078 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1079 AliceInp << setw(10) << 0.0; // ignored
1080 AliceInp << setw(10) << 0.0; // ignored
1081 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1082 AliceInp << setw(10) << setprecision(2);
1083 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1084 AliceInp << setprecision(1);
1085 AliceInp << setw(10) << 1.0; // step length in assigning indices
1088 else if (iProcessValue[i] == 1 || iProcessValue[i] == 2 || iProcessValue[i] == 3) {
1091 AliceInp << "*Kinetic energy threshold (GeV) for delta ray production";
1093 AliceInp << "*Generated from call: SetProcess('DRAY',flag), flag=1,2,3";
1095 AliceInp << "*Delta ray production by muons switched on";
1097 AliceInp << "*Energy threshold set by call SetCut('DCUTM',cut) or set to 0.";
1099 AliceInp << setw(10) << "DELTARAY ";
1100 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1102 for (j=0; j<iNbOfCut; j++) {
1103 if (strncmp(&sCutFlag[j][0],"DCUTM",5) == 0) fCut = fCutValue[j];
1105 AliceInp << setw(10) << fCut; // kinetic energy threshold (GeV) for delta ray production (discrete energy transfer)
1106 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1107 AliceInp << setw(10) << 0.0; // ignored
1108 AliceInp << setw(10) << 0.0; // ignored
1109 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1110 AliceInp << setw(10) << setprecision(2);
1111 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1112 AliceInp << setprecision(1);
1113 AliceInp << setw(10) << 1.0; // step length in assigning indices
1119 AliceInp << "*Illegal flag value in SetProcess('DRAY',?) call.";
1121 AliceInp << "*No FLUKA card generated";
1124 } // end of else if (strncmp(&sProcessFlag[i][0],"DRAY",4) == 0)
1127 // G3 default value: 1
1128 // G4 processes: all defined by TG4PhysicsConstructorHadron
1130 // Particles: hadrons
1132 // flag = 0 no multiple scattering
1133 // flag = 1 hadronic interactions, secondaries processed
1134 // flag = 2 hadronic interactions, no secondaries stored
1135 // gMC ->SetProcess("HADR",1); // ??? hadronic process
1136 //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3) ?????
1137 else if (strncmp(&sProcessFlag[i][0],"HADR",4) == 0) {
1138 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
1141 AliceInp << "*Hadronic interaction is ON by default in FLUKA";
1143 AliceInp << "*No FLUKA card generated";
1146 else if (iProcessValue[i] == 0) {
1149 AliceInp << "*Hadronic interaction is set OFF";
1151 AliceInp << "*Generated from call: SetProcess('HADR',0);";
1153 AliceInp << setw(10) << "MULSOPT ";
1154 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1155 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1156 AliceInp << setw(10) << 0.0; // ignored
1157 AliceInp << setw(10) << 3.0; // multiple scattering for hadrons and muons is completely suppressed
1158 AliceInp << setw(10) << 0.0; // no spin-relativistic corrections
1159 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1160 AliceInp << setprecision(2);
1161 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1168 AliceInp << "*Illegal flag value in SetProcess('HADR',?) call.";
1170 AliceInp << "*No FLUKA card generated";
1173 } // end of else if (strncmp(&sProcessFlag[i][0],"HADR",4) == 0)
1177 // G3 default value: 2
1178 // G4 processes: G4eIonisation/G4IeIonization,
1179 // G4MuIonisation/G4IMuIonization,
1180 // G4hIonisation/G4IhIonisation
1182 // Particles: charged
1184 // flag=0 no energy loss
1185 // flag=1 restricted energy loss fluctuations
1186 // flag=2 complete energy loss fluctuations
1188 // flag=4 no energy loss fluctuations
1189 // If the value ILOSS is changed, then (in G3) cross-sections and energy
1190 // loss tables must be recomputed via the command 'PHYSI'
1191 // gMC ->SetProcess("LOSS",2); // ??? IONFLUCT ? energy loss
1192 else if (strncmp(&sProcessFlag[i][0],"LOSS",4) == 0) {
1193 if (iProcessValue[i] == 2) { // complete energy loss fluctuations
1196 AliceInp << "*Complete energy loss fluctuations do not exist in FLUKA";
1198 AliceInp << "*Generated from call: SetProcess('LOSS',2);";
1200 AliceInp << "*flag=2=complete energy loss fluctuations";
1202 AliceInp << "*No input card generated";
1205 else if (iProcessValue[i] == 1 || iProcessValue[i] == 3) { // restricted energy loss fluctuations
1208 AliceInp << "*Restricted energy loss fluctuations";
1210 AliceInp << "*Generated from call: SetProcess('LOSS',1) or SetProcess('LOSS',3)";
1212 AliceInp << setw(10) << "IONFLUCT ";
1213 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1214 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1215 AliceInp << setw(10) << 1.0; // restricted energy loss fluctuations (for hadrons and muons) switched on
1216 AliceInp << setw(10) << 1.0; // restricted energy loss fluctuations (for e+ and e-) switched on
1217 AliceInp << setw(10) << 1.0; // minimal accuracy
1218 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1219 AliceInp << setprecision(2);
1220 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1223 else if (iProcessValue[i] == 4) { // no energy loss fluctuations
1226 AliceInp << "*No energy loss fluctuations";
1228 AliceInp << "*Generated from call: SetProcess('LOSS',4)";
1230 AliceInp << setw(10) << -1.0; // restricted energy loss fluctuations (for hadrons and muons) switched off
1231 AliceInp << setw(10) << -1.0; // restricted energy loss fluctuations (for e+ and e-) switched off
1232 AliceInp << setw(10) << 1.0; // minimal accuracy
1233 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1234 AliceInp << setprecision(2);
1235 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1241 AliceInp << "*Illegal flag value in SetProcess('LOSS',?) call.";
1243 AliceInp << "*No FLUKA card generated";
1246 } // end of else if (strncmp(&sProcessFlag[i][0],"LOSS",4) == 0)
1249 // multiple scattering
1250 // G3 default value: 1
1251 // G4 process: G4MultipleScattering/G4IMultipleScattering
1253 // Particles: charged
1255 // flag = 0 no multiple scattering
1256 // flag = 1 Moliere or Coulomb scattering
1257 // flag = 2 Moliere or Coulomb scattering
1258 // flag = 3 Gaussian scattering
1259 // gMC ->SetProcess("MULS",1); // MULSOPT multiple scattering
1260 else if (strncmp(&sProcessFlag[i][0],"MULS",4) == 0) {
1261 if (iProcessValue[i] == 1 || iProcessValue[i] == 2 || iProcessValue[i] == 3) {
1264 AliceInp << "*Multiple scattering is ON by default for e+e- and for hadrons/muons";
1266 AliceInp << "*No FLUKA card generated";
1269 else if (iProcessValue[i] == 0) {
1272 AliceInp << "*Multiple scattering is set OFF";
1274 AliceInp << "*Generated from call: SetProcess('MULS',0);";
1276 AliceInp << setw(10) << "MULSOPT ";
1277 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1278 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1279 AliceInp << setw(10) << 0.0; // ignored
1280 AliceInp << setw(10) << 3.0; // multiple scattering for hadrons and muons is completely suppressed
1281 AliceInp << setw(10) << 3.0; // multiple scattering for e+ and e- is completely suppressed
1282 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1283 AliceInp << setprecision(2);
1284 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1290 AliceInp << "*Illegal flag value in SetProcess('MULS',?) call.";
1292 AliceInp << "*No FLUKA card generated";
1295 } // end of else if (strncmp(&sProcessFlag[i][0],"MULS",4) == 0)
1298 // muon nuclear interaction
1299 // G3 default value: 0
1300 // G4 processes: G4MuNuclearInteraction,
1301 // G4MuonMinusCaptureAtRest
1305 // flag = 0 no muon-nuclear interaction
1306 // flag = 1 nuclear interaction, secondaries processed
1307 // flag = 2 nuclear interaction, secondaries not processed
1308 // gMC ->SetProcess("MUNU",1); // MUPHOTON 1. 0. 0. 3. lastmat
1309 else if (strncmp(&sProcessFlag[i][0],"MUNU",4) == 0) {
1310 if (iProcessValue[i] == 1) {
1313 AliceInp << "*Muon nuclear interactions with production of secondary hadrons";
1315 AliceInp << "*Generated from call: SetProcess('MUNU',1);";
1317 AliceInp << setw(10) << "MUPHOTON ";
1318 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1319 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1320 AliceInp << setw(10) << 1.0; // full simulation of muon nuclear interactions and production of secondary hadrons
1321 AliceInp << setw(10) << 0.0; // ratio of longitudinal to transverse virtual photon cross-section - Default = 0.25.
1322 AliceInp << setw(10) << 0.0; // fraction of rho-like interactions ( must be < 1) - Default = 0.75.
1323 AliceInp << setprecision(1);
1324 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1325 AliceInp << setprecision(2);
1326 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1329 else if (iProcessValue[i] == 2) {
1332 AliceInp << "*Muon nuclear interactions without production of secondary hadrons";
1334 AliceInp << "*Generated from call: SetProcess('MUNU',2);";
1336 AliceInp << setw(10) << "MUPHOTON ";
1337 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1338 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1339 AliceInp << setw(10) << 2.0; // full simulation of muon nuclear interactions and production of secondary hadrons
1340 AliceInp << setw(10) << 0.0; // ratio of longitudinal to transverse virtual photon cross-section - Default = 0.25.
1341 AliceInp << setw(10) << 0.0; // fraction of rho-like interactions ( must be < 1) - Default = 0.75.
1342 AliceInp << setprecision(1);
1343 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1344 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1347 else if (iProcessValue[i] == 0) {
1350 AliceInp << "*No muon nuclear interaction - no FLUKA card generated";
1352 AliceInp << "*Generated from call: SetProcess('MUNU',0)";
1358 AliceInp << "*Illegal flag value in SetProcess('MUNU',?) call.";
1360 AliceInp << "*No FLUKA card generated";
1363 } // end of else if (strncmp(&sProcessFlag[i][0],"MUNU",4) == 0)
1367 // G3 default value: 0
1372 // gMC ->SetProcess("PFIS",0); // PHOTONUC -1. 0. 0. 3. lastmat 0.
1373 // flag = 0 no photon fission
1374 // flag = 1 photon fission, secondaries processed
1375 // flag = 2 photon fission, no secondaries stored
1376 else if (strncmp(&sProcessFlag[i][0],"PFIS",4) == 0) {
1377 if (iProcessValue[i] == 0) {
1380 AliceInp << "*No photonuclear interactions";
1382 AliceInp << "*Generated from call: SetProcess('PFIS',0);";
1384 AliceInp << setw(10) << "PHOTONUC ";
1385 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1386 AliceInp << setw(10) << -1.0; // no photonuclear interactions
1387 AliceInp << setw(10) << 0.0; // not used
1388 AliceInp << setw(10) << 0.0; // not used
1389 AliceInp << setw(10) << 3.0; // upper bound of the material indices in which the respective thresholds apply
1390 AliceInp << setprecision(2);
1391 AliceInp << setw(10) << fLastMaterial;
1392 AliceInp << setprecision(1); // upper bound of the material indices in which the respective thresholds apply
1393 AliceInp << setprecision(1);
1394 AliceInp << setw(10) << 1.0; // step length in assigning indices
1397 else if (iProcessValue[i] == 1) {
1400 AliceInp << "*Photon nuclear interactions are activated at all energies";
1402 AliceInp << "*Generated from call: SetProcess('PFIS',1);";
1404 AliceInp << setw(10) << "PHOTONUC ";
1405 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1406 AliceInp << setw(10) << 1.0; // photonuclear interactions are activated at all energies
1407 AliceInp << setw(10) << 0.0; // not used
1408 AliceInp << setw(10) << 0.0; // not used
1409 AliceInp << setprecision(2);
1410 AliceInp << setw(10) << 3.0; // upper bound of the material indices in which the respective thresholds apply
1411 AliceInp << setw(10) << fLastMaterial;
1412 AliceInp << setprecision(1); // upper bound of the material indices in which the respective thresholds apply
1413 AliceInp << setprecision(1);
1414 AliceInp << setw(10) << 1.0; // step length in assigning indices
1417 else if (iProcessValue[i] == 0) {
1420 AliceInp << "*No photofission - no FLUKA card generated";
1422 AliceInp << "*Generated from call: SetProcess('PFIS',0)";
1428 AliceInp << "*Illegal flag value in SetProcess('PFIS',?) call.";
1430 AliceInp << "*No FLUKA card generated";
1436 // photo electric effect
1437 // G3 default value: 1
1438 // G4 processes: G4PhotoElectricEffect
1439 // G4LowEnergyPhotoElectric
1442 // flag = 0 no photo electric effect
1443 // flag = 1 photo electric effect, electron processed
1444 // flag = 2 photo electric effect, no electron stored
1445 // gMC ->SetProcess("PHOT",1); // EMFCUT 0. -1. 0. 3. lastmat 0. PHOT-THR
1446 else if (strncmp(&sProcessFlag[i][0],"PHOT",4) == 0) {
1447 if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
1450 AliceInp << "*Photo electric effect is activated";
1452 AliceInp << "*Generated from call: SetProcess('PHOT',1);";
1454 AliceInp << setw(10) << "EMFCUT ";
1455 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1456 AliceInp << setw(10) << 0.0; // ignored
1457 AliceInp << setw(10) << -1.0; // resets to default=0.
1458 AliceInp << setw(10) << 0.0; // ignored
1459 AliceInp << setw(10) << 3.0; // upper bound of the material indices in which the respective thresholds apply
1460 AliceInp << setprecision(2);
1461 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1462 AliceInp << setprecision(1);
1463 AliceInp << setw(10) << 1.0; // step length in assigning indices
1464 AliceInp << setw(8) << "PHOT-THR";
1467 else if (iProcessValue[i] == 0) {
1470 AliceInp << "*No photo electric effect - no FLUKA card generated";
1472 AliceInp << "*Generated from call: SetProcess('PHOT',0)";
1478 AliceInp << "*Illegal flag value in SetProcess('PHOT',?) call.";
1480 AliceInp << "*No FLUKA card generated";
1483 } // else if (strncmp(&sProcessFlag[i][0],"PHOT",4) == 0)
1486 // Rayleigh scattering
1487 // G3 default value: 0
1488 // G4 process: G4OpRayleigh
1490 // Particles: optical photon
1492 // flag = 0 Rayleigh scattering off
1493 // flag = 1 Rayleigh scattering on
1494 //xx gMC ->SetProcess("RAYL",1);
1495 else if (strncmp(&sProcessFlag[i][0],"RAYL",4) == 0) {
1496 if (iProcessValue[i] == 1) {
1499 AliceInp << "*Rayleigh scattering is ON by default in FLUKA";
1501 AliceInp << "*No FLUKA card generated";
1504 else if (iProcessValue[i] == 0) {
1507 AliceInp << "*Rayleigh scattering is set OFF";
1509 AliceInp << "*Generated from call: SetProcess('RAYL',0);";
1511 AliceInp << setw(10) << "EMFRAY ";
1512 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1513 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1514 AliceInp << setw(10) << -1.0; // no Rayleigh scattering and no binding corrections for Compton
1515 AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1516 AliceInp << setprecision(2);
1517 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1523 AliceInp << "*Illegal flag value in SetProcess('RAYL',?) call.";
1525 AliceInp << "*No FLUKA card generated";
1528 } // end of else if (strncmp(&sProcessFlag[i][0],"RAYL",4) == 0)
1531 // synchrotron radiation in magnetic field
1532 // G3 default value: 0
1533 // G4 process: G4SynchrotronRadiation
1537 // flag = 0 no synchrotron radiation
1538 // flag = 1 synchrotron radiation
1539 //xx gMC ->SetProcess("SYNC",1); // synchrotron radiation generation
1540 else if (strncmp(&sProcessFlag[i][0],"SYNC",4) == 0) {
1543 AliceInp << "*Synchrotron radiation generation is NOT implemented in FLUKA";
1545 AliceInp << "*No FLUKA card generated";
1550 // Automatic calculation of tracking medium parameters
1551 // flag = 0 no automatic calculation
1552 // flag = 1 automatic calculation
1553 //xx gMC ->SetProcess("AUTO",1); // ??? automatic computation of the tracking medium parameters
1554 else if (strncmp(&sProcessFlag[i][0],"AUTO",4) == 0) {
1557 AliceInp << "*Automatic calculation of tracking medium parameters is always ON in FLUKA";
1559 AliceInp << "*No FLUKA card generated";
1564 // To control energy loss fluctuation model
1565 // flag = 0 Urban model
1566 // flag = 1 PAI model
1567 // flag = 2 PAI+ASHO model (not active at the moment)
1568 //xx gMC ->SetProcess("STRA",1); // ??? energy fluctuation model
1569 else if (strncmp(&sProcessFlag[i][0],"STRA",4) == 0) {
1570 if (iProcessValue[i] == 0 || iProcessValue[i] == 2 || iProcessValue[i] == 3) {
1573 AliceInp << "*Ionization energy losses calculation is activated";
1575 AliceInp << "*Generated from call: SetProcess('STRA',n);, n=0,1,2";
1577 AliceInp << setw(10) << "IONFLUCT ";
1578 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1579 AliceInp << setw(10) << 1.0; // restricted energy loss fluctuations
1580 // (for hadrons and muons) switched on
1581 AliceInp << setw(10) << 1.0; // restricted energy loss fluctuations
1582 // (for e+ and e-) switched on
1583 AliceInp << setw(10) << 1.0; // minimal accuracy
1584 AliceInp << setw(10) << 3.0; // upper bound of the material indices in
1585 // which the respective thresholds apply
1586 AliceInp << setprecision(2);
1587 AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1588 AliceInp << setprecision(1);
1589 AliceInp << setw(10) << 1.0; // step length in assigning indices
1595 AliceInp << "*Illegal flag value in SetProcess('STRA',?) call.";
1597 AliceInp << "*No FLUKA card generated";
1600 } // else if (strncmp(&sProcessFlag[i][0],"STRA",4) == 0)
1605 else { // processes not yet treated
1607 // light photon absorption (Cerenkov photons)
1608 // it is turned on when Cerenkov process is turned on
1609 // G3 default value: 0
1610 // G4 process: G4OpAbsorption, G4OpBoundaryProcess
1612 // Particles: optical photon
1614 // flag = 0 no absorption of Cerenkov photons
1615 // flag = 1 absorption of Cerenkov photons
1616 // gMC ->SetProcess("LABS",2); // ??? Cerenkov light absorption
1620 cout << "SetProcess for flag=" << &sProcessFlag[i][0] << " value=" << iProcessValue[i] << " not yet implemented!" << endl;
1622 } //end of loop number of SetProcess calls
1625 // Loop over number of SetCut calls
1626 for (Int_t i=0; i<iNbOfCut; i++) {
1628 // cuts used in SetProcess calls
1629 if (strncmp(&sCutFlag[i][0],"BCUTM",5) == 0) continue;
1630 else if (strncmp(&sCutFlag[i][0],"BCUTE",5) == 0) continue;
1631 else if (strncmp(&sCutFlag[i][0],"DCUTM",5) == 0) continue;
1632 else if (strncmp(&sCutFlag[i][0],"PPCUTM",6) == 0) continue;
1635 // G4 particles: "gamma"
1636 // G3 default value: 0.001 GeV
1637 //gMC ->SetCut("CUTGAM",cut); // cut for gammas
1638 else if (strncmp(&sCutFlag[i][0],"CUTGAM",6) == 0) {
1641 AliceInp << "*Cut for gamma";
1643 AliceInp << "*Generated from call: SetCut('CUTGAM',cut);";
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(1);
1649 AliceInp << setw(10) << 7.0;
1654 // G4 particles: "e-"
1656 // G3 default value: 0.001 GeV
1657 //gMC ->SetCut("CUTELE",cut); // cut for e+,e-
1658 else if (strncmp(&sCutFlag[i][0],"CUTELE",6) == 0) {
1661 AliceInp << "*Cut for electrons";
1663 AliceInp << "*Generated from call: SetCut('CUTELE',cut);";
1665 AliceInp << setw(10) << "PART-THR ";
1666 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1667 AliceInp << setw(10) << -fCutValue[i];
1668 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1669 AliceInp << setw(10) << 3.0;
1670 AliceInp << setw(10) << 4.0;
1671 AliceInp << setw(10) << 1.0;
1676 // G4 particles: of type "baryon", "meson", "nucleus" with zero charge
1677 // G3 default value: 0.01 GeV
1678 //gMC ->SetCut("CUTNEU",cut); // cut for neutral hadrons
1679 else if (strncmp(&sCutFlag[i][0],"CUTNEU",6) == 0) {
1682 AliceInp << "*Cut for neutral hadrons";
1684 AliceInp << "*Generated from call: SetCut('CUTNEU',cut);";
1686 AliceInp << setw(10) << "PART-THR ";
1687 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1688 AliceInp << setw(10) << -fCutValue[i];
1689 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1690 AliceInp << setw(10) << 8.0; // Neutron
1691 AliceInp << setw(10) << 9.0; // Antineutron
1693 AliceInp << setw(10) << "PART-THR ";
1694 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1695 AliceInp << setw(10) << -fCutValue[i];
1696 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1697 AliceInp << setw(10) << 12.0; // Kaon zero long
1698 AliceInp << setw(10) << 12.0; // Kaon zero long
1700 AliceInp << setw(10) << "PART-THR ";
1701 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1702 AliceInp << setw(10) << -fCutValue[i];
1703 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1704 AliceInp << setw(10) << 17.0; // Lambda, 18=Antilambda
1705 AliceInp << setw(10) << 19.0; // Kaon zero short
1707 AliceInp << setw(10) << "PART-THR ";
1708 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1709 AliceInp << setw(10) << -fCutValue[i];
1710 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1711 AliceInp << setw(10) << 22.0; // Sigma zero, Pion zero, Kaon zero
1712 AliceInp << setw(10) << 25.0; // Antikaon zero
1714 AliceInp << setw(10) << "PART-THR ";
1715 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1716 AliceInp << setw(10) << -fCutValue[i];
1717 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1718 AliceInp << setw(10) << 32.0; // Antisigma zero
1719 AliceInp << setw(10) << 32.0; // Antisigma zero
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(2);
1725 AliceInp << setw(10) << 34.0; // Xi zero
1726 AliceInp << setw(10) << 35.0; // AntiXi zero
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) << 47.0; // D zero
1733 AliceInp << setw(10) << 48.0; // AntiD zero
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) << 53.0; // Xi_c zero
1740 AliceInp << setw(10) << 53.0; // Xi_c zero
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) << 55.0; // Xi'_c zero
1747 AliceInp << setw(10) << 56.0; // Omega_c zero
1749 AliceInp << setw(10) << "PART-THR ";
1750 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1751 AliceInp << setw(10) << -fCutValue[i];
1752 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1753 AliceInp << setw(10) << 59.0; // AntiXi_c zero
1754 AliceInp << setw(10) << 59.0; // AntiXi_c zero
1756 AliceInp << setw(10) << "PART-THR ";
1757 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1758 AliceInp << setw(10) << -fCutValue[i];
1759 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1760 AliceInp << setw(10) << 61.0; // AntiXi'_c zero
1761 AliceInp << setw(10) << 62.0; // AntiOmega_c zero
1766 // G4 particles: of type "baryon", "meson", "nucleus" with non-zero charge
1767 // G3 default value: 0.01 GeV
1768 //gMC ->SetCut("CUTHAD",cut); // cut for charged hadrons
1769 else if (strncmp(&sCutFlag[i][0],"CUTHAD",6) == 0) {
1772 AliceInp << "*Cut for charged hadrons";
1774 AliceInp << "*Generated from call: SetCut('CUTHAD',cut);";
1776 AliceInp << setw(10) << "PART-THR ";
1777 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1778 AliceInp << setw(10) << -fCutValue[i];
1779 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1780 AliceInp << setw(10) << 1.0; // Proton
1781 AliceInp << setw(10) << 2.0; // Antiproton
1783 AliceInp << setw(10) << "PART-THR ";
1784 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1785 AliceInp << setw(10) << -fCutValue[i];
1786 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1787 AliceInp << setw(10) << 13.0; // Positive Pion, Negative Pion, Positive Kaon
1788 AliceInp << setw(10) << 16.0; // Negative Kaon
1790 AliceInp << setw(10) << "PART-THR ";
1791 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1792 AliceInp << setw(10) << -fCutValue[i];
1793 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1794 AliceInp << setw(10) << 20.0; // Negative Sigma
1795 AliceInp << setw(10) << 16.0; // Positive Sigma
1797 AliceInp << setw(10) << "PART-THR ";
1798 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1799 AliceInp << setw(10) << -fCutValue[i];
1800 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1801 AliceInp << setw(10) << 31.0; // Antisigma minus
1802 AliceInp << setw(10) << 33.0; // Antisigma plus
1803 AliceInp << setprecision(1);
1804 AliceInp << setw(10) << 2.0; // step length
1806 AliceInp << setw(10) << "PART-THR ";
1807 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1808 AliceInp << setw(10) << -fCutValue[i];
1809 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1810 AliceInp << setw(10) << 36.0; // Negative Xi, Positive Xi, Omega minus
1811 AliceInp << setw(10) << 39.0; // Antiomega
1813 AliceInp << setw(10) << "PART-THR ";
1814 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1815 AliceInp << setw(10) << -fCutValue[i];
1816 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1817 AliceInp << setw(10) << 45.0; // D plus
1818 AliceInp << setw(10) << 46.0; // D minus
1820 AliceInp << setw(10) << "PART-THR ";
1821 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1822 AliceInp << setw(10) << -fCutValue[i];
1823 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1824 AliceInp << setw(10) << 49.0; // D_s plus, D_s minus, Lambda_c plus
1825 AliceInp << setw(10) << 52.0; // Xi_c plus
1827 AliceInp << setw(10) << "PART-THR ";
1828 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1829 AliceInp << setw(10) << -fCutValue[i];
1830 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1831 AliceInp << setw(10) << 54.0; // Xi'_c plus
1832 AliceInp << setw(10) << 60.0; // AntiXi'_c minus
1833 AliceInp << setprecision(1);
1834 AliceInp << setw(10) << 6.0; // step length
1836 AliceInp << setw(10) << "PART-THR ";
1837 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1838 AliceInp << setw(10) << -fCutValue[i];
1839 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1840 AliceInp << setw(10) << 57.0; // Antilambda_c minus
1841 AliceInp << setw(10) << 58.0; // AntiXi_c minus
1846 // G4 particles: "mu+", "mu-"
1847 // G3 default value: 0.01 GeV
1848 //gMC ->SetCut("CUTMUO",cut); // cut for mu+, mu-
1849 else if (strncmp(&sCutFlag[i][0],"CUTMUO",6) == 0) {
1852 AliceInp << "*Cut for muons";
1854 AliceInp << "*Generated from call: SetCut('CUTMUO',cut);";
1856 AliceInp << setw(10) << "PART-THR ";
1857 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1858 AliceInp << setw(10) << -fCutValue[i];
1859 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1860 AliceInp << setprecision(2);
1861 AliceInp << setw(10) << 10.0;
1862 AliceInp << setw(10) << 11.0;
1865 // delta-rays by electrons
1866 // G4 particles: "e-"
1867 // G3 default value: 10**4 GeV
1868 // gMC ->SetCut("DCUTE",cut); // cut for deltarays by electrons ???????????????
1869 else if (strncmp(&sCutFlag[i][0],"DCUTE",5) == 0) {
1872 AliceInp << "*Cut for delta rays by electrons ????????????";
1874 AliceInp << "*Generated from call: SetCut('DCUTE',cut);";
1876 AliceInp << setw(10) << "EMFCUT ";
1877 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1878 AliceInp << setw(10) << -fCutValue[i];
1879 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1880 AliceInp << setw(10) << 0.0;
1881 AliceInp << setw(10) << 0.0;
1882 AliceInp << setw(10) << 3.0;
1883 AliceInp << setprecision(2);
1884 AliceInp << setw(10) << fLastMaterial;
1885 AliceInp << setprecision(1);
1886 AliceInp << setw(10) << 1.0;
1891 // time of flight cut in seconds
1892 // G4 particles: all
1893 // G3 default value: 0.01 GeV
1894 //gMC ->SetCut("TOFMAX",tofmax); // time of flight cuts in seconds
1895 else if (strncmp(&sCutFlag[i][0],"TOFMAX",6) == 0) {
1898 AliceInp << "*Time of flight cuts in seconds";
1900 AliceInp << "*Generated from call: SetCut('TOFMAX',tofmax);";
1902 AliceInp << setw(10) << "TIME-CUT ";
1903 AliceInp << setiosflags(ios::scientific) << setprecision(5);
1904 AliceInp << setw(10) << fCutValue[i]*1.e9;
1905 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1906 AliceInp << setw(10) << 0.0;
1907 AliceInp << setw(10) << 0.0;
1908 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
1909 AliceInp << setprecision(2);
1910 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
1911 AliceInp << setprecision(1);
1912 AliceInp << setw(10) << 1.0; // step length in assigning numbers
1917 cout << "SetCut for flag=" << &sCutFlag[i][0] << " value=" << fCutValue[i] << " not yet implemented!" << endl;
1919 } //end of loop over SeCut calls
1921 // Add START and STOP card
1922 AliceInp << setw(10) << "START ";
1923 AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint);
1924 AliceInp << setw(10) << fEventsPerRun;
1926 AliceInp << setw(10) << "STOP ";
1929 } // end of InitPhysics
1932 //______________________________________________________________________________
1933 void TFluka::SetMaxStep(Double_t)
1935 // SetMaxStep is dummy procedure in TFluka !
1936 if (fVerbosityLevel >=3)
1937 cout << "SetMaxStep is dummy procedure in TFluka !" << endl;
1940 //______________________________________________________________________________
1941 void TFluka::SetMaxNStep(Int_t)
1943 // SetMaxNStep is dummy procedure in TFluka !
1944 if (fVerbosityLevel >=3)
1945 cout << "SetMaxNStep is dummy procedure in TFluka !" << endl;
1948 //______________________________________________________________________________
1949 void TFluka::SetUserDecay(Int_t)
1951 // SetUserDecay is dummy procedure in TFluka !
1952 if (fVerbosityLevel >=3)
1953 cout << "SetUserDecay is dummy procedure in TFluka !" << endl;
1957 // dynamic properties
1959 //______________________________________________________________________________
1960 void TFluka::TrackPosition(TLorentzVector& position) const
1962 // Return the current position in the master reference frame of the
1963 // track being transported
1964 // TRACKR.atrack = age of the particle
1965 // TRACKR.xtrack = x-position of the last point
1966 // TRACKR.ytrack = y-position of the last point
1967 // TRACKR.ztrack = z-position of the last point
1968 Int_t caller = GetCaller();
1969 if (caller == 1 || caller == 3 || caller == 6 || caller == 11 || caller == 12) { //bxdraw,endraw,usdraw
1970 position.SetX(GetXsco());
1971 position.SetY(GetYsco());
1972 position.SetZ(GetZsco());
1973 position.SetT(TRACKR.atrack);
1975 else if (caller == 4) { // mgdraw
1976 position.SetX(TRACKR.xtrack[TRACKR.ntrack]);
1977 position.SetY(TRACKR.ytrack[TRACKR.ntrack]);
1978 position.SetZ(TRACKR.ztrack[TRACKR.ntrack]);
1979 position.SetT(TRACKR.atrack);
1981 else if (caller == 5) { // sodraw
1982 position.SetX(TRACKR.xtrack[TRACKR.ntrack]);
1983 position.SetY(TRACKR.ytrack[TRACKR.ntrack]);
1984 position.SetZ(TRACKR.ztrack[TRACKR.ntrack]);
1988 Warning("TrackPosition","position not available");
1991 //______________________________________________________________________________
1992 void TFluka::TrackPosition(Double_t& x, Double_t& y, Double_t& z) const
1994 // Return the current position in the master reference frame of the
1995 // track being transported
1996 // TRACKR.atrack = age of the particle
1997 // TRACKR.xtrack = x-position of the last point
1998 // TRACKR.ytrack = y-position of the last point
1999 // TRACKR.ztrack = z-position of the last point
2000 Int_t caller = GetCaller();
2001 if (caller == 1 || caller == 3 || caller == 6 || caller == 11 || caller == 12) { //bxdraw,endraw,usdraw
2006 else if (caller == 4) { // mgdraw
2007 x = TRACKR.xtrack[TRACKR.ntrack];
2008 y = TRACKR.ytrack[TRACKR.ntrack];
2009 z = TRACKR.ztrack[TRACKR.ntrack];
2011 else if (caller == 5) { // sodraw
2012 x = TRACKR.xtrack[TRACKR.ntrack];
2013 y = TRACKR.ytrack[TRACKR.ntrack];
2014 z = TRACKR.ztrack[TRACKR.ntrack];
2017 Warning("TrackPosition","position not available");
2020 //______________________________________________________________________________
2021 void TFluka::TrackMomentum(TLorentzVector& momentum) const
2023 // Return the direction and the momentum (GeV/c) of the track
2024 // currently being transported
2025 // TRACKR.ptrack = momentum of the particle (not always defined, if
2026 // < 0 must be obtained from etrack)
2027 // TRACKR.cx,y,ztrck = direction cosines of the current particle
2028 // TRACKR.etrack = total energy of the particle
2029 // TRACKR.jtrack = identity number of the particle
2030 // PAPROP.am[TRACKR.jtrack] = particle mass in gev
2031 Int_t caller = GetCaller();
2032 if (caller != 2) { // not eedraw
2033 if (TRACKR.ptrack >= 0) {
2034 momentum.SetPx(TRACKR.ptrack*TRACKR.cxtrck);
2035 momentum.SetPy(TRACKR.ptrack*TRACKR.cytrck);
2036 momentum.SetPz(TRACKR.ptrack*TRACKR.cztrck);
2037 momentum.SetE(TRACKR.etrack);
2041 Double_t p = sqrt(TRACKR.etrack*TRACKR.etrack - PAPROP.am[TRACKR.jtrack+6]*PAPROP.am[TRACKR.jtrack+6]);
2042 momentum.SetPx(p*TRACKR.cxtrck);
2043 momentum.SetPy(p*TRACKR.cytrck);
2044 momentum.SetPz(p*TRACKR.cztrck);
2045 momentum.SetE(TRACKR.etrack);
2050 Warning("TrackMomentum","momentum not available");
2053 //______________________________________________________________________________
2054 void TFluka::TrackMomentum(Double_t& px, Double_t& py, Double_t& pz, Double_t& e) const
2056 // Return the direction and the momentum (GeV/c) of the track
2057 // currently being transported
2058 // TRACKR.ptrack = momentum of the particle (not always defined, if
2059 // < 0 must be obtained from etrack)
2060 // TRACKR.cx,y,ztrck = direction cosines of the current particle
2061 // TRACKR.etrack = total energy of the particle
2062 // TRACKR.jtrack = identity number of the particle
2063 // PAPROP.am[TRACKR.jtrack] = particle mass in gev
2064 Int_t caller = GetCaller();
2065 if (caller != 2) { // not eedraw
2066 if (TRACKR.ptrack >= 0) {
2067 px = TRACKR.ptrack*TRACKR.cxtrck;
2068 py = TRACKR.ptrack*TRACKR.cytrck;
2069 pz = TRACKR.ptrack*TRACKR.cztrck;
2074 Double_t p = sqrt(TRACKR.etrack*TRACKR.etrack - PAPROP.am[TRACKR.jtrack+6]*PAPROP.am[TRACKR.jtrack+6]);
2075 px = p*TRACKR.cxtrck;
2076 py = p*TRACKR.cytrck;
2077 pz = p*TRACKR.cztrck;
2083 Warning("TrackMomentum","momentum not available");
2086 //______________________________________________________________________________
2087 Double_t TFluka::TrackStep() const
2089 // Return the length in centimeters of the current step
2090 // TRACKR.ctrack = total curved path
2091 Int_t caller = GetCaller();
2092 if (caller == 1 || caller == 3 || caller == 6) //bxdraw,endraw,usdraw
2094 else if (caller == 4) //mgdraw
2095 return TRACKR.ctrack;
2100 //______________________________________________________________________________
2101 Double_t TFluka::TrackLength() const
2103 // TRACKR.cmtrck = cumulative curved path since particle birth
2104 Int_t caller = GetCaller();
2105 if (caller == 1 || caller == 3 || caller == 4 || caller == 6) //bxdraw,endraw,mgdraw,usdraw
2106 return TRACKR.cmtrck;
2111 //______________________________________________________________________________
2112 Double_t TFluka::TrackTime() const
2114 // Return the current time of flight of the track being transported
2115 // TRACKR.atrack = age of the particle
2116 Int_t caller = GetCaller();
2117 if (caller == 1 || caller == 3 || caller == 4 || caller == 6) //bxdraw,endraw,mgdraw,usdraw
2118 return TRACKR.atrack;
2123 //______________________________________________________________________________
2124 Double_t TFluka::Edep() const
2126 // Energy deposition
2127 // if TRACKR.ntrack = 0, TRACKR.mtrack = 0:
2128 // -->local energy deposition (the value and the point are not recorded in TRACKR)
2129 // but in the variable "rull" of the procedure "endraw.cxx"
2130 // if TRACKR.ntrack > 0, TRACKR.mtrack = 0:
2131 // -->no energy loss along the track
2132 // if TRACKR.ntrack > 0, TRACKR.mtrack > 0:
2133 // -->energy loss distributed along the track
2134 // TRACKR.dtrack = energy deposition of the jth deposition even
2136 for ( Int_t j=0;j<TRACKR.mtrack;j++) {
2137 sum +=TRACKR.dtrack[j];
2139 if (TRACKR.ntrack == 0 && TRACKR.mtrack == 0)
2146 //______________________________________________________________________________
2147 Int_t TFluka::TrackPid() const
2149 // Return the id of the particle transported
2150 // TRACKR.jtrack = identity number of the particle
2151 Int_t caller = GetCaller();
2152 if (caller != 2) // not eedraw
2153 return PDGFromId(TRACKR.jtrack);
2158 //______________________________________________________________________________
2159 Double_t TFluka::TrackCharge() const
2161 // Return charge of the track currently transported
2162 // PAPROP.ichrge = electric charge of the particle
2163 // TRACKR.jtrack = identity number of the particle
2164 Int_t caller = GetCaller();
2165 if (caller != 2) // not eedraw
2166 return PAPROP.ichrge[TRACKR.jtrack+6];
2171 //______________________________________________________________________________
2172 Double_t TFluka::TrackMass() const
2174 // PAPROP.am = particle mass in GeV
2175 // TRACKR.jtrack = identity number of the particle
2176 Int_t caller = GetCaller();
2177 if (caller != 2) // not eedraw
2178 return PAPROP.am[TRACKR.jtrack+6];
2183 //______________________________________________________________________________
2184 Double_t TFluka::Etot() const
2186 // TRACKR.etrack = total energy of the particle
2187 Int_t caller = GetCaller();
2188 if (caller != 2) // not eedraw
2189 return TRACKR.etrack;
2197 //______________________________________________________________________________
2198 Bool_t TFluka::IsNewTrack() const
2200 // True if the track has positive cummulative length
2201 Int_t caller = GetCaller();
2202 if (caller != 2) { // not eedraw
2203 if (TRACKR.cmtrck > 0.0)
2212 //______________________________________________________________________________
2213 Bool_t TFluka::IsTrackInside() const
2215 // True if the track is not at the boundary of the current volume
2216 // In Fluka a step is always inside one kind of material
2217 // If the step would go behind the region of one material,
2218 // it will be shortened to reach only the boundary.
2219 // Therefore IsTrackInside() is always true.
2220 Int_t caller = GetCaller();
2221 if (caller == 1) // bxdraw
2227 //______________________________________________________________________________
2228 Bool_t TFluka::IsTrackEntering() const
2230 // True if this is the first step of the track in the current volume
2232 Int_t caller = GetCaller();
2233 if (caller == 11) // bxdraw entering
2238 //______________________________________________________________________________
2239 Bool_t TFluka::IsTrackExiting() const
2241 Int_t caller = GetCaller();
2242 if (caller == 12) // bxdraw exiting
2247 //______________________________________________________________________________
2248 Bool_t TFluka::IsTrackOut() const
2250 // True if the track is out of the setup
2252 // Icode = 14: escape - call from Kaskad
2253 // Icode = 23: escape - call from Emfsco
2254 // Icode = 32: escape - call from Kasneu
2255 // Icode = 40: escape - call from Kashea
2256 // Icode = 51: escape - call from Kasoph
2261 fIcode == 51) return 1;
2265 //______________________________________________________________________________
2266 Bool_t TFluka::IsTrackDisappeared() const
2268 // means all inelastic interactions and decays
2269 // fIcode from usdraw
2270 if (fIcode == 101 || // inelastic interaction
2271 fIcode == 102 || // particle decay
2272 fIcode == 214 || // in-flight annihilation
2273 fIcode == 215 || // annihilation at rest
2274 fIcode == 217 || // pair production
2275 fIcode == 221) return 1;
2279 //______________________________________________________________________________
2280 Bool_t TFluka::IsTrackStop() const
2282 // True if the track energy has fallen below the threshold
2283 // means stopped by signal or below energy threshold
2284 // Icode = 12: stopping particle - call from Kaskad
2285 // Icode = 15: time kill - call from Kaskad
2286 // Icode = 21: below threshold, iarg=1 - call from Emfsco
2287 // Icode = 22: below threshold, iarg=2 - call from Emfsco
2288 // Icode = 24: time kill - call from Emfsco
2289 // Icode = 31: below threshold - call from Kasneu
2290 // Icode = 33: time kill - call from Kasneu
2291 // Icode = 41: time kill - call from Kashea
2292 // Icode = 52: time kill - call from Kasoph
2301 fIcode == 52) return 1;
2305 //______________________________________________________________________________
2306 Bool_t TFluka::IsTrackAlive() const
2308 // means not disappeared or not out
2309 if (IsTrackDisappeared() || IsTrackOut() ) return 0;
2317 //______________________________________________________________________________
2318 Int_t TFluka::NSecondaries() const
2319 // Number of secondary particles generated in the current step
2320 // FINUC.np = number of secondaries except light and heavy ions
2321 // FHEAVY.npheav = number of secondaries for light and heavy secondary ions
2323 Int_t caller = GetCaller();
2324 if (caller == 6) // valid only after usdraw
2325 return FINUC.np + FHEAVY.npheav;
2328 } // end of NSecondaries
2330 //______________________________________________________________________________
2331 void TFluka::GetSecondary(Int_t isec, Int_t& particleId,
2332 TLorentzVector& position, TLorentzVector& momentum)
2334 Int_t caller = GetCaller();
2335 if (caller == 6) { // valid only after usdraw
2336 if (isec >= 0 && isec < FINUC.np) {
2337 particleId = PDGFromId(FINUC.kpart[isec]);
2338 position.SetX(fXsco);
2339 position.SetY(fYsco);
2340 position.SetZ(fZsco);
2341 position.SetT(TRACKR.atrack);
2342 // position.SetT(TRACKR.atrack+FINUC.agesec[isec]); //not yet implem.
2343 momentum.SetPx(FINUC.plr[isec]*FINUC.cxr[isec]);
2344 momentum.SetPy(FINUC.plr[isec]*FINUC.cyr[isec]);
2345 momentum.SetPz(FINUC.plr[isec]*FINUC.czr[isec]);
2346 momentum.SetE(FINUC.tki[isec] + PAPROP.am[FINUC.kpart[isec]+6]);
2348 else if (isec >= FINUC.np && isec < FINUC.np + FHEAVY.npheav) {
2349 Int_t jsec = isec - FINUC.np;
2350 particleId = FHEAVY.kheavy[jsec]; // this is Fluka id !!!
2351 position.SetX(fXsco);
2352 position.SetY(fYsco);
2353 position.SetZ(fZsco);
2354 position.SetT(TRACKR.atrack);
2355 // position.SetT(TRACKR.atrack+FHEAVY.agheav[jsec]); //not yet implem.
2356 momentum.SetPx(FHEAVY.pheavy[jsec]*FHEAVY.cxheav[jsec]);
2357 momentum.SetPy(FHEAVY.pheavy[jsec]*FHEAVY.cyheav[jsec]);
2358 momentum.SetPz(FHEAVY.pheavy[jsec]*FHEAVY.czheav[jsec]);
2359 if (FHEAVY.tkheav[jsec] >= 3 && FHEAVY.tkheav[jsec] <= 6)
2360 momentum.SetE(FHEAVY.tkheav[jsec] + PAPROP.am[jsec+6]);
2361 else if (FHEAVY.tkheav[jsec] > 6)
2362 momentum.SetE(FHEAVY.tkheav[jsec] + FHEAVY.amnhea[jsec]); // to be checked !!!
2365 Warning("GetSecondary","isec out of range");
2368 Warning("GetSecondary","no secondaries available");
2369 } // end of GetSecondary
2371 //______________________________________________________________________________
2372 TMCProcess TFluka::ProdProcess(Int_t) const
2373 // Name of the process that has produced the secondary particles
2374 // in the current step
2376 const TMCProcess kIpNoProc = kPNoProcess;
2377 const TMCProcess kIpPDecay = kPDecay;
2378 const TMCProcess kIpPPair = kPPair;
2379 // const TMCProcess kIpPPairFromPhoton = kPPairFromPhoton;
2380 // const TMCProcess kIpPPairFromVirtualPhoton = kPPairFromVirtualPhoton;
2381 const TMCProcess kIpPCompton = kPCompton;
2382 const TMCProcess kIpPPhotoelectric = kPPhotoelectric;
2383 const TMCProcess kIpPBrem = kPBrem;
2384 // const TMCProcess kIpPBremFromHeavy = kPBremFromHeavy;
2385 // const TMCProcess kIpPBremFromElectronOrPositron = kPBremFromElectronOrPositron;
2386 const TMCProcess kIpPDeltaRay = kPDeltaRay;
2387 // const TMCProcess kIpPMoller = kPMoller;
2388 // const TMCProcess kIpPBhabha = kPBhabha;
2389 const TMCProcess kIpPAnnihilation = kPAnnihilation;
2390 // const TMCProcess kIpPAnnihilInFlight = kPAnnihilInFlight;
2391 // const TMCProcess kIpPAnnihilAtRest = kPAnnihilAtRest;
2392 const TMCProcess kIpPHadronic = kPHadronic;
2393 const TMCProcess kIpPMuonNuclear = kPMuonNuclear;
2394 const TMCProcess kIpPPhotoFission = kPPhotoFission;
2395 const TMCProcess kIpPRayleigh = kPRayleigh;
2396 // const TMCProcess kIpPCerenkov = kPCerenkov;
2397 // const TMCProcess kIpPSynchrotron = kPSynchrotron;
2399 Int_t mugamma = TRACKR.jtrack == 7 || TRACKR.jtrack == 10 || TRACKR.jtrack == 11;
2400 if (fIcode == 102) return kIpPDecay;
2401 else if (fIcode == 104 || fIcode == 217) return kIpPPair;
2402 // else if (fIcode == 104) return kIpPairFromPhoton;
2403 // else if (fIcode == 217) return kIpPPairFromVirtualPhoton;
2404 else if (fIcode == 219) return kIpPCompton;
2405 else if (fIcode == 221) return kIpPPhotoelectric;
2406 else if (fIcode == 105 || fIcode == 208) return kIpPBrem;
2407 // else if (fIcode == 105) return kIpPBremFromHeavy;
2408 // else if (fIcode == 208) return kPBremFromElectronOrPositron;
2409 else if (fIcode == 103 || fIcode == 400) return kIpPDeltaRay;
2410 else if (fIcode == 210 || fIcode == 212) return kIpPDeltaRay;
2411 // else if (fIcode == 210) return kIpPMoller;
2412 // else if (fIcode == 212) return kIpPBhabha;
2413 else if (fIcode == 214 || fIcode == 215) return kIpPAnnihilation;
2414 // else if (fIcode == 214) return kIpPAnnihilInFlight;
2415 // else if (fIcode == 215) return kIpPAnnihilAtRest;
2416 else if (fIcode == 101) return kIpPHadronic;
2417 else if (fIcode == 101) {
2418 if (!mugamma) return kIpPHadronic;
2419 else if (TRACKR.jtrack == 7) return kIpPPhotoFission;
2420 else return kIpPMuonNuclear;
2422 else if (fIcode == 225) return kIpPRayleigh;
2423 // Fluka codes 100, 300 and 400 still to be investigasted
2424 else return kIpNoProc;
2427 //Int_t StepProcesses(TArrayI &proc) const
2428 // Return processes active in the current step
2430 //ck = total energy of the particl ????????????????
2434 //______________________________________________________________________________
2435 Int_t TFluka::VolId2Mate(Int_t id) const
2438 // Returns the material number for a given volume ID
2440 return fGeom->VolId2Mate(id);
2443 //______________________________________________________________________________
2444 const char* TFluka::VolName(Int_t id) const
2447 // Returns the volume name for a given volume ID
2449 return fGeom->VolName(id);
2452 //______________________________________________________________________________
2453 Int_t TFluka::VolId(const Text_t* volName) const
2456 // Converts from volume name to volume ID.
2457 // Time consuming. (Only used during set-up)
2458 // Could be replaced by hash-table
2460 return fGeom->VolId(volName);
2463 //______________________________________________________________________________
2464 Int_t TFluka::CurrentVolID(Int_t& copyNo) const
2467 // Return the logical id and copy number corresponding to the current fluka region
2469 return fGeom->CurrentVolID(copyNo);
2472 //______________________________________________________________________________
2473 Int_t TFluka::CurrentVolOffID(Int_t off, Int_t& copyNo) const
2476 // Return the logical id and copy number of off'th mother
2477 // corresponding to the current fluka region
2479 return fGeom->CurrentVolOffID(off, copyNo);
2482 //______________________________________________________________________________
2483 const char* TFluka::CurrentVolName() const
2486 // Return the current volume name
2488 return fGeom->CurrentVolName();
2491 //______________________________________________________________________________
2492 const char* TFluka::CurrentVolOffName(Int_t off) const
2495 // Return the volume name of the off'th mother of the current volume
2497 return fGeom->CurrentVolOffName(off);
2500 //______________________________________________________________________________
2501 Int_t TFluka::CurrentMaterial(Float_t & /*a*/, Float_t & /*z*/,
2502 Float_t & /*dens*/, Float_t & /*radl*/, Float_t & /*absl*/) const
2505 // Return the current medium number ??? what about material properties
2508 Int_t id = TFluka::CurrentVolID(copy);
2509 Int_t med = TFluka::VolId2Mate(id);
2513 //______________________________________________________________________________
2514 void TFluka::Gmtod(Float_t* xm, Float_t* xd, Int_t iflag)
2516 // Transforms a position from the world reference frame
2517 // to the current volume reference frame.
2519 // Geant3 desription:
2520 // ==================
2521 // Computes coordinates XD (in DRS)
2522 // from known coordinates XM in MRS
2523 // The local reference system can be initialized by
2524 // - the tracking routines and GMTOD used in GUSTEP
2525 // - a call to GMEDIA(XM,NUMED)
2526 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2527 // (inverse routine is GDTOM)
2529 // If IFLAG=1 convert coordinates
2530 // IFLAG=2 convert direction cosinus
2533 fGeom->Gmtod(xm,xd,iflag);
2536 //______________________________________________________________________________
2537 void TFluka::Gmtod(Double_t* xm, Double_t* xd, Int_t iflag)
2539 // Transforms a position from the world reference frame
2540 // to the current volume reference frame.
2542 // Geant3 desription:
2543 // ==================
2544 // Computes coordinates XD (in DRS)
2545 // from known coordinates XM in MRS
2546 // The local reference system can be initialized by
2547 // - the tracking routines and GMTOD used in GUSTEP
2548 // - a call to GMEDIA(XM,NUMED)
2549 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2550 // (inverse routine is GDTOM)
2552 // If IFLAG=1 convert coordinates
2553 // IFLAG=2 convert direction cosinus
2556 fGeom->Gmtod(xm,xd,iflag);
2559 //______________________________________________________________________________
2560 void TFluka::Gdtom(Float_t* xd, Float_t* xm, Int_t iflag)
2562 // Transforms a position from the current volume reference frame
2563 // to the world reference frame.
2565 // Geant3 desription:
2566 // ==================
2567 // Computes coordinates XM (Master Reference System
2568 // knowing the coordinates XD (Detector Ref System)
2569 // The local reference system can be initialized by
2570 // - the tracking routines and GDTOM used in GUSTEP
2571 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2572 // (inverse routine is GMTOD)
2574 // If IFLAG=1 convert coordinates
2575 // IFLAG=2 convert direction cosinus
2578 fGeom->Gdtom(xd,xm,iflag);
2581 //______________________________________________________________________________
2582 void TFluka::Gdtom(Double_t* xd, Double_t* xm, Int_t iflag)
2584 // Transforms a position from the current volume reference frame
2585 // to the world reference frame.
2587 // Geant3 desription:
2588 // ==================
2589 // Computes coordinates XM (Master Reference System
2590 // knowing the coordinates XD (Detector Ref System)
2591 // The local reference system can be initialized by
2592 // - the tracking routines and GDTOM used in GUSTEP
2593 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2594 // (inverse routine is GMTOD)
2596 // If IFLAG=1 convert coordinates
2597 // IFLAG=2 convert direction cosinus
2600 fGeom->Gdtom(xd,xm,iflag);
2603 // ===============================================================
2604 void TFluka::FutoTest()
2606 Int_t icode, mreg, newreg, particleId;
2607 Double_t rull, xsco, ysco, zsco;
2608 TLorentzVector position, momentum;
2611 if (fVerbosityLevel >=3)
2612 cout << " icode=" << icode << endl;
2613 } else if (icode > 0 && icode <= 5) {
2616 if (fVerbosityLevel >=3)
2617 cout << " icode=" << icode
2620 TrackPosition(position);
2621 TrackMomentum(momentum);
2622 if (fVerbosityLevel >=3) {
2623 cout << "TLorentzVector positionX=" << position.X()
2624 << "positionY=" << position.Y()
2625 << "positionZ=" << position.Z()
2626 << "timeT=" << position.T() << endl;
2627 cout << "TLorentzVector momentumX=" << momentum.X()
2628 << "momentumY=" << momentum.Y()
2629 << "momentumZ=" << momentum.Z()
2630 << "energyE=" << momentum.E() << endl;
2631 cout << "TrackStep=" << TrackStep() << endl;
2632 cout << "TrackLength=" << TrackLength() << endl;
2633 cout << "TrackTime=" << TrackTime() << endl;
2634 cout << "Edep=" << Edep() << endl;
2635 cout << "TrackPid=" << TrackPid() << endl;
2636 cout << "TrackCharge=" << TrackCharge() << endl;
2637 cout << "TrackMass=" << TrackMass() << endl;
2638 cout << "Etot=" << Etot() << endl;
2639 cout << "IsNewTrack=" << IsNewTrack() << endl;
2640 cout << "IsTrackInside=" << IsTrackInside() << endl;
2641 cout << "IsTrackEntering=" << IsTrackEntering() << endl;
2642 cout << "IsTrackExiting=" << IsTrackExiting() << endl;
2643 cout << "IsTrackOut=" << IsTrackOut() << endl;
2644 cout << "IsTrackDisappeared=" << IsTrackDisappeared() << endl;
2645 cout << "IsTrackAlive=" << IsTrackAlive() << endl;
2648 Float_t x = position.X();
2649 Float_t y = position.Y();
2650 Float_t z = position.Z();
2653 xm[0] = x; xm[1] = y; xm[2] = z;
2654 if (fVerbosityLevel >= 3)
2655 printf("Global trackPosition: %f %f %f \n", x, y, z);
2657 if (fVerbosityLevel >= 3)
2658 printf("Local trackPosition: %f %f %f \n", xd[0], xd[1], xd[2]);
2660 if (fVerbosityLevel >= 3)
2661 printf("New trackPosition: %f %f %f \n", xm[0], xm[1], xm[2]);
2662 } else if((icode >= 10 && icode <= 15) ||
2663 (icode >= 20 && icode <= 24) ||
2664 (icode >= 30 && icode <= 33) ||
2665 (icode >= 40 && icode <= 41) ||
2666 (icode >= 50 && icode <= 52)) {
2674 if (fVerbosityLevel >=3) {
2675 cout << " icode=" << icode
2680 << " zsco=" << zsco << endl;
2682 TrackPosition(position);
2683 TrackMomentum(momentum);
2684 if (fVerbosityLevel >=3) {
2685 cout << "Edep=" << Edep() << endl;
2686 cout << "Etot=" << Etot() << endl;
2687 cout << "TrackPid=" << TrackPid() << endl;
2688 cout << "TrackCharge=" << TrackCharge() << endl;
2689 cout << "TrackMass=" << TrackMass() << endl;
2690 cout << "IsTrackOut=" << IsTrackOut() << endl;
2691 cout << "IsTrackDisappeared=" << IsTrackDisappeared() << endl;
2692 cout << "IsTrackStop=" << IsTrackStop() << endl;
2693 cout << "IsTrackAlive=" << IsTrackAlive() << endl;
2695 } else if((icode >= 100 && icode <= 105) ||
2699 (icode >= 214 && icode <= 215) ||
2712 if (fVerbosityLevel >=3) {
2713 cout << " icode=" << icode
2717 << " zsco=" << zsco << endl;
2718 cout << "TrackPid=" << TrackPid() << endl;
2719 cout << "NSecondaries=" << NSecondaries() << endl;
2722 for (Int_t isec=0; isec< NSecondaries(); isec++) {
2723 TFluka::GetSecondary(isec, particleId, position, momentum);
2724 if (fVerbosityLevel >=3) {
2725 cout << "TLorentzVector positionX=" << position.X()
2726 << "positionY=" << position.Y()
2727 << "positionZ=" << position.Z()
2728 << "timeT=" << position.T() << endl;
2729 cout << "TLorentzVector momentumX=" << momentum.X()
2730 << "momentumY=" << momentum.Y()
2731 << "momentumZ=" << momentum.Z()
2732 << "energyE=" << momentum.E() << endl;
2733 cout << "TrackPid=" << particleId << endl;
2736 } else if((icode == 19) ||
2742 newreg = GetNewreg();
2746 if (fVerbosityLevel >=3) {
2747 cout << " icode=" << icode
2749 << " newreg=" << newreg
2752 << " zsco=" << zsco << endl;
2755 } // end of FutoTest