]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TFluka/TFluka.cxx
b25cc8bd39b4bfdc84facc9c99cf91332881699e
[u/mrichter/AliRoot.git] / TFluka / TFluka.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 #include <Riostream.h>
19
20 #include "TClonesArray.h"
21 #include "TFluka.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
32
33 #include "TVirtualMC.h"
34 #include "TG4GeometryManager.h" //For the geometry management
35 #include "TG4DetConstruction.h" //For the detector construction
36
37 #include "FGeometryInit.hh"
38 #include "TLorentzVector.h"
39 #include "FlukaVolume.h"
40
41 // Fluka methods that may be needed.
42 #ifndef WIN32 
43 # define flukam  flukam_
44 # define fluka_openinp fluka_openinp_
45 # define fluka_closeinp fluka_closeinp_
46 # define mcihad mcihad_
47 # define mpdgha mpdgha_
48 #else 
49 # define flukam  FLUKAM
50 # define fluka_openinp FLUKA_OPENINP
51 # define fluka_closeinp FLUKA_CLOSEINP
52 # define mcihad MCIHAD
53 # define mpdgha MPDGHA
54 #endif
55
56 extern "C" 
57 {
58   //
59   // Prototypes for FLUKA functions
60   //
61   void type_of_call flukam(const int&);
62   void type_of_call fluka_openinp(const int&, DEFCHARA);
63   void type_of_call fluka_closeinp(const int&);
64   int  type_of_call mcihad(const int&);
65   int  type_of_call mpdgha(const int&);
66 }
67
68 //
69 // Class implementation for ROOT
70 //
71 ClassImp(TFluka)
72
73 //
74 //----------------------------------------------------------------------------
75 // TFluka constructors and destructors.
76 //____________________________________________________________________________ 
77 TFluka::TFluka()
78   :TVirtualMC(),
79    fVerbosityLevel(0),
80    sInputFileName(""),
81    fDetector(0),
82    fCurrentFlukaRegion(-1)
83
84   //
85   // Default constructor
86   //
87
88  
89 TFluka::TFluka(const char *title, Int_t verbosity, Bool_t isRootGeometrySupported)
90   :TVirtualMC("TFluka",title, isRootGeometrySupported),
91    fVerbosityLevel(verbosity),
92    sInputFileName(""),
93    fTrackIsEntering(0),
94    fTrackIsExiting(0),
95    fTrackIsNew(0),
96    fDetector(0),
97    fCurrentFlukaRegion(-1)
98 {
99   if (fVerbosityLevel >=3)
100     cout << "==> TFluka::TFluka(" << title << ") constructor called." << endl;
101
102   
103   // create geometry manager
104   if (fVerbosityLevel >=2)
105     cout << "\t* Creating G4 Geometry manager..." << endl;
106   fGeometryManager = new TG4GeometryManager();
107   if (fVerbosityLevel >=2)
108     cout << "\t* Creating G4 Detector..." << endl;
109   fDetector = new TG4DetConstruction();
110   FGeometryInit* geominit = FGeometryInit::GetInstance();
111   if (geominit)
112     geominit->setDetConstruction(fDetector);
113   else {
114     cerr << "ERROR: Could not create FGeometryInit!" << endl;
115     cerr << "       Exiting!!!" << endl;
116     abort();
117   }
118
119   if (fVerbosityLevel >=3)
120     cout << "<== TFluka::TFluka(" << title << ") constructor called." << endl;
121
122   fVolumeMediaMap = new TClonesArray("FlukaVolume",1000);
123   fNVolumes      = 0;
124   fMediaByRegion = 0;
125 }
126
127 TFluka::~TFluka() {
128   if (fVerbosityLevel >=3)
129     cout << "==> TFluka::~TFluka() destructor called." << endl;
130
131   delete fGeometryManager;
132   fVolumeMediaMap->Delete();
133   delete  fVolumeMediaMap;
134   
135
136   if (fVerbosityLevel >=3)
137     cout << "<== TFluka::~TFluka() destructor called." << endl;
138 }
139
140 //
141 //_____________________________________________________________________________
142 // TFluka control methods
143 //____________________________________________________________________________ 
144 void TFluka::Init() {
145
146   FGeometryInit* geominit = FGeometryInit::GetInstance();
147   if (fVerbosityLevel >=3)
148     cout << "==> TFluka::Init() called." << endl;
149
150   cout << "\t* InitPhysics() - Prepare input file to be called" << endl; 
151   geominit->Init();
152   // now we have G4 geometry created and we have to patch alice.inp
153   // with the material mapping file FlukaMat.inp
154   InitPhysics(); // prepare input file with the current physics settings
155     cout << "\t* InitPhysics() - Prepare input file was called" << endl; 
156
157   if (fVerbosityLevel >=2)
158     cout << "\t* Changing lfdrtr = (" << (GLOBAL.lfdrtr?'T':'F')
159          << ") in fluka..." << endl;
160   GLOBAL.lfdrtr = true;
161
162   if (fVerbosityLevel >=2)
163     cout << "\t* Opening file " << sInputFileName << endl;
164   const char* fname = sInputFileName;
165   fluka_openinp(lunin, PASSCHARA(fname));
166
167   if (fVerbosityLevel >=2)
168     cout << "\t* Calling flukam..." << endl;
169   flukam(1);
170
171   if (fVerbosityLevel >=2)
172     cout << "\t* Closing file " << sInputFileName << endl;
173   fluka_closeinp(lunin);
174
175   FinishGeometry();
176
177   if (fVerbosityLevel >=3)
178     cout << "<== TFluka::Init() called." << endl;
179
180 }
181
182 void TFluka::FinishGeometry() {
183 //
184 // Build-up table with region to medium correspondance
185 //
186     char tmp[5];
187     
188   if (fVerbosityLevel >=3)
189     cout << "==> TFluka::FinishGeometry() called." << endl;
190
191 //  fGeometryManager->Ggclos();
192
193   FGeometryInit* flugg = FGeometryInit::GetInstance();  
194   
195   fMediaByRegion = new Int_t[fNVolumes+2];
196   for (Int_t i = 0; i < fNVolumes; i++)
197   {
198       FlukaVolume* vol = dynamic_cast<FlukaVolume*>((*fVolumeMediaMap)[i]);
199       TString volName = vol->GetName();
200       Int_t   media   = vol->GetMedium();
201       if (fVerbosityLevel >= 3)
202       printf("Finish Geometry: volName, media %d %s %d \n", i, volName.Data(), media);
203       strcpy(tmp, volName.Data());
204       tmp[4] = '\0';
205       flugg->SetMediumFromName(tmp, media, i+1);
206       fMediaByRegion[i] = media;
207   }
208
209   flugg->BuildMediaMap();
210   
211   if (fVerbosityLevel >=3)
212     cout << "<== TFluka::FinishGeometry() called." << endl;
213
214
215 void TFluka::BuildPhysics() {
216   if (fVerbosityLevel >=3)
217     cout << "==> TFluka::BuildPhysics() called." << endl;
218
219
220   if (fVerbosityLevel >=3)
221     cout << "<== TFluka::BuildPhysics() called." << endl;
222 }  
223
224 void TFluka::ProcessEvent() {
225   if (fVerbosityLevel >=3)
226     cout << "==> TFluka::ProcessEvent() called." << endl;
227   fApplication->GeneratePrimaries();
228   EPISOR.lsouit = true;
229   flukam(1);
230   if (fVerbosityLevel >=3)
231     cout << "<== TFluka::ProcessEvent() called." << endl;
232 }
233
234
235 void TFluka::ProcessRun(Int_t nevent) {
236   if (fVerbosityLevel >=3)
237     cout << "==> TFluka::ProcessRun(" << nevent << ") called." 
238          << endl;
239
240   if (fVerbosityLevel >=2) {
241     cout << "\t* GLOBAL.fdrtr = " << (GLOBAL.lfdrtr?'T':'F') << endl;
242     cout << "\t* Calling flukam again..." << endl;
243   }
244   fApplication->InitGeometry();
245   fApplication->BeginEvent();
246   ProcessEvent();
247   fApplication->FinishEvent();
248   if (fVerbosityLevel >=3)
249     cout << "<== TFluka::ProcessRun(" << nevent << ") called." 
250          << endl;
251
252 }
253
254 //_____________________________________________________________________________
255 // methods for building/management of geometry
256 //____________________________________________________________________________ 
257 // functions from GCONS 
258 void TFluka::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,  
259                     Float_t &dens, Float_t &radl, Float_t &absl,
260                     Float_t* ubuf, Int_t& nbuf) {
261 //
262   fGeometryManager->Gfmate(imat, name, a, z, dens, radl, absl, ubuf, nbuf);
263
264
265 void TFluka::Gfmate(Int_t imat, char *name, Double_t &a, Double_t &z,  
266                     Double_t &dens, Double_t &radl, Double_t &absl,
267                     Double_t* ubuf, Int_t& nbuf) {
268 //
269   fGeometryManager->Gfmate(imat, name, a, z, dens, radl, absl, ubuf, nbuf);
270
271
272 // detector composition
273 void TFluka::Material(Int_t& kmat, const char* name, Double_t a, 
274                       Double_t z, Double_t dens, Double_t radl, Double_t absl,
275                       Float_t* buf, Int_t nwbuf) {
276 //
277   fGeometryManager
278     ->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf); 
279
280 void TFluka::Material(Int_t& kmat, const char* name, Double_t a, 
281                       Double_t z, Double_t dens, Double_t radl, Double_t absl,
282                       Double_t* buf, Int_t nwbuf) {
283 //
284   fGeometryManager
285     ->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf); 
286
287
288 void TFluka::Mixture(Int_t& kmat, const char *name, Float_t *a, 
289                      Float_t *z, Double_t dens, Int_t nlmat, Float_t *wmat) {
290 //
291    fGeometryManager
292      ->Mixture(kmat, name, a, z, dens, nlmat, wmat); 
293
294 void TFluka::Mixture(Int_t& kmat, const char *name, Double_t *a, 
295                      Double_t *z, Double_t dens, Int_t nlmat, Double_t *wmat) {
296 //
297    fGeometryManager
298      ->Mixture(kmat, name, a, z, dens, nlmat, wmat); 
299
300
301 void TFluka::Medium(Int_t& kmed, const char *name, Int_t nmat, 
302                     Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd, 
303                     Double_t stemax, Double_t deemax, Double_t epsil, 
304                     Double_t stmin, Float_t* ubuf, Int_t nbuf) { 
305   //
306   fGeometryManager
307     ->Medium(kmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax, 
308              epsil, stmin, ubuf, nbuf);
309
310 void TFluka::Medium(Int_t& kmed, const char *name, Int_t nmat, 
311                     Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd, 
312                     Double_t stemax, Double_t deemax, Double_t epsil, 
313                     Double_t stmin, Double_t* ubuf, Int_t nbuf) { 
314   //
315   fGeometryManager
316     ->Medium(kmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax, 
317              epsil, stmin, ubuf, nbuf);
318
319
320 void TFluka::Matrix(Int_t& krot, Double_t thetaX, Double_t phiX, 
321                     Double_t thetaY, Double_t phiY, Double_t thetaZ, 
322                     Double_t phiZ) {
323 //                   
324   fGeometryManager
325     ->Matrix(krot, thetaX, phiX, thetaY, phiY, thetaZ, phiZ); 
326
327
328 void TFluka::Gstpar(Int_t itmed, const char *param, Double_t parval) {
329 //
330   fGeometryManager->Gstpar(itmed, param, parval); 
331 }    
332
333 // functions from GGEOM 
334 Int_t TFluka::Gsvolu(const char *name, const char *shape, Int_t nmed,  
335                      Float_t *upar, Int_t np)  {
336 //
337 //  fVolumeMediaMap[TString(name)] = nmed;
338     if (fVerbosityLevel >= 3)
339     printf("TFluka::Gsvolu() name = %s, nmed = %d\n", name, nmed);
340     
341     TClonesArray &lvols = *fVolumeMediaMap;
342     new(lvols[fNVolumes++]) 
343         FlukaVolume(name, nmed);
344     return fGeometryManager->Gsvolu(name, shape, nmed, upar, np); 
345 }
346 Int_t TFluka::Gsvolu(const char *name, const char *shape, Int_t nmed,  
347                      Double_t *upar, Int_t np)  {
348 //
349     TClonesArray &lvols = *fVolumeMediaMap;
350     new(lvols[fNVolumes++]) 
351         FlukaVolume(name, nmed);
352
353     return fGeometryManager->Gsvolu(name, shape, nmed, upar, np); 
354 }
355  
356 void TFluka::Gsdvn(const char *name, const char *mother, Int_t ndiv, 
357                    Int_t iaxis) {
358 //
359 //  The medium of the daughter is the one of the mother
360     Int_t volid = TFluka::VolId(mother);
361     Int_t med   = TFluka::VolId2Mate(volid);
362     TClonesArray &lvols = *fVolumeMediaMap;
363     new(lvols[fNVolumes++]) 
364         FlukaVolume(name, med);
365     fGeometryManager->Gsdvn(name, mother, ndiv, iaxis); 
366
367
368 void TFluka::Gsdvn2(const char *name, const char *mother, Int_t ndiv, 
369                     Int_t iaxis, Double_t c0i, Int_t numed) {
370 //
371     TClonesArray &lvols = *fVolumeMediaMap;
372     new(lvols[fNVolumes++]) 
373         FlukaVolume(name, numed);
374     fGeometryManager->Gsdvn2(name, mother, ndiv, iaxis, c0i, numed); 
375
376
377 void TFluka::Gsdvt(const char *name, const char *mother, Double_t step, 
378                    Int_t iaxis, Int_t numed, Int_t ndvmx) {
379 //      
380     TClonesArray &lvols = *fVolumeMediaMap;
381     new(lvols[fNVolumes++]) 
382         FlukaVolume(name, numed);               
383     fGeometryManager->Gsdvt(name, mother, step, iaxis, numed, ndvmx); 
384
385
386 void TFluka::Gsdvt2(const char *name, const char *mother, Double_t step, 
387                     Int_t iaxis, Double_t c0, Int_t numed, Int_t ndvmx) { 
388 //
389     TClonesArray &lvols = *fVolumeMediaMap;
390     new(lvols[fNVolumes++]) 
391         FlukaVolume(name, numed);
392     fGeometryManager->Gsdvt2(name, mother, step, iaxis, c0, numed, ndvmx); 
393
394
395 void TFluka::Gsord(const char *name, Int_t iax) {
396 //
397   fGeometryManager->Gsord(name, iax); 
398
399
400 void TFluka::Gspos(const char *name, Int_t nr, const char *mother,  
401                    Double_t x, Double_t y, Double_t z, Int_t irot, 
402                    const char *konly) {
403 //
404   fGeometryManager->Gspos(name, nr, mother, x, y, z, irot, konly); 
405
406
407 void TFluka::Gsposp(const char *name, Int_t nr, const char *mother,  
408                     Double_t x, Double_t y, Double_t z, Int_t irot,
409                     const char *konly, Float_t *upar, Int_t np)  {
410   //
411   fGeometryManager->Gsposp(name, nr, mother, x, y, z, irot, konly, upar, np); 
412
413 void TFluka::Gsposp(const char *name, Int_t nr, const char *mother,  
414                     Double_t x, Double_t y, Double_t z, Int_t irot,
415                     const char *konly, Double_t *upar, Int_t np)  {
416   //
417   fGeometryManager->Gsposp(name, nr, mother, x, y, z, irot, konly, upar, np); 
418
419
420 void TFluka::Gsbool(const char* onlyVolName, const char* manyVolName) {
421 //
422   fGeometryManager->Gsbool(onlyVolName, manyVolName);
423 }
424
425 void TFluka::SetCerenkov(Int_t itmed, Int_t npckov, Float_t *ppckov,
426                          Float_t *absco, Float_t *effic, Float_t *rindex) {
427 //
428   fGeometryManager->SetCerenkov(itmed, npckov, ppckov, absco, effic, rindex);
429 }  
430 void TFluka::SetCerenkov(Int_t itmed, Int_t npckov, Double_t *ppckov,
431                          Double_t *absco, Double_t *effic, Double_t *rindex) {
432 //
433   fGeometryManager->SetCerenkov(itmed, npckov, ppckov, absco, effic, rindex);
434 }  
435     
436 // Euclid
437 void TFluka::WriteEuclid(const char* fileName, const char* topVol, 
438                           Int_t number, Int_t nlevel) {
439 //
440   fGeometryManager->WriteEuclid(fileName, topVol, number, nlevel); 
441
442
443
444
445 //_____________________________________________________________________________
446 // methods needed by the stepping
447 //____________________________________________________________________________ 
448
449 Int_t TFluka::GetMedium() const {
450 //
451 //  Get the medium number for the current fluka region
452 //
453     FGeometryInit* flugg = FGeometryInit::GetInstance();  
454     return flugg->GetMedium(fCurrentFlukaRegion);
455 }
456
457
458
459 //____________________________________________________________________________ 
460 // particle table usage
461 // ID <--> PDG transformations
462 //_____________________________________________________________________________
463 Int_t TFluka::IdFromPDG(Int_t pdg) const 
464 {
465     //
466     // Return Fluka code from PDG and pseudo ENDF code
467     
468     // Catch the feedback photons
469     if (pdg == 50000051) return (-1);
470     // MCIHAD() goes from pdg to fluka internal.
471     Int_t intfluka = mcihad(pdg);
472     // KPTOIP array goes from internal to official
473     return GetFlukaKPTOIP(intfluka);
474 }
475
476 Int_t TFluka::PDGFromId(Int_t id) const 
477 {
478   //
479   // Return PDG code and pseudo ENDF code from Fluka code
480
481   // IPTOKP array goes from official to internal
482
483     if (id == -1) {
484 // Cerenkov photon
485         if (fVerbosityLevel >= 1)
486             printf("\n PDGFromId: Cerenkov Photon \n");
487         return  50000050;
488     }
489 // Error id    
490     if (id == 0) {
491         if (fVerbosityLevel >= 1)
492             printf("PDGFromId: Error id = 0\n");
493         return -1;
494     }
495 // Good id    
496     Int_t intfluka = GetFlukaIPTOKP(id);
497     if (intfluka == 0) {
498         if (fVerbosityLevel >= 1)
499             printf("PDGFromId: Error intfluka = 0: %d\n", id);
500         return -1;
501     } else if (intfluka < 0) {
502         if (fVerbosityLevel >= 1)
503             printf("PDGFromId: Error intfluka < 0: %d\n", id);
504         return -1;
505     }
506     if (fVerbosityLevel >= 3)
507         printf("mpdgha called with %d %d \n", id, intfluka);
508     // MPDGHA() goes from fluka internal to pdg.
509     return mpdgha(intfluka);
510 }
511
512 //_____________________________________________________________________________
513 // methods for physics management
514 //____________________________________________________________________________ 
515 //
516 // set methods
517 //
518
519 void TFluka::SetProcess(const char* flagName, Int_t flagValue)
520 {
521   Int_t i;
522   if (iNbOfProc < 100) {
523     for (i=0; i<iNbOfProc; i++) {
524       if (strcmp(&sProcessFlag[i][0],flagName) == 0) {
525         iProcessValue[iNbOfProc] = flagValue;
526         goto fin;
527       }
528     }
529     strcpy(&sProcessFlag[iNbOfProc][0],flagName);
530     iProcessValue[iNbOfProc++] = flagValue;
531   }
532   else
533     cout << "Nb of SetProcess calls exceeds 100 - ignored" << endl;
534 fin:
535   iNbOfProc = iNbOfProc;
536 }
537
538 void TFluka::SetCut(const char* cutName, Double_t cutValue)
539 {
540   Int_t i;
541   if (iNbOfCut < 100) {
542     for (i=0; i<iNbOfCut; i++) {
543       if (strcmp(&sCutFlag[i][0],cutName) == 0) {
544         fCutValue[iNbOfCut] = cutValue;
545         goto fin;
546       }
547     }
548     strcpy(&sCutFlag[iNbOfCut][0],cutName);
549     fCutValue[iNbOfCut++] = cutValue;
550   }
551   else
552     cout << "Nb of SetCut calls exceeds 100 - ignored" << endl;
553 fin:
554   iNbOfCut = iNbOfCut;
555 }
556
557 Double_t TFluka::Xsec(char*, Double_t, Int_t, Int_t)
558 {
559   printf("WARNING: Xsec not yet implemented !\n"); return -1.;
560 }
561
562
563 void TFluka::InitPhysics()
564 {
565 // Last material number taken from the "corealice.inp" file, presently 31
566 // !!! it should be available from Flugg !!!
567   Int_t i, j, k;
568   Double_t fCut;
569   FGeometryInit* geominit = FGeometryInit::GetInstance();
570   Float_t fLastMaterial = geominit->GetLastMaterialIndex();
571   printf("   last FLUKA material is %g\n", fLastMaterial);
572  
573 // construct file names
574   TString sAliceCoreInp = getenv("ALICE_ROOT");
575   sAliceCoreInp +="/TFluka/input/";
576   TString sAliceTmp = "flukaMat.inp";
577   TString sAliceInp = GetInputFileName();
578   sAliceCoreInp += GetCoreInputFileName();
579   ifstream AliceCoreInp(sAliceCoreInp.Data());
580   ifstream AliceFlukaMat(sAliceTmp.Data());
581   ofstream AliceInp(sAliceInp.Data());
582
583 // copy core input file 
584   Char_t sLine[255];
585   Float_t fEventsPerRun;
586
587   while (AliceCoreInp.getline(sLine,255)) {
588     if (strncmp(sLine,"GEOEND",6) != 0)
589       AliceInp << sLine << endl; // copy until GEOEND card
590     else {
591       AliceInp << "GEOEND" << endl; // add GEOEND card
592       goto flukamat;
593     }
594   } // end of while until GEOEND card
595
596 flukamat:
597   while (AliceFlukaMat.getline(sLine,255)) { // copy flukaMat.inp file
598         AliceInp << sLine << endl;
599   }
600
601   while (AliceCoreInp.getline(sLine,255)) {
602     if (strncmp(sLine,"START",5) != 0)
603       AliceInp << sLine << endl;
604     else {
605       sscanf(sLine+10,"%10f",&fEventsPerRun);
606       goto fin;
607     }
608   } //end of while until START card
609
610 fin:
611 // in G3 the process control values meaning can be different for
612 // different processes, but for most of them is:
613 //   0  process is not activated
614 //   1  process is activated WITH generation of secondaries
615 //   2  process is activated WITHOUT generation of secondaries
616 // if process does not generate secondaries => 1 same as 2
617 //
618 // Exceptions:
619 //   MULS:  also 3
620 //   LOSS:  also 3, 4
621 //   RAYL:  only 0,1
622 //   HADR:  may be > 2
623 //
624  
625 // Loop over number of SetProcess calls  
626   AliceInp << "*----------------------------------------------------------------------------- "; 
627   AliceInp << endl;
628   AliceInp << "*----- The following data are generated from SetProcess and SetCut calls ----- "; 
629   AliceInp << endl;
630   AliceInp << "*----------------------------------------------------------------------------- "; 
631     AliceInp << endl;
632   for (i=0; i<iNbOfProc; i++) {
633
634     // annihilation
635     // G3 default value: 1
636     // G4 processes: G4eplusAnnihilation/G4IeplusAnnihilation
637     // Particles: e+
638     // Physics:   EM
639     // flag = 0 no annihilation
640     // flag = 1 annihilation, decays processed
641     // flag = 2 annihilation, no decay product stored
642     // gMC ->SetProcess("ANNI",1); // EMFCUT   -1.   0.  0. 3. lastmat 0. ANNH-THR
643     if (strncmp(&sProcessFlag[i][0],"ANNI",4) == 0) {
644       if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
645         AliceInp << "*"; 
646         AliceInp << endl;
647         AliceInp << "*Kinetic energy threshold (GeV) for e+ annihilation - resets to default=0."; 
648         AliceInp << endl;
649         AliceInp << "*Generated from call: SetProcess('ANNI',1) or SetProcess('ANNI',2)"; 
650         AliceInp << endl;
651         AliceInp << setw(10) << "EMFCUT    "; 
652         AliceInp << setiosflags(ios::scientific) << setprecision(5);
653         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
654         AliceInp << setw(10) << -1.0; // kinetic energy threshold (GeV) for e+ annihilation (resets to default=0)
655         AliceInp << setw(10) << 0.0;  // not used
656         AliceInp << setw(10) << 0.0;  // not used
657         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
658         AliceInp << setw(10) << setprecision(2);
659         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
660         AliceInp << setprecision(1);
661         AliceInp << setw(10) << 1.0;  // step length in assigning indices
662         AliceInp << setw(8)  << "ANNH-THR"; 
663         AliceInp << endl;
664       }
665       else if (iProcessValue[i] == 0) {
666         AliceInp << "*"; 
667         AliceInp << endl;
668         AliceInp << "*No annihilation - no FLUKA card generated"; 
669         AliceInp << endl;
670         AliceInp << "*Generated from call: SetProcess('ANNI',0)"; 
671         AliceInp << endl;
672       }
673       else  {
674         AliceInp << "*"; 
675         AliceInp << endl;
676         AliceInp << "*Illegal flag value in SetProcess('ANNI',?) call."; 
677         AliceInp << endl;
678         AliceInp << "*No FLUKA card generated"; 
679         AliceInp << endl;
680       }
681     }
682     
683     // bremsstrahlung and pair production are both activated
684     // G3 default value: 1
685     // G4 processes: G4eBremsstrahlung/G4IeBremsstrahlung,
686     //               G4MuBremsstrahlung/G4IMuBremsstrahlung,
687     //               G4LowEnergyBremstrahlung
688     // Particles: e-/e+; mu+/mu-
689     // Physics:   EM
690     // flag = 0 no bremsstrahlung
691     // flag = 1 bremsstrahlung, photon processed
692     // flag = 2 bremsstrahlung, no photon stored
693     // gMC ->SetProcess("BREM",1); // PAIRBREM  2.   0.  0. 3. lastmat
694                                  // EMFCUT   -1.   0.  0. 3. lastmat 0. ELPO-THR
695     // G3 default value: 1
696     // G4 processes: G4GammaConversion,
697     //               G4MuPairProduction/G4IMuPairProduction
698     //               G4LowEnergyGammaConversion
699     // Particles: gamma, mu
700     // Physics:   EM
701     // flag = 0 no delta rays
702     // flag = 1 delta rays, secondaries processed
703     // flag = 2 delta rays, no secondaries stored
704     // gMC ->SetProcess("PAIR",1); // PAIRBREM  1.   0.  0. 3. lastmat
705                                  // EMFCUT    0.   0. -1. 3. lastmat 0. PHOT-THR
706     else if ((strncmp(&sProcessFlag[i][0],"PAIR",4) == 0) && (iProcessValue[i] == 1 || iProcessValue[i] == 2)) {
707       for (j=0; j<iNbOfProc; j++) {
708         if ((strncmp(&sProcessFlag[j][0],"BREM",4) == 0) && (iProcessValue[j] == 1 || iProcessValue[j] == 2)) {
709           AliceInp << "*"; 
710           AliceInp << endl;
711           AliceInp << "*Bremsstrahlung and pair production by muons and charged hadrons both activated"; 
712           AliceInp << endl;
713           AliceInp << "*Generated from call: SetProcess('BREM',1) and SetProcess('PAIR',1)"; 
714           AliceInp << endl;
715           AliceInp << "*Energy threshold set by call SetCut('BCUTM',cut) or set to 0."; 
716           AliceInp << endl;
717           AliceInp << "*Energy threshold set by call SetCut('PPCUTM',cut) or set to 0."; 
718           AliceInp << endl;
719           AliceInp << setw(10) << "PAIRBREM  "; 
720           AliceInp << setiosflags(ios::scientific) << setprecision(5);
721           AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
722           AliceInp << setw(10) << 3.0; // bremsstrahlung and pair production by muons and charged hadrons both are activated
723           // direct pair production by muons
724           // G4 particles: "e-", "e+"
725           // G3 default value: 0.01 GeV
726           //gMC ->SetCut("PPCUTM",cut); // total energy cut for direct pair prod. by muons
727           fCut = 0.0;
728           for (k=0; k<iNbOfCut; k++) {
729             if (strncmp(&sCutFlag[k][0],"PPCUTM",6) == 0) fCut = fCutValue[k];
730           }
731           AliceInp << setiosflags(ios::scientific) << setprecision(5);
732           AliceInp << setw(10) << fCut; // e+, e- kinetic energy threshold (in GeV) for explicit pair production.
733           // muon and hadron bremsstrahlung
734           // G4 particles: "gamma"
735           // G3 default value: CUTGAM=0.001 GeV
736           //gMC ->SetCut("BCUTM",cut);  // cut for muon and hadron bremsstrahlung
737           fCut = 0.0;
738           for (k=0; k<iNbOfCut; k++) {
739             if (strncmp(&sCutFlag[k][0],"BCUTM",5) == 0) fCut = fCutValue[k];
740           }
741           AliceInp << setiosflags(ios::scientific) << setprecision(5);
742           AliceInp << setw(10) << fCut; // photon energy threshold (GeV) for explicit bremsstrahlung production
743           AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
744           AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
745           AliceInp << setw(10) << setprecision(2);
746           AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
747           AliceInp << endl;
748           
749           // for e+ and e-
750           AliceInp << "*"; 
751           AliceInp << endl;
752           AliceInp << "*Kinetic energy threshold (GeV) for e+/e- bremsstrahlung - resets to default=0."; 
753           AliceInp << endl;
754           AliceInp << "*Generated from call: SetProcess('BREM',1);"; 
755           AliceInp << endl;
756           AliceInp << setw(10) << "EMFCUT    "; 
757           fCut = -1.0;
758           for (k=0; k<iNbOfCut; k++) {
759             if (strncmp(&sCutFlag[k][0],"BCUTE",5) == 0) fCut = fCutValue[k];
760           }
761           AliceInp << setiosflags(ios::scientific) << setprecision(5);
762           AliceInp << setw(10) << fCut; // kinetic energy threshold (GeV) for e+/e- bremsstrahlung (resets to default=0)
763           AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
764           AliceInp << setw(10) << 0.0;  // not used
765           AliceInp << setw(10) << 0.0;  // not used
766           AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
767           AliceInp << setw(10) << setprecision(2);
768           AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
769           AliceInp << setprecision(1);
770           AliceInp << setw(10) << 1.0; // step length in assigning indices
771           AliceInp << setw(8)  << "ELPO-THR"; 
772           AliceInp << endl;
773       
774           // for e+ and e-
775           AliceInp << "*"; 
776           AliceInp << endl;
777           AliceInp << "*Pair production by electrons is activated"; 
778           AliceInp << endl;
779           AliceInp << "*Generated from call: SetProcess('PAIR',1);"; 
780           AliceInp << endl;
781           AliceInp << setw(10) << "EMFCUT    "; 
782           AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
783           AliceInp << setw(10) << 0.0;  // energy threshold (GeV) for Compton scattering (= 0.0 : ignored)
784           AliceInp << setw(10) << 0.0;  // energy threshold (GeV) for Photoelectric (= 0.0 : ignored)
785           fCut = -1.0;
786           for (j=0; j<iNbOfCut; j++) {
787             if (strncmp(&sCutFlag[j][0],"CUTGAM",6) == 0) fCut = fCutValue[j];
788           }
789           AliceInp << setiosflags(ios::scientific) << setprecision(5);
790           AliceInp << setw(10) << fCut; // energy threshold (GeV) for gamma pair production (< 0.0 : resets to default, = 0.0 : ignored)
791           AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
792           AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
793           AliceInp << setprecision(2);
794           AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
795           AliceInp << setprecision(1);
796           AliceInp << setw(10) << 1.0;  // step length in assigning indices
797           AliceInp << setw(8) << "PHOT-THR"; 
798           AliceInp << endl;
799           goto BOTH;
800         } // end of if for BREM
801       } // end of loop for BREM
802
803       // only pair production by muons and charged hadrons is activated
804       AliceInp << "*"; 
805       AliceInp << endl;
806       AliceInp << "*Pair production by muons and charged hadrons is activated"; 
807       AliceInp << endl;
808       AliceInp << "*Generated from call: SetProcess('PAIR',1) or SetProcess('PAIR',2)"; 
809       AliceInp << endl;
810       AliceInp << "*Energy threshold set by call SetCut('PPCUTM',cut) or set to 0."; 
811       AliceInp << endl;
812       AliceInp << setw(10) << "PAIRBREM  "; 
813       AliceInp << setiosflags(ios::scientific) << setprecision(5);
814       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
815       AliceInp << setw(10) << 1.0; // pair production by muons and charged hadrons is activated
816       // direct pair production by muons
817       // G4 particles: "e-", "e+"
818       // G3 default value: 0.01 GeV
819       //gMC ->SetCut("PPCUTM",cut); // total energy cut for direct pair prod. by muons
820       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
821       AliceInp << setw(10) << 0.0; // e+, e- kinetic energy threshold (in GeV) for explicit pair production.
822       AliceInp << setw(10) << 0.0; // no explicit bremsstrahlung production is simulated
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 << endl;
827       
828       // for e+ and e-
829       AliceInp << "*"; 
830       AliceInp << endl;
831       AliceInp << "*Pair production by electrons is activated"; 
832       AliceInp << endl;
833       AliceInp << "*Generated from call: SetProcess('PAIR',1) or SetProcess('PAIR',2)"; 
834       AliceInp << endl;
835       AliceInp << setw(10) << "EMFCUT    "; 
836       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
837       AliceInp << setw(10) << 0.0;  // energy threshold (GeV) for Compton scattering (= 0.0 : ignored)
838       AliceInp << setw(10) << 0.0;  // energy threshold (GeV) for Photoelectric (= 0.0 : ignored)
839
840       fCut = -1.0;
841       for (j=0; j<iNbOfCut; j++) {
842         if (strncmp(&sCutFlag[j][0],"CUTGAM",6) == 0) fCut = fCutValue[j];
843       }
844       AliceInp << setiosflags(ios::scientific) << setprecision(5);
845       AliceInp << setw(10) << fCut; // energy threshold (GeV) for gamma pair production (< 0.0 : resets to default, = 0.0 : ignored)
846       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
847       AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
848       AliceInp << setprecision(2);
849       AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
850       AliceInp << setprecision(1);
851       AliceInp << setw(10) << 1.0;  // step length in assigning indices
852       AliceInp << setw(8) << "PHOT-THR"; 
853       AliceInp << endl;
854
855 BOTH:
856     k = 0;
857     } // end of if for PAIR
858
859
860
861     // bremsstrahlung
862     // G3 default value: 1
863     // G4 processes: G4eBremsstrahlung/G4IeBremsstrahlung,
864     //               G4MuBremsstrahlung/G4IMuBremsstrahlung,
865     //               G4LowEnergyBremstrahlung
866     // Particles: e-/e+; mu+/mu-
867     // Physics:   EM
868     // flag = 0 no bremsstrahlung
869     // flag = 1 bremsstrahlung, photon processed
870     // flag = 2 bremsstrahlung, no photon stored
871     // gMC ->SetProcess("BREM",1); // PAIRBREM  2.   0.  0. 3. lastmat
872                                  // EMFCUT   -1.   0.  0. 3. lastmat 0. ELPO-THR
873     else if (strncmp(&sProcessFlag[i][0],"BREM",4) == 0) {
874       for (j=0; j<iNbOfProc; j++) {
875         if ((strncmp(&sProcessFlag[j][0],"PAIR",4) == 0) && iProcessValue[j] == 1) goto NOBREM;
876       }
877       if (iProcessValue[i] == 1 || iProcessValue[i] == 2) { 
878         AliceInp << "*"; 
879         AliceInp << endl;
880         AliceInp << "*Bremsstrahlung by muons and charged hadrons is activated"; 
881         AliceInp << endl;
882         AliceInp << "*Generated from call: SetProcess('BREM',1) or SetProcess('BREM',2)"; 
883         AliceInp << endl;
884         AliceInp << "*Energy threshold set by call SetCut('BCUTM',cut) or set to 0."; 
885         AliceInp << endl;
886         AliceInp << setw(10) << "PAIRBREM  "; 
887         AliceInp << setiosflags(ios::scientific) << setprecision(5);
888         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
889         AliceInp << setw(10) << 2.0; // bremsstrahlung by muons and charged hadrons is activated
890         AliceInp << setw(10) << 0.0; // no meaning
891         // muon and hadron bremsstrahlung
892         // G4 particles: "gamma"
893         // G3 default value: CUTGAM=0.001 GeV
894         //gMC ->SetCut("BCUTM",cut);  // cut for muon and hadron bremsstrahlung
895         fCut = 0.0;
896         for (j=0; j<iNbOfCut; j++) {
897           if (strncmp(&sCutFlag[j][0],"BCUTM",5) == 0) fCut = fCutValue[j];
898         }
899         AliceInp << setw(10) << fCut; // photon energy threshold (GeV) for explicit bremsstrahlung production
900         AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
901         AliceInp << setw(10) << setprecision(2);
902         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
903         AliceInp << endl;
904  
905         // for e+ and e-
906         AliceInp << "*"; 
907         AliceInp << endl;
908         AliceInp << "*Kinetic energy threshold (GeV) for e+/e- bremsstrahlung - resets to default=0."; 
909         AliceInp << endl;
910         AliceInp << "*Generated from call: SetProcess('BREM',1);"; 
911         AliceInp << endl;
912         AliceInp << setw(10) << "EMFCUT    "; 
913         AliceInp << setiosflags(ios::scientific) << setprecision(5);
914         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
915         AliceInp << setw(10) << -1.0; // kinetic energy threshold (GeV) for e+/e- bremsstrahlung (resets to default=0)
916         AliceInp << setw(10) << 0.0;  // not used
917         AliceInp << setw(10) << 0.0;  // not used
918         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
919         AliceInp << setw(10) << setprecision(2);
920         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
921         AliceInp << setprecision(1);
922         AliceInp << setw(10) << 1.0; // step length in assigning indices
923         AliceInp << setw(8)  << "ELPO-THR"; 
924         AliceInp << endl;
925       }
926       else if (iProcessValue[i] == 0) {
927         AliceInp << "*"; 
928         AliceInp << endl;
929         AliceInp << "*No bremsstrahlung - no FLUKA card generated"; 
930         AliceInp << endl;
931         AliceInp << "*Generated from call: SetProcess('BREM',0)"; 
932         AliceInp << endl;
933       }
934       else  {
935         AliceInp << "*"; 
936         AliceInp << endl;
937         AliceInp << "*Illegal flag value in SetProcess('BREM',?) call."; 
938         AliceInp << endl;
939         AliceInp << "*No FLUKA card generated"; 
940         AliceInp << endl;
941       }
942 NOBREM:
943       j = 0;
944     } // end of else if (strncmp(&sProcessFlag[i][0],"BREM",4) == 0)
945
946     
947     // Cerenkov photon generation
948     // G3 default value: 0
949     // G4 process: G4Cerenkov
950     // 
951     // Particles: charged
952     // Physics:   Optical
953     // flag = 0 no Cerenkov photon generation
954     // flag = 1 Cerenkov photon generation
955     // flag = 2 Cerenkov photon generation with primary stopped at each step
956     //xx gMC ->SetProcess("CKOV",1); // ??? Cerenkov photon generation
957     else if (strncmp(&sProcessFlag[i][0],"CKOV",4) == 0) {
958       if (iProcessValue[i] == 1 || iProcessValue[i] == 2) { 
959         AliceInp << "*"; 
960         AliceInp << endl;
961         AliceInp << "*Cerenkov photon generation"; 
962         AliceInp << endl;
963         AliceInp << "*Generated from call: SetProcess('CKOV',1) or SetProcess('CKOV',2)"; 
964         AliceInp << endl;
965         AliceInp << setw(10) << "OPT-PROD  "; 
966         AliceInp << setiosflags(ios::scientific) << setprecision(5);
967         AliceInp << setw(10) <<  2.07e-9 ; //  minimum Cerenkov photon emission energy (in GeV!). Default: 2.07E-9 GeV (corresponding to 600 nm)
968         AliceInp << setw(10) << 4.96e-9;  // maximum Cerenkov photon emission energy (in GeV!). Default: 4.96E-9 GeV (corresponding to 250 nm)
969         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
970         AliceInp << setw(10) << 0.0;  // not used
971         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
972         AliceInp << setprecision(2);
973         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
974         AliceInp << setprecision(1);
975         AliceInp << setw(10) << 1.0; // step length in assigning indices
976         AliceInp << setw(8) << "CERENKOV"; 
977         AliceInp << endl;
978       }
979       else if (iProcessValue[i] == 0) {
980         AliceInp << "*"; 
981         AliceInp << endl;
982         AliceInp << "*No Cerenkov photon generation"; 
983         AliceInp << endl;
984         AliceInp << "*Generated from call: SetProcess('CKOV',0)"; 
985         AliceInp << endl;
986         AliceInp << setw(10) << "OPT-PROD  "; 
987         AliceInp << setiosflags(ios::scientific) << setprecision(5);
988         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
989         AliceInp << setw(10) << 0.0;  // not used
990         AliceInp << setw(10) << 0.0;  // not used
991         AliceInp << setw(10) << 0.0;  // not used
992         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
993         AliceInp << setprecision(2);
994         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
995         AliceInp << setprecision(1);
996         AliceInp << setw(10) << 1.0; // step length in assigning indices
997         AliceInp << setw(8) << "CERE-OFF"; 
998         AliceInp << endl;
999       }
1000       else  {
1001         AliceInp << "*"; 
1002         AliceInp << endl;
1003         AliceInp << "*Illegal flag value in SetProcess('CKOV',?) call."; 
1004         AliceInp << endl;
1005         AliceInp << "*No FLUKA card generated"; 
1006         AliceInp << endl;
1007       }
1008     } // end of else if (strncmp(&sProcessFlag[i][0],"CKOV",4) == 0)
1009             
1010
1011     // Compton scattering
1012     // G3 default value: 1
1013     // G4 processes: G4ComptonScattering,
1014     //               G4LowEnergyCompton,
1015     //               G4PolarizedComptonScattering
1016     // Particles: gamma
1017     // Physics:   EM
1018     // flag = 0 no Compton scattering
1019     // flag = 1 Compton scattering, electron processed
1020     // flag = 2 Compton scattering, no electron stored
1021     // gMC ->SetProcess("COMP",1); // EMFCUT   -1.   0.  0. 3. lastmat 0. PHOT-THR
1022     else if (strncmp(&sProcessFlag[i][0],"COMP",4) == 0) {
1023       if (iProcessValue[i] == 1 || iProcessValue[i] == 2) { 
1024         AliceInp << "*"; 
1025         AliceInp << endl;
1026         AliceInp << "*Energy threshold (GeV) for Compton scattering - resets to default=0."; 
1027         AliceInp << endl;
1028         AliceInp << "*Generated from call: SetProcess('COMP',1);"; 
1029         AliceInp << endl;
1030         AliceInp << setw(10) << "EMFCUT    "; 
1031         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1032         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1033         AliceInp << setw(10) << -1.0; // energy threshold (GeV) for Compton scattering - resets to default=0.
1034         AliceInp << setw(10) << 0.0;  // not used
1035         AliceInp << setw(10) << 0.0;  // not used
1036         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
1037         AliceInp << setprecision(2);
1038         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1039         AliceInp << setprecision(1);
1040         AliceInp << setw(10) << 1.0; // step length in assigning indices
1041         AliceInp << setw(8) << "PHOT-THR"; 
1042         AliceInp << endl;
1043       }
1044       else if (iProcessValue[i] == 0) {
1045         AliceInp << "*"; 
1046         AliceInp << endl;
1047         AliceInp << "*No Compton scattering - no FLUKA card generated"; 
1048         AliceInp << endl;
1049         AliceInp << "*Generated from call: SetProcess('COMP',0)"; 
1050         AliceInp << endl;
1051       }
1052       else  {
1053         AliceInp << "*"; 
1054         AliceInp << endl;
1055         AliceInp << "*Illegal flag value in SetProcess('COMP',?) call."; 
1056         AliceInp << endl;
1057         AliceInp << "*No FLUKA card generated"; 
1058         AliceInp << endl;
1059       }
1060     } // end of else if (strncmp(&sProcessFlag[i][0],"COMP",4) == 0)
1061
1062     // decay
1063     // G3 default value: 1
1064     // G4 process: G4Decay
1065     // 
1066     // Particles: all which decay is applicable for
1067     // Physics:   General
1068     // flag = 0 no decays
1069     // flag = 1 decays, secondaries processed
1070     // flag = 2 decays, no secondaries stored
1071     //gMC ->SetProcess("DCAY",1); // not available
1072     else if ((strncmp(&sProcessFlag[i][0],"DCAY",4) == 0) && iProcessValue[i] == 1) 
1073       cout << "SetProcess for flag=" << &sProcessFlag[i][0] << " value=" << iProcessValue[i] << " not avaliable!" << endl;
1074       
1075     // delta-ray
1076     // G3 default value: 2
1077     // !! G4 treats delta rays in different way
1078     // G4 processes: G4eIonisation/G4IeIonization,
1079     //               G4MuIonisation/G4IMuIonization,
1080     //               G4hIonisation/G4IhIonisation
1081     // Particles: charged
1082     // Physics:   EM
1083     // flag = 0 no energy loss
1084     // flag = 1 restricted energy loss fluctuations
1085     // flag = 2 complete energy loss fluctuations
1086     // flag = 3 same as 1
1087     // flag = 4 no energy loss fluctuations
1088     // gMC ->SetProcess("DRAY",0); // DELTARAY 1.E+6 0.  0. 3. lastmat 0.
1089     else if (strncmp(&sProcessFlag[i][0],"DRAY",4) == 0) {
1090       if (iProcessValue[i] == 0 || iProcessValue[i] == 4) {
1091         AliceInp << "*"; 
1092         AliceInp << endl;
1093         AliceInp << "*Kinetic energy threshold (GeV) for delta ray production"; 
1094         AliceInp << endl;
1095         AliceInp << "*Generated from call: SetProcess('DRAY',0) or SetProcess('DRAY',4)"; 
1096         AliceInp << endl;
1097         AliceInp << "*No delta ray production by muons - threshold set artificially high"; 
1098         AliceInp << endl;
1099         AliceInp << setw(10) << "DELTARAY  "; 
1100         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1101         AliceInp << setw(10) << 1.0e+6; // kinetic energy threshold (GeV) for delta ray production (discrete energy transfer)
1102         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1103         AliceInp << setw(10) << 0.0; // ignored
1104         AliceInp << setw(10) << 0.0; // ignored
1105         AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1106         AliceInp << setw(10) << setprecision(2);
1107         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1108         AliceInp << setprecision(1);
1109         AliceInp << setw(10) << 1.0; // step length in assigning indices
1110         AliceInp << endl;
1111       }
1112       else if (iProcessValue[i] == 1 || iProcessValue[i] == 2 || iProcessValue[i] == 3) {
1113         AliceInp << "*"; 
1114         AliceInp << endl;
1115         AliceInp << "*Kinetic energy threshold (GeV) for delta ray production"; 
1116         AliceInp << endl;
1117         AliceInp << "*Generated from call: SetProcess('DRAY',flag), flag=1,2,3"; 
1118         AliceInp << endl;
1119         AliceInp << "*Delta ray production by muons switched on"; 
1120         AliceInp << endl;
1121         AliceInp << "*Energy threshold set by call SetCut('DCUTM',cut) or set to 0."; 
1122         AliceInp << endl;
1123         AliceInp << setw(10) << "DELTARAY  "; 
1124         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1125         fCut = 1.0e+6;
1126         for (j=0; j<iNbOfCut; j++) {
1127           if (strncmp(&sCutFlag[j][0],"DCUTM",5) == 0) fCut = fCutValue[j];
1128         }
1129         AliceInp << setw(10) << fCut; // kinetic energy threshold (GeV) for delta ray production (discrete energy transfer)
1130         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1131         AliceInp << setw(10) << 0.0; // ignored
1132         AliceInp << setw(10) << 0.0; // ignored
1133         AliceInp << setw(10) << 3.0; // lower bound of the material indices in which the respective thresholds apply
1134         AliceInp << setw(10) << setprecision(2);
1135         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1136         AliceInp << setprecision(1);
1137         AliceInp << setw(10) << 1.0; // step length in assigning indices
1138         AliceInp << endl;
1139       }
1140       else  {
1141         AliceInp << "*"; 
1142         AliceInp << endl;
1143         AliceInp << "*Illegal flag value in SetProcess('DRAY',?) call."; 
1144         AliceInp << endl;
1145         AliceInp << "*No FLUKA card generated"; 
1146         AliceInp << endl;
1147       }
1148     } // end of else if (strncmp(&sProcessFlag[i][0],"DRAY",4) == 0)
1149     
1150     // hadronic process
1151     // G3 default value: 1
1152     // G4 processes: all defined by TG4PhysicsConstructorHadron
1153     //  
1154     // Particles: hadrons
1155     // Physics:   Hadron
1156     // flag = 0 no multiple scattering
1157     // flag = 1 hadronic interactions, secondaries processed
1158     // flag = 2 hadronic interactions, no secondaries stored
1159     // gMC ->SetProcess("HADR",1); // ??? hadronic process
1160     //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3) ?????
1161     else if (strncmp(&sProcessFlag[i][0],"HADR",4) == 0) {
1162       if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
1163         AliceInp << "*"; 
1164         AliceInp << endl;
1165         AliceInp << "*Hadronic interaction is ON by default in FLUKA"; 
1166         AliceInp << endl;
1167         AliceInp << "*No FLUKA card generated"; 
1168         AliceInp << endl;
1169       }
1170       else if (iProcessValue[i] == 0) {
1171         AliceInp << "*"; 
1172         AliceInp << endl;
1173         AliceInp << "*Hadronic interaction is set OFF"; 
1174         AliceInp << endl;
1175         AliceInp << "*Generated from call: SetProcess('HADR',0);"; 
1176         AliceInp << endl;
1177         AliceInp << setw(10) << "MULSOPT  "; 
1178         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1179         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1180         AliceInp << setw(10) << 0.0;  // ignored
1181         AliceInp << setw(10) << 3.0;  // multiple scattering for hadrons and muons is completely suppressed
1182         AliceInp << setw(10) << 0.0;  // no spin-relativistic corrections
1183         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
1184         AliceInp << setprecision(2);
1185         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1186         AliceInp << endl;
1187
1188       }
1189       else  {
1190         AliceInp << "*"; 
1191         AliceInp << endl;
1192         AliceInp << "*Illegal flag value in SetProcess('HADR',?) call."; 
1193         AliceInp << endl;
1194         AliceInp << "*No FLUKA card generated"; 
1195         AliceInp << endl;
1196       }
1197     } // end of else if (strncmp(&sProcessFlag[i][0],"HADR",4) == 0)
1198
1199
1200     // energy loss
1201     // G3 default value: 2
1202     // G4 processes: G4eIonisation/G4IeIonization,
1203     //               G4MuIonisation/G4IMuIonization,
1204     //               G4hIonisation/G4IhIonisation
1205     // 
1206     // Particles: charged
1207     // Physics:   EM
1208     // flag=0 no energy loss
1209     // flag=1 restricted energy loss fluctuations
1210     // flag=2 complete energy loss fluctuations
1211     // flag=3 same as 1
1212     // flag=4 no energy loss fluctuations
1213     // If the value ILOSS is changed, then (in G3) cross-sections and energy
1214     // loss tables must be recomputed via the command 'PHYSI'
1215     // gMC ->SetProcess("LOSS",2); // ??? IONFLUCT ? energy loss
1216     else if (strncmp(&sProcessFlag[i][0],"LOSS",4) == 0) {
1217       if (iProcessValue[i] == 2) { // complete energy loss fluctuations
1218         AliceInp << "*"; 
1219         AliceInp << endl;
1220         AliceInp << "*Complete energy loss fluctuations do not exist in FLUKA"; 
1221         AliceInp << endl;
1222         AliceInp << "*Generated from call: SetProcess('LOSS',2);"; 
1223         AliceInp << endl;
1224         AliceInp << "*flag=2=complete energy loss fluctuations"; 
1225         AliceInp << endl;
1226         AliceInp << "*No input card generated"; 
1227         AliceInp << endl;
1228       }
1229       else if (iProcessValue[i] == 1 || iProcessValue[i] == 3) { // restricted energy loss fluctuations
1230         AliceInp << "*"; 
1231         AliceInp << endl;
1232         AliceInp << "*Restricted energy loss fluctuations"; 
1233         AliceInp << endl;
1234         AliceInp << "*Generated from call: SetProcess('LOSS',1) or SetProcess('LOSS',3)"; 
1235         AliceInp << endl;
1236         AliceInp << setw(10) << "IONFLUCT  "; 
1237         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1238         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1239         AliceInp << setw(10) << 1.0;  // restricted energy loss fluctuations (for hadrons and muons) switched on
1240         AliceInp << setw(10) << 1.0;  // restricted energy loss fluctuations (for e+ and e-) switched on
1241         AliceInp << setw(10) << 1.0;  // minimal accuracy
1242         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
1243         AliceInp << setprecision(2);
1244         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1245         AliceInp << endl;
1246       }
1247       else if (iProcessValue[i] == 4) { // no energy loss fluctuations
1248         AliceInp << "*"; 
1249         AliceInp << endl;
1250         AliceInp << "*No energy loss fluctuations"; 
1251         AliceInp << endl;
1252         AliceInp << "*Generated from call: SetProcess('LOSS',4)"; 
1253         AliceInp << endl;
1254         AliceInp << setw(10) << -1.0;  // restricted energy loss fluctuations (for hadrons and muons) switched off
1255         AliceInp << setw(10) << -1.0;  // restricted energy loss fluctuations (for e+ and e-) switched off
1256         AliceInp << setw(10) << 1.0;  // minimal accuracy
1257         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
1258         AliceInp << setprecision(2);
1259         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1260         AliceInp << endl;
1261       }
1262       else  {
1263         AliceInp << "*"; 
1264         AliceInp << endl;
1265         AliceInp << "*Illegal flag value in SetProcess('LOSS',?) call."; 
1266         AliceInp << endl;
1267         AliceInp << "*No FLUKA card generated"; 
1268         AliceInp << endl;
1269       }
1270     } // end of else if (strncmp(&sProcessFlag[i][0],"LOSS",4) == 0)
1271         
1272  
1273     // multiple scattering
1274     // G3 default value: 1
1275     // G4 process: G4MultipleScattering/G4IMultipleScattering
1276     // 
1277     // Particles: charged
1278     // Physics:   EM
1279     // flag = 0 no multiple scattering
1280     // flag = 1 Moliere or Coulomb scattering
1281     // flag = 2 Moliere or Coulomb scattering
1282     // flag = 3 Gaussian scattering
1283     // gMC ->SetProcess("MULS",1); // MULSOPT multiple scattering
1284     else if (strncmp(&sProcessFlag[i][0],"MULS",4) == 0) {
1285       if (iProcessValue[i] == 1 || iProcessValue[i] == 2 || iProcessValue[i] == 3) {
1286         AliceInp << "*"; 
1287         AliceInp << endl;
1288         AliceInp << "*Multiple scattering is ON by default for e+e- and for hadrons/muons"; 
1289         AliceInp << endl;
1290         AliceInp << "*No FLUKA card generated"; 
1291         AliceInp << endl;
1292       }
1293       else if (iProcessValue[i] == 0) {
1294         AliceInp << "*"; 
1295         AliceInp << endl;
1296         AliceInp << "*Multiple scattering is set OFF"; 
1297         AliceInp << endl;
1298         AliceInp << "*Generated from call: SetProcess('MULS',0);"; 
1299         AliceInp << endl;
1300         AliceInp << setw(10) << "MULSOPT  "; 
1301         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1302         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1303         AliceInp << setw(10) << 0.0;  // ignored
1304         AliceInp << setw(10) << 3.0;  // multiple scattering for hadrons and muons is completely suppressed
1305         AliceInp << setw(10) << 3.0;  // multiple scattering for e+ and e- is completely suppressed
1306         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
1307         AliceInp << setprecision(2);
1308         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1309         AliceInp << endl;
1310       }
1311       else  {
1312         AliceInp << "*"; 
1313         AliceInp << endl;
1314         AliceInp << "*Illegal flag value in SetProcess('MULS',?) call."; 
1315         AliceInp << endl;
1316         AliceInp << "*No FLUKA card generated"; 
1317         AliceInp << endl;
1318       }
1319     } // end of else if (strncmp(&sProcessFlag[i][0],"MULS",4) == 0)
1320
1321
1322     // muon nuclear interaction
1323     // G3 default value: 0
1324     // G4 processes: G4MuNuclearInteraction,
1325     // G4MuonMinusCaptureAtRest
1326     // 
1327     // Particles: mu
1328     // Physics:   Not set
1329     // flag = 0 no muon-nuclear interaction
1330     // flag = 1 nuclear interaction, secondaries processed
1331     // flag = 2 nuclear interaction, secondaries not processed
1332     // gMC ->SetProcess("MUNU",1); // MUPHOTON  1.   0.  0. 3. lastmat
1333     else if (strncmp(&sProcessFlag[i][0],"MUNU",4) == 0) {
1334       if (iProcessValue[i] == 1) {
1335         AliceInp << "*"; 
1336         AliceInp << endl;
1337         AliceInp << "*Muon nuclear interactions with production of secondary hadrons"; 
1338         AliceInp << endl;
1339         AliceInp << "*Generated from call: SetProcess('MUNU',1);"; 
1340         AliceInp << endl;
1341         AliceInp << setw(10) << "MUPHOTON  "; 
1342         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1343         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1344         AliceInp << setw(10) << 1.0;  // full simulation of muon nuclear interactions and production of secondary hadrons
1345         AliceInp << setw(10) << 0.0; // ratio of longitudinal to transverse virtual photon cross-section - Default = 0.25.
1346         AliceInp << setw(10) << 0.0; // fraction of rho-like interactions ( must be < 1) - Default = 0.75.
1347         AliceInp << setprecision(1);
1348         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
1349         AliceInp << setprecision(2);
1350         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1351         AliceInp << endl;
1352       }
1353       else if (iProcessValue[i] == 2) {
1354         AliceInp << "*"; 
1355         AliceInp << endl;
1356         AliceInp << "*Muon nuclear interactions without production of secondary hadrons"; 
1357         AliceInp << endl;
1358         AliceInp << "*Generated from call: SetProcess('MUNU',2);"; 
1359         AliceInp << endl;
1360         AliceInp << setw(10) << "MUPHOTON  "; 
1361         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1362         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1363         AliceInp << setw(10) << 2.0; // full simulation of muon nuclear interactions and production of secondary hadrons
1364         AliceInp << setw(10) << 0.0; // ratio of longitudinal to transverse virtual photon cross-section - Default = 0.25.
1365         AliceInp << setw(10) << 0.0; // fraction of rho-like interactions ( must be < 1) - Default = 0.75.
1366         AliceInp << setprecision(1);
1367         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
1368         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1369         AliceInp << endl;
1370       }
1371       else if (iProcessValue[i] == 0) {
1372         AliceInp << "*"; 
1373         AliceInp << endl;
1374         AliceInp << "*No muon nuclear interaction - no FLUKA card generated"; 
1375         AliceInp << endl;
1376         AliceInp << "*Generated from call: SetProcess('MUNU',0)"; 
1377         AliceInp << endl;
1378       }
1379       else  {
1380         AliceInp << "*"; 
1381         AliceInp << endl;
1382         AliceInp << "*Illegal flag value in SetProcess('MUNU',?) call."; 
1383         AliceInp << endl;
1384         AliceInp << "*No FLUKA card generated"; 
1385         AliceInp << endl;
1386       }
1387     } // end of else if (strncmp(&sProcessFlag[i][0],"MUNU",4) == 0)
1388
1389
1390     // photofission
1391     // G3 default value: 0
1392     // G4 process: ??
1393     //
1394     // Particles: gamma
1395     // Physics:   ??
1396     // gMC ->SetProcess("PFIS",0); // PHOTONUC -1.   0.  0. 3. lastmat 0.
1397     // flag = 0 no photon fission
1398     // flag = 1 photon fission, secondaries processed
1399     // flag = 2 photon fission, no secondaries stored
1400     else if (strncmp(&sProcessFlag[i][0],"PFIS",4) == 0) {
1401       if (iProcessValue[i] == 0) {
1402         AliceInp << "*"; 
1403         AliceInp << endl;
1404         AliceInp << "*No photonuclear interactions";
1405         AliceInp << endl;
1406         AliceInp << "*Generated from call: SetProcess('PFIS',0);"; 
1407         AliceInp << endl;
1408         AliceInp << setw(10) << "PHOTONUC  "; 
1409         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1410         AliceInp << setw(10) << -1.0; // no photonuclear interactions
1411         AliceInp << setw(10) << 0.0;  // not used
1412         AliceInp << setw(10) << 0.0;  // not used
1413         AliceInp << setw(10) << 3.0;  // upper bound of the material indices in which the respective thresholds apply
1414         AliceInp << setprecision(2); 
1415         AliceInp << setw(10) << fLastMaterial;
1416         AliceInp << setprecision(1);  // upper bound of the material indices in which the respective thresholds apply
1417         AliceInp << setprecision(1);
1418         AliceInp << setw(10) << 1.0;  // step length in assigning indices
1419         AliceInp << endl;
1420       }
1421       else if (iProcessValue[i] == 1) {
1422         AliceInp << "*"; 
1423         AliceInp << endl;
1424         AliceInp << "*Photon nuclear interactions are activated at all energies";
1425         AliceInp << endl;
1426         AliceInp << "*Generated from call: SetProcess('PFIS',1);"; 
1427         AliceInp << endl;
1428         AliceInp << setw(10) << "PHOTONUC  "; 
1429         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1430         AliceInp << setw(10) << 1.0; // photonuclear interactions are activated at all energies
1431         AliceInp << setw(10) << 0.0; // not used
1432         AliceInp << setw(10) << 0.0; // not used
1433         AliceInp << setprecision(2); 
1434         AliceInp << setw(10) << 3.0; // upper bound of the material indices in which the respective thresholds apply
1435         AliceInp << setw(10) << fLastMaterial;
1436         AliceInp << setprecision(1); // upper bound of the material indices in which the respective thresholds apply
1437         AliceInp << setprecision(1);
1438         AliceInp << setw(10) << 1.0; // step length in assigning indices
1439         AliceInp << endl;
1440       }
1441       else if (iProcessValue[i] == 0) {
1442         AliceInp << "*"; 
1443         AliceInp << endl;
1444         AliceInp << "*No photofission - no FLUKA card generated"; 
1445         AliceInp << endl;
1446         AliceInp << "*Generated from call: SetProcess('PFIS',0)"; 
1447         AliceInp << endl;
1448       }
1449       else {
1450         AliceInp << "*"; 
1451         AliceInp << endl;
1452         AliceInp << "*Illegal flag value in SetProcess('PFIS',?) call."; 
1453         AliceInp << endl;
1454         AliceInp << "*No FLUKA card generated"; 
1455         AliceInp << endl;
1456       }
1457     }
1458
1459  
1460     // photo electric effect
1461     // G3 default value: 1
1462     // G4 processes: G4PhotoElectricEffect
1463     //               G4LowEnergyPhotoElectric
1464     // Particles: gamma
1465     // Physics:   EM
1466     // flag = 0 no photo electric effect
1467     // flag = 1 photo electric effect, electron processed
1468     // flag = 2 photo electric effect, no electron stored
1469     // gMC ->SetProcess("PHOT",1); // EMFCUT    0.  -1.  0. 3. lastmat 0. PHOT-THR
1470     else if (strncmp(&sProcessFlag[i][0],"PHOT",4) == 0) {
1471       if (iProcessValue[i] == 1 || iProcessValue[i] == 2) {
1472         AliceInp << "*"; 
1473         AliceInp << endl;
1474         AliceInp << "*Photo electric effect is activated"; 
1475         AliceInp << endl;
1476         AliceInp << "*Generated from call: SetProcess('PHOT',1);"; 
1477         AliceInp << endl;
1478         AliceInp << setw(10) << "EMFCUT    "; 
1479         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1480         AliceInp << setw(10) << 0.0;  // ignored
1481         AliceInp << setw(10) << -1.0; // resets to default=0.
1482         AliceInp << setw(10) << 0.0;  // ignored
1483         AliceInp << setw(10) << 3.0;  // upper bound of the material indices in which the respective thresholds apply
1484         AliceInp << setprecision(2);
1485         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1486         AliceInp << setprecision(1);
1487         AliceInp << setw(10) << 1.0;  // step length in assigning indices
1488         AliceInp << setw(8) << "PHOT-THR"; 
1489         AliceInp << endl;
1490       }
1491       else if (iProcessValue[i] == 0) {
1492         AliceInp << "*"; 
1493         AliceInp << endl;
1494         AliceInp << "*No photo electric effect - no FLUKA card generated"; 
1495         AliceInp << endl;
1496         AliceInp << "*Generated from call: SetProcess('PHOT',0)"; 
1497         AliceInp << endl;
1498       }
1499       else {
1500         AliceInp << "*"; 
1501         AliceInp << endl;
1502         AliceInp << "*Illegal flag value in SetProcess('PHOT',?) call."; 
1503         AliceInp << endl;
1504         AliceInp << "*No FLUKA card generated"; 
1505         AliceInp << endl;
1506       }
1507     } // else if (strncmp(&sProcessFlag[i][0],"PHOT",4) == 0)
1508
1509
1510     // Rayleigh scattering
1511     // G3 default value: 0
1512     // G4 process: G4OpRayleigh
1513     // 
1514     // Particles: optical photon
1515     // Physics:   Optical
1516     // flag = 0 Rayleigh scattering off
1517     // flag = 1 Rayleigh scattering on
1518     //xx gMC ->SetProcess("RAYL",1);
1519     else if (strncmp(&sProcessFlag[i][0],"RAYL",4) == 0) {
1520       if (iProcessValue[i] == 1) {
1521         AliceInp << "*"; 
1522         AliceInp << endl;
1523         AliceInp << "*Rayleigh scattering is ON by default in FLUKA"; 
1524         AliceInp << endl;
1525         AliceInp << "*No FLUKA card generated"; 
1526         AliceInp << endl;
1527       }
1528       else if (iProcessValue[i] == 0) {
1529         AliceInp << "*"; 
1530         AliceInp << endl;
1531         AliceInp << "*Rayleigh scattering is set OFF"; 
1532         AliceInp << endl;
1533         AliceInp << "*Generated from call: SetProcess('RAYL',0);"; 
1534         AliceInp << endl;
1535         AliceInp << setw(10) << "EMFRAY    "; 
1536         AliceInp << setiosflags(ios::scientific) << setprecision(5);
1537         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1538         AliceInp << setw(10) << -1.0;  // no Rayleigh scattering and no binding corrections for Compton
1539         AliceInp << setw(10) << 3.0;  // lower bound of the material indices in which the respective thresholds apply
1540         AliceInp << setprecision(2);
1541         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1542         AliceInp << endl;
1543       }
1544       else  {
1545         AliceInp << "*"; 
1546         AliceInp << endl;
1547         AliceInp << "*Illegal flag value in SetProcess('RAYL',?) call."; 
1548         AliceInp << endl;
1549         AliceInp << "*No FLUKA card generated"; 
1550         AliceInp << endl;
1551       }
1552     } // end of else if (strncmp(&sProcessFlag[i][0],"RAYL",4) == 0)
1553
1554
1555     // synchrotron radiation in magnetic field
1556     // G3 default value: 0
1557     // G4 process: G4SynchrotronRadiation
1558     // 
1559     // Particles: ??
1560     // Physics:   Not set
1561     // flag = 0 no synchrotron radiation
1562     // flag = 1 synchrotron radiation
1563     //xx gMC ->SetProcess("SYNC",1); // synchrotron radiation generation
1564     else if (strncmp(&sProcessFlag[i][0],"SYNC",4) == 0) {
1565         AliceInp << "*"; 
1566         AliceInp << endl;
1567         AliceInp << "*Synchrotron radiation generation is NOT implemented in FLUKA"; 
1568         AliceInp << endl;
1569         AliceInp << "*No FLUKA card generated"; 
1570         AliceInp << endl;
1571     }
1572         
1573
1574     // Automatic calculation of tracking medium parameters
1575     // flag = 0 no automatic calculation
1576     // flag = 1 automatic calculation
1577     //xx gMC ->SetProcess("AUTO",1); // ??? automatic computation of the tracking medium parameters
1578     else if (strncmp(&sProcessFlag[i][0],"AUTO",4) == 0) {
1579         AliceInp << "*"; 
1580         AliceInp << endl;
1581         AliceInp << "*Automatic calculation of tracking medium parameters is always ON in FLUKA"; 
1582         AliceInp << endl;
1583         AliceInp << "*No FLUKA card generated"; 
1584         AliceInp << endl;
1585     }
1586         
1587
1588     // To control energy loss fluctuation model
1589     // flag = 0 Urban model
1590     // flag = 1 PAI model
1591     // flag = 2 PAI+ASHO model (not active at the moment)
1592     //xx gMC ->SetProcess("STRA",1); // ??? energy fluctuation model
1593     else if (strncmp(&sProcessFlag[i][0],"STRA",4) == 0) {
1594       if (iProcessValue[i] == 0 || iProcessValue[i] == 2 || iProcessValue[i] == 3) {
1595         AliceInp << "*";
1596         AliceInp << endl;
1597         AliceInp << "*Ionization energy losses calculation is activated";
1598         AliceInp << endl;
1599         AliceInp << "*Generated from call: SetProcess('STRA',n);, n=0,1,2";
1600         AliceInp << endl;
1601         AliceInp << setw(10) << "IONFLUCT  ";
1602         AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1603         AliceInp << setw(10) << 1.0;  // restricted energy loss fluctuations
1604                                       // (for hadrons and muons) switched on
1605         AliceInp << setw(10) << 1.0;  // restricted energy loss fluctuations
1606                                       // (for e+ and e-) switched on
1607         AliceInp << setw(10) << 1.0;  // minimal accuracy
1608         AliceInp << setw(10) << 3.0;  // upper bound of the material indices in
1609                                       // which the respective thresholds apply
1610         AliceInp << setprecision(2);
1611         AliceInp << setw(10) << fLastMaterial; // upper bound of the material indices in which the respective thresholds apply
1612         AliceInp << setprecision(1);
1613         AliceInp << setw(10) << 1.0;  // step length in assigning indices
1614         AliceInp << endl;
1615       }
1616       else {
1617         AliceInp << "*";
1618         AliceInp << endl;
1619         AliceInp << "*Illegal flag value in SetProcess('STRA',?) call.";
1620         AliceInp << endl;
1621         AliceInp << "*No FLUKA card generated";
1622         AliceInp << endl;
1623       }
1624     } // else if (strncmp(&sProcessFlag[i][0],"STRA",4) == 0)
1625
1626
1627
1628
1629     else { // processes not yet treated
1630
1631     // light photon absorption (Cerenkov photons)
1632     // it is turned on when Cerenkov process is turned on
1633     // G3 default value: 0
1634     // G4 process: G4OpAbsorption, G4OpBoundaryProcess
1635     // 
1636     // Particles: optical photon
1637     // Physics:   Optical
1638     // flag = 0 no absorption of Cerenkov photons
1639     // flag = 1 absorption of Cerenkov photons
1640     // gMC ->SetProcess("LABS",2); // ??? Cerenkov light absorption
1641
1642
1643
1644       cout << "SetProcess for flag=" << &sProcessFlag[i][0] << " value=" << iProcessValue[i] << " not yet implemented!" << endl;
1645     }
1646   } //end of loop number of SetProcess calls
1647
1648  
1649 // Loop over number of SetCut calls  
1650   for (Int_t i=0; i<iNbOfCut; i++) {
1651
1652     // cuts used in SetProcess calls
1653     if (strncmp(&sCutFlag[i][0],"BCUTM",5) == 0) continue;
1654     else if (strncmp(&sCutFlag[i][0],"BCUTE",5) == 0) continue;
1655     else if (strncmp(&sCutFlag[i][0],"DCUTM",5) == 0) continue;
1656     else if (strncmp(&sCutFlag[i][0],"PPCUTM",6) == 0) continue;
1657
1658     // gammas
1659     // G4 particles: "gamma"
1660     // G3 default value: 0.001 GeV
1661     //gMC ->SetCut("CUTGAM",cut); // cut for gammas
1662     else if (strncmp(&sCutFlag[i][0],"CUTGAM",6) == 0) {
1663       AliceInp << "*"; 
1664       AliceInp << endl;
1665       AliceInp << "*Cut for gamma"; 
1666       AliceInp << endl;
1667       AliceInp << "*Generated from call: SetCut('CUTGAM',cut);"; 
1668       AliceInp << endl;
1669       AliceInp << setw(10) << "PART-THR  "; 
1670       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1671       AliceInp << setw(10) << -fCutValue[i];
1672       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1673       AliceInp << setw(10) << 7.0;
1674       AliceInp << endl;
1675     }
1676
1677     // electrons
1678     // G4 particles: "e-"
1679     // ?? positrons
1680     // G3 default value: 0.001 GeV
1681     //gMC ->SetCut("CUTELE",cut); // cut for e+,e-
1682     else if (strncmp(&sCutFlag[i][0],"CUTELE",6) == 0) {
1683       AliceInp << "*"; 
1684       AliceInp << endl;
1685       AliceInp << "*Cut for electrons"; 
1686       AliceInp << endl;
1687       AliceInp << "*Generated from call: SetCut('CUTELE',cut);"; 
1688       AliceInp << endl;
1689       AliceInp << setw(10) << "PART-THR  "; 
1690       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1691       AliceInp << setw(10) << -fCutValue[i];
1692       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1693       AliceInp << setw(10) << 3.0;
1694       AliceInp << setw(10) << 4.0;
1695       AliceInp << setw(10) << 1.0;
1696       AliceInp << endl;
1697     }
1698
1699     // neutral hadrons
1700     // G4 particles: of type "baryon", "meson", "nucleus" with zero charge
1701     // G3 default value: 0.01 GeV
1702     //gMC ->SetCut("CUTNEU",cut); // cut for neutral hadrons
1703     else if (strncmp(&sCutFlag[i][0],"CUTNEU",6) == 0) {
1704       AliceInp << "*"; 
1705       AliceInp << endl;
1706       AliceInp << "*Cut for neutral hadrons"; 
1707       AliceInp << endl;
1708       AliceInp << "*Generated from call: SetCut('CUTNEU',cut);"; 
1709       AliceInp << endl;
1710       AliceInp << setw(10) << "PART-THR  "; 
1711       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1712       AliceInp << setw(10) << -fCutValue[i];
1713       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1714       AliceInp << setw(10) << 8.0; // Neutron
1715       AliceInp << setw(10) << 9.0; // Antineutron
1716       AliceInp << endl;
1717       AliceInp << setw(10) << "PART-THR  "; 
1718       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1719       AliceInp << setw(10) << -fCutValue[i];
1720       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1721       AliceInp << setw(10) << 12.0; // Kaon zero long
1722       AliceInp << setw(10) << 12.0; // Kaon zero long
1723       AliceInp << endl;
1724       AliceInp << setw(10) << "PART-THR  "; 
1725       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1726       AliceInp << setw(10) << -fCutValue[i];
1727       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1728       AliceInp << setw(10) << 17.0; // Lambda, 18=Antilambda
1729       AliceInp << setw(10) << 19.0; // Kaon zero short
1730       AliceInp << endl;
1731       AliceInp << setw(10) << "PART-THR  "; 
1732       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1733       AliceInp << setw(10) << -fCutValue[i];
1734       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1735       AliceInp << setw(10) << 22.0; // Sigma zero, Pion zero, Kaon zero
1736       AliceInp << setw(10) << 25.0; // Antikaon zero
1737       AliceInp << endl;
1738       AliceInp << setw(10) << "PART-THR  "; 
1739       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1740       AliceInp << setw(10) << -fCutValue[i];
1741       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1742       AliceInp << setw(10) << 32.0; // Antisigma zero
1743       AliceInp << setw(10) << 32.0; // Antisigma zero
1744       AliceInp << endl;
1745       AliceInp << setw(10) << "PART-THR  "; 
1746       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1747       AliceInp << setw(10) << -fCutValue[i];
1748       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1749       AliceInp << setw(10) << 34.0; // Xi zero
1750       AliceInp << setw(10) << 35.0; // AntiXi zero
1751       AliceInp << endl;
1752       AliceInp << setw(10) << "PART-THR  "; 
1753       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1754       AliceInp << setw(10) << -fCutValue[i];
1755       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1756       AliceInp << setw(10) << 47.0; // D zero
1757       AliceInp << setw(10) << 48.0; // AntiD zero
1758       AliceInp << endl;
1759       AliceInp << setw(10) << "PART-THR  "; 
1760       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1761       AliceInp << setw(10) << -fCutValue[i];
1762       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1763       AliceInp << setw(10) << 53.0; // Xi_c zero
1764       AliceInp << setw(10) << 53.0; // Xi_c zero
1765       AliceInp << endl;
1766       AliceInp << setw(10) << "PART-THR  "; 
1767       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1768       AliceInp << setw(10) << -fCutValue[i];
1769       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1770       AliceInp << setw(10) << 55.0; // Xi'_c zero
1771       AliceInp << setw(10) << 56.0; // Omega_c zero
1772       AliceInp << endl;
1773       AliceInp << setw(10) << "PART-THR  "; 
1774       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1775       AliceInp << setw(10) << -fCutValue[i];
1776       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1777       AliceInp << setw(10) << 59.0; // AntiXi_c zero
1778       AliceInp << setw(10) << 59.0; // AntiXi_c zero
1779       AliceInp << endl;
1780       AliceInp << setw(10) << "PART-THR  "; 
1781       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1782       AliceInp << setw(10) << -fCutValue[i];
1783       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1784       AliceInp << setw(10) << 61.0; // AntiXi'_c zero
1785       AliceInp << setw(10) << 62.0; // AntiOmega_c zero
1786       AliceInp << endl;
1787     }
1788
1789     // charged hadrons
1790     // G4 particles: of type "baryon", "meson", "nucleus" with non-zero charge
1791     // G3 default value: 0.01 GeV
1792     //gMC ->SetCut("CUTHAD",cut); // cut for charged hadrons
1793     else if (strncmp(&sCutFlag[i][0],"CUTHAD",6) == 0) {
1794       AliceInp << "*"; 
1795       AliceInp << endl;
1796       AliceInp << "*Cut for charged hadrons"; 
1797       AliceInp << endl;
1798       AliceInp << "*Generated from call: SetCut('CUTHAD',cut);"; 
1799       AliceInp << endl;
1800       AliceInp << setw(10) << "PART-THR  "; 
1801       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1802       AliceInp << setw(10) << -fCutValue[i];
1803       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1804       AliceInp << setw(10) << 1.0; // Proton
1805       AliceInp << setw(10) << 2.0; // Antiproton
1806       AliceInp << endl;
1807       AliceInp << setw(10) << "PART-THR  "; 
1808       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1809       AliceInp << setw(10) << -fCutValue[i];
1810       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1811       AliceInp << setw(10) << 13.0; // Positive Pion, Negative Pion, Positive Kaon
1812       AliceInp << setw(10) << 16.0; // Negative Kaon
1813       AliceInp << endl;
1814       AliceInp << setw(10) << "PART-THR  "; 
1815       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1816       AliceInp << setw(10) << -fCutValue[i];
1817       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1818       AliceInp << setw(10) << 20.0; // Negative Sigma
1819       AliceInp << setw(10) << 16.0; // Positive Sigma
1820       AliceInp << endl;
1821       AliceInp << setw(10) << "PART-THR  "; 
1822       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1823       AliceInp << setw(10) << -fCutValue[i];
1824       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1825       AliceInp << setw(10) << 31.0; // Antisigma minus
1826       AliceInp << setw(10) << 33.0; // Antisigma plus
1827       AliceInp << setprecision(1);
1828       AliceInp << setw(10) << 2.0;  // step length
1829       AliceInp << endl;
1830       AliceInp << setw(10) << "PART-THR  "; 
1831       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1832       AliceInp << setw(10) << -fCutValue[i];
1833       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1834       AliceInp << setw(10) << 36.0; // Negative Xi, Positive Xi, Omega minus
1835       AliceInp << setw(10) << 39.0; // Antiomega
1836       AliceInp << endl;
1837       AliceInp << setw(10) << "PART-THR  "; 
1838       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1839       AliceInp << setw(10) << -fCutValue[i];
1840       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1841       AliceInp << setw(10) << 45.0; // D plus
1842       AliceInp << setw(10) << 46.0; // D minus
1843       AliceInp << endl;
1844       AliceInp << setw(10) << "PART-THR  "; 
1845       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1846       AliceInp << setw(10) << -fCutValue[i];
1847       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1848       AliceInp << setw(10) << 49.0; // D_s plus, D_s minus, Lambda_c plus
1849       AliceInp << setw(10) << 52.0; // Xi_c plus
1850       AliceInp << endl;
1851       AliceInp << setw(10) << "PART-THR  "; 
1852       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1853       AliceInp << setw(10) << -fCutValue[i];
1854       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1855       AliceInp << setw(10) << 54.0; // Xi'_c plus
1856       AliceInp << setw(10) << 60.0; // AntiXi'_c minus
1857       AliceInp << setprecision(1);
1858       AliceInp << setw(10) << 6.0;  // step length
1859       AliceInp << endl;
1860       AliceInp << setw(10) << "PART-THR  "; 
1861       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1862       AliceInp << setw(10) << -fCutValue[i];
1863       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(2);
1864       AliceInp << setw(10) << 57.0; // Antilambda_c minus
1865       AliceInp << setw(10) << 58.0; // AntiXi_c minus
1866       AliceInp << endl;
1867     }
1868
1869     // muons
1870     // G4 particles: "mu+", "mu-"
1871     // G3 default value: 0.01 GeV
1872     //gMC ->SetCut("CUTMUO",cut); // cut for mu+, mu-
1873     else if (strncmp(&sCutFlag[i][0],"CUTMUO",6) == 0) {
1874       AliceInp << "*"; 
1875       AliceInp << endl;
1876       AliceInp << "*Cut for muons"; 
1877       AliceInp << endl;
1878       AliceInp << "*Generated from call: SetCut('CUTMUO',cut);"; 
1879       AliceInp << endl;
1880       AliceInp << setw(10) << "PART-THR  "; 
1881       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1882       AliceInp << setw(10) << -fCutValue[i];
1883       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1884       AliceInp << setprecision(2);
1885       AliceInp << setw(10) << 10.0;
1886       AliceInp << setw(10) << 11.0;
1887       AliceInp << endl;
1888     }
1889     // delta-rays by electrons
1890     // G4 particles: "e-"
1891     // G3 default value: 10**4 GeV
1892     // gMC ->SetCut("DCUTE",cut);  // cut for deltarays by electrons ???????????????
1893     else if (strncmp(&sCutFlag[i][0],"DCUTE",5) == 0) {
1894       AliceInp << "*"; 
1895       AliceInp << endl;
1896       AliceInp << "*Cut for delta rays by electrons ????????????";
1897       AliceInp << endl;
1898       AliceInp << "*Generated from call: SetCut('DCUTE',cut);";
1899       AliceInp << endl;
1900       AliceInp << setw(10) << "EMFCUT    ";
1901       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1902       AliceInp << setw(10) << -fCutValue[i];
1903       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1904       AliceInp << setw(10) << 0.0;
1905       AliceInp << setw(10) << 0.0;
1906       AliceInp << setw(10) << 3.0;
1907       AliceInp << setprecision(2);
1908       AliceInp << setw(10) << fLastMaterial;
1909       AliceInp << setprecision(1);
1910       AliceInp << setw(10) << 1.0;
1911       AliceInp << endl;
1912     }
1913     
1914     //
1915     // time of flight cut in seconds
1916     // G4 particles: all
1917     // G3 default value: 0.01 GeV
1918     //gMC ->SetCut("TOFMAX",tofmax); // time of flight cuts in seconds
1919     else if (strncmp(&sCutFlag[i][0],"TOFMAX",6) == 0) {
1920       AliceInp << "*"; 
1921       AliceInp << endl;
1922       AliceInp << "*Time of flight cuts in seconds"; 
1923       AliceInp << endl;
1924       AliceInp << "*Generated from call: SetCut('TOFMAX',tofmax);"; 
1925       AliceInp << endl;
1926       AliceInp << setw(10) << "TIME-CUT  "; 
1927       AliceInp << setiosflags(ios::scientific) << setprecision(5);
1928       AliceInp << setw(10) << fCutValue[i]*1.e9;
1929       AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint) << setprecision(1);
1930       AliceInp << setw(10) << 0.0;
1931       AliceInp << setw(10) << 0.0;
1932       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
1933       AliceInp << setprecision(2);
1934       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
1935       AliceInp << setprecision(1);
1936       AliceInp << setw(10) << 1.0; // step length in assigning numbers
1937       AliceInp << endl;
1938     }
1939
1940     else {
1941       cout << "SetCut for flag=" << &sCutFlag[i][0] << " value=" << fCutValue[i] << " not yet implemented!" << endl;
1942     }
1943   } //end of loop over SeCut calls
1944     
1945 // Add START and STOP card
1946    AliceInp << setw(10) << "START     "; 
1947    AliceInp << setiosflags(ios::fixed) << setiosflags(ios::showpoint);
1948    AliceInp << setw(10) << fEventsPerRun;
1949    AliceInp << endl;
1950    AliceInp << setw(10) << "STOP      "; 
1951    AliceInp << endl;
1952
1953 } // end of InitPhysics
1954
1955
1956 void TFluka::SetMaxStep(Double_t)
1957 {
1958 // SetMaxStep is dummy procedure in TFluka !
1959   if (fVerbosityLevel >=3)
1960   cout << "SetMaxStep is dummy procedure in TFluka !" << endl;
1961 }
1962
1963 void TFluka::SetMaxNStep(Int_t)
1964 {
1965 // SetMaxNStep is dummy procedure in TFluka !
1966   if (fVerbosityLevel >=3)
1967   cout << "SetMaxNStep is dummy procedure in TFluka !" << endl;
1968 }
1969
1970 void TFluka::SetUserDecay(Int_t)
1971 {
1972 // SetUserDecay is dummy procedure in TFluka !
1973   if (fVerbosityLevel >=3)
1974   cout << "SetUserDecay is dummy procedure in TFluka !" << endl;
1975 }
1976
1977 //
1978 // dynamic properties
1979 //
1980 void TFluka::TrackPosition(TLorentzVector& position) const
1981 {
1982 // Return the current position in the master reference frame of the
1983 // track being transported
1984 // TRACKR.atrack = age of the particle
1985 // TRACKR.xtrack = x-position of the last point
1986 // TRACKR.ytrack = y-position of the last point
1987 // TRACKR.ztrack = z-position of the last point
1988   Int_t caller = GetCaller();
1989   if (caller == 3 || caller == 6 || caller == 11 || caller == 12) { //bxdraw,endraw,usdraw
1990     position.SetX(GetXsco());
1991     position.SetY(GetYsco());
1992     position.SetZ(GetZsco());
1993     position.SetT(TRACKR.atrack);
1994   }
1995   else if (caller == 4) { // mgdraw
1996     position.SetX(TRACKR.xtrack[TRACKR.ntrack]);
1997     position.SetY(TRACKR.ytrack[TRACKR.ntrack]);
1998     position.SetZ(TRACKR.ztrack[TRACKR.ntrack]);
1999     position.SetT(TRACKR.atrack);
2000   }
2001   else if (caller == 5) { // sodraw
2002     position.SetX(TRACKR.xtrack[TRACKR.ntrack]);
2003     position.SetY(TRACKR.ytrack[TRACKR.ntrack]);
2004     position.SetZ(TRACKR.ztrack[TRACKR.ntrack]);
2005     position.SetT(0);
2006   }
2007   else
2008     Warning("TrackPosition","position not available");
2009 }
2010
2011 //
2012 void TFluka::TrackPosition(Double_t& x, Double_t& y, Double_t& z) const
2013 {
2014 // Return the current position in the master reference frame of the
2015 // track being transported
2016 // TRACKR.atrack = age of the particle
2017 // TRACKR.xtrack = x-position of the last point
2018 // TRACKR.ytrack = y-position of the last point
2019 // TRACKR.ztrack = z-position of the last point
2020   Int_t caller = GetCaller();
2021   if (caller == 3 || caller == 6 || caller == 11 || caller == 12) { //bxdraw,endraw,usdraw
2022     x = GetXsco();
2023     y = GetYsco();
2024     z = GetZsco();
2025   }
2026   else if (caller == 4) { // mgdraw
2027     x = TRACKR.xtrack[TRACKR.ntrack];
2028     y = TRACKR.ytrack[TRACKR.ntrack];
2029     z = TRACKR.ztrack[TRACKR.ntrack];
2030   }
2031   else if (caller == 5) { // sodraw
2032     x = TRACKR.xtrack[TRACKR.ntrack];
2033     y = TRACKR.ytrack[TRACKR.ntrack];
2034     z = TRACKR.ztrack[TRACKR.ntrack];
2035   }
2036   else
2037     Warning("TrackPosition","position not available");
2038 }
2039
2040 void TFluka::TrackMomentum(TLorentzVector& momentum) const
2041 {
2042 // Return the direction and the momentum (GeV/c) of the track
2043 // currently being transported
2044 // TRACKR.ptrack = momentum of the particle (not always defined, if
2045 //               < 0 must be obtained from etrack) 
2046 // TRACKR.cx,y,ztrck = direction cosines of the current particle
2047 // TRACKR.etrack = total energy of the particle
2048 // TRACKR.jtrack = identity number of the particle
2049 // PAPROP.am[TRACKR.jtrack] = particle mass in gev
2050   Int_t caller = GetCaller();
2051   if (caller != 2) { // not eedraw 
2052     if (TRACKR.ptrack >= 0) {
2053       momentum.SetPx(TRACKR.ptrack*TRACKR.cxtrck);
2054       momentum.SetPy(TRACKR.ptrack*TRACKR.cytrck);
2055       momentum.SetPz(TRACKR.ptrack*TRACKR.cztrck);
2056       momentum.SetE(TRACKR.etrack);
2057       return;
2058     }
2059     else {
2060       Double_t p = sqrt(TRACKR.etrack*TRACKR.etrack - PAPROP.am[TRACKR.jtrack+6]*PAPROP.am[TRACKR.jtrack+6]);
2061       momentum.SetPx(p*TRACKR.cxtrck);
2062       momentum.SetPy(p*TRACKR.cytrck);
2063       momentum.SetPz(p*TRACKR.cztrck);
2064       momentum.SetE(TRACKR.etrack);
2065       return;
2066     }
2067   }
2068   else
2069     Warning("TrackMomentum","momentum not available");
2070 }
2071
2072 void TFluka::TrackMomentum(Double_t& px, Double_t& py, Double_t& pz, Double_t& e) const
2073 {
2074 // Return the direction and the momentum (GeV/c) of the track
2075 // currently being transported
2076 // TRACKR.ptrack = momentum of the particle (not always defined, if
2077 //               < 0 must be obtained from etrack) 
2078 // TRACKR.cx,y,ztrck = direction cosines of the current particle
2079 // TRACKR.etrack = total energy of the particle
2080 // TRACKR.jtrack = identity number of the particle
2081 // PAPROP.am[TRACKR.jtrack] = particle mass in gev
2082   Int_t caller = GetCaller();
2083   if (caller != 2) { // not eedraw 
2084     if (TRACKR.ptrack >= 0) {
2085       px = TRACKR.ptrack*TRACKR.cxtrck;
2086       py = TRACKR.ptrack*TRACKR.cytrck;
2087       pz = TRACKR.ptrack*TRACKR.cztrck;
2088       e = TRACKR.etrack;
2089       return;
2090     }
2091     else {
2092       Double_t p = sqrt(TRACKR.etrack*TRACKR.etrack - PAPROP.am[TRACKR.jtrack+6]*PAPROP.am[TRACKR.jtrack+6]);
2093       px = p*TRACKR.cxtrck;
2094       py = p*TRACKR.cytrck;
2095       pz = p*TRACKR.cztrck;
2096       e = TRACKR.etrack;
2097       return;
2098     }
2099   }
2100   else
2101     Warning("TrackMomentum","momentum not available");
2102 }
2103
2104 Double_t TFluka::TrackStep() const
2105 {
2106 // Return the length in centimeters of the current step
2107 // TRACKR.ctrack = total curved path
2108   Int_t caller = GetCaller();
2109   if (caller == 11 || caller==12 || caller == 3 || caller == 6) //bxdraw,endraw,usdraw
2110     return 0.0;
2111   else if (caller == 4) //mgdraw
2112     return TRACKR.ctrack;
2113   else
2114     return -1.0;
2115 }
2116
2117 Double_t TFluka::TrackLength() const
2118 {
2119 // TRACKR.cmtrck = cumulative curved path since particle birth
2120   Int_t caller = GetCaller();
2121   if (caller == 111 || caller==12 || caller == 3 || caller == 4 || caller == 6) //bxdraw,endraw,mgdraw,usdraw
2122     return TRACKR.cmtrck;
2123   else 
2124     return -1.0;
2125 }
2126
2127 Double_t TFluka::TrackTime() const
2128 {
2129 // Return the current time of flight of the track being transported
2130 // TRACKR.atrack = age of the particle
2131   Int_t caller = GetCaller();
2132   if (caller == 11 || caller==12 || caller == 3 || caller == 4 || caller == 6) //bxdraw,endraw,mgdraw,usdraw
2133     return TRACKR.atrack;
2134   else 
2135     return -1;
2136 }
2137
2138 Double_t TFluka::Edep() const
2139 {
2140 // Energy deposition
2141 // if TRACKR.ntrack = 0, TRACKR.mtrack = 0:
2142 // -->local energy deposition (the value and the point are not recorded in TRACKR)
2143 //    but in the variable "rull" of the procedure "endraw.cxx"
2144 // if TRACKR.ntrack > 0, TRACKR.mtrack = 0:
2145 // -->no energy loss along the track
2146 // if TRACKR.ntrack > 0, TRACKR.mtrack > 0:
2147 // -->energy loss distributed along the track
2148 // TRACKR.dtrack = energy deposition of the jth deposition even
2149
2150   // If coming from bxdraw we have 2 steps of 0 length and 0 edep
2151   Int_t caller = GetCaller();
2152   if (caller == 11 || caller==12) return 0.0;
2153   Double_t sum = 0;
2154   for ( Int_t j=0;j<TRACKR.mtrack;j++) {
2155     sum +=TRACKR.dtrack[j];  
2156   }
2157   if (TRACKR.ntrack == 0 && TRACKR.mtrack == 0)
2158     return fRull + sum;
2159   else {
2160     return sum;
2161   }
2162 }
2163
2164 Int_t TFluka::TrackPid() const
2165 {
2166 // Return the id of the particle transported
2167 // TRACKR.jtrack = identity number of the particle
2168   Int_t caller = GetCaller();
2169   if (caller != 2)  // not eedraw 
2170     return PDGFromId(TRACKR.jtrack);
2171   else
2172     return -1000;
2173 }
2174
2175 Double_t TFluka::TrackCharge() const
2176 {
2177 // Return charge of the track currently transported
2178 // PAPROP.ichrge = electric charge of the particle
2179 // TRACKR.jtrack = identity number of the particle
2180   Int_t caller = GetCaller();
2181   if (caller != 2)  // not eedraw 
2182     return PAPROP.ichrge[TRACKR.jtrack+6];
2183   else
2184     return -1000.0;
2185 }
2186
2187 Double_t TFluka::TrackMass() const
2188 {
2189 // PAPROP.am = particle mass in GeV
2190 // TRACKR.jtrack = identity number of the particle
2191   Int_t caller = GetCaller();
2192   if (caller != 2) {  // not eedraw 
2193 //    cout << "JTRACK=" << TRACKR.jtrack << " mass=" << PAPROP.am[TRACKR.jtrack+6] << endl;
2194     return PAPROP.am[TRACKR.jtrack+6];
2195   }
2196   else
2197     return -1000.0;
2198 }
2199
2200 Double_t TFluka::Etot() const
2201 {
2202 // TRACKR.etrack = total energy of the particle
2203   Int_t caller = GetCaller();
2204   if (caller != 2)  // not eedraw
2205     return TRACKR.etrack;
2206   else
2207     return -1000.0;
2208 }
2209
2210 //
2211 // track status
2212 //
2213 Bool_t   TFluka::IsNewTrack() const
2214 {
2215 // Return true for the first call of Stepping()
2216 /*
2217 // True if the track has positive cummulative length
2218   Int_t caller = GetCaller();
2219   if (caller != 2) { // not eedraw
2220     if (TRACKR.cmtrck > 0.0)
2221       return 1; 
2222     else
2223       return 0; 
2224   }
2225   else
2226     return 0;
2227 */    
2228    return fTrackIsNew;
2229 }
2230
2231 Bool_t   TFluka::IsTrackInside() const
2232 {
2233 // True if the track is not at the boundary of the current volume
2234 // In Fluka a step is always inside one kind of material
2235 // If the step would go behind the region of one material,
2236 // it will be shortened to reach only the boundary.
2237 // Therefore IsTrackInside() is always true.
2238   Int_t caller = GetCaller();
2239   if (caller == 11 || caller==12)  // bxdraw
2240     return 0;
2241   else
2242     return 1;
2243 }
2244
2245 Bool_t   TFluka::IsTrackEntering() const
2246 {
2247 // True if this is the first step of the track in the current volume
2248
2249   Int_t caller = GetCaller();
2250   if (caller == 11)  // bxdraw entering
2251     return 1;
2252   else return 0;
2253 }
2254
2255 Bool_t   TFluka::IsTrackExiting() const
2256 {
2257   Int_t caller = GetCaller();
2258   if (caller == 12)  // bxdraw exiting
2259     return 1;
2260   else return 0;
2261 }
2262
2263 Bool_t   TFluka::IsTrackOut() const
2264 {
2265 // True if the track is out of the setup
2266 // means escape
2267 // Icode = 14: escape - call from Kaskad
2268 // Icode = 23: escape - call from Emfsco
2269 // Icode = 32: escape - call from Kasneu
2270 // Icode = 40: escape - call from Kashea
2271 // Icode = 51: escape - call from Kasoph
2272   if (fIcode == 14 ||
2273       fIcode == 23 ||
2274       fIcode == 32 ||
2275       fIcode == 40 ||
2276       fIcode == 51) return 1;
2277   else return 0;
2278 }
2279
2280 Bool_t   TFluka::IsTrackDisappeared() const
2281 {
2282 // means all inelastic interactions and decays
2283 // fIcode from usdraw
2284   if (fIcode == 101 || // inelastic interaction
2285       fIcode == 102 || // particle decay
2286       fIcode == 214 || // in-flight annihilation
2287       fIcode == 215 || // annihilation at rest
2288       fIcode == 217 || // pair production
2289       fIcode == 221) return 1;
2290   else return 0;
2291 }
2292
2293 Bool_t   TFluka::IsTrackStop() const
2294 {
2295 // True if the track energy has fallen below the threshold
2296 // means stopped by signal or below energy threshold
2297 // Icode = 12: stopping particle       - call from Kaskad
2298 // Icode = 15: time kill               - call from Kaskad
2299 // Icode = 21: below threshold, iarg=1 - call from Emfsco
2300 // Icode = 22: below threshold, iarg=2 - call from Emfsco
2301 // Icode = 24: time kill               - call from Emfsco
2302 // Icode = 31: below threshold         - call from Kasneu
2303 // Icode = 33: time kill               - call from Kasneu
2304 // Icode = 41: time kill               - call from Kashea
2305 // Icode = 52: time kill               - call from Kasoph
2306   if (fIcode == 12 ||
2307       fIcode == 15 ||
2308       fIcode == 21 ||
2309       fIcode == 22 ||
2310       fIcode == 24 ||
2311       fIcode == 31 ||
2312       fIcode == 33 ||
2313       fIcode == 41 ||
2314       fIcode == 52) return 1;
2315   else return 0;
2316 }
2317
2318 Bool_t   TFluka::IsTrackAlive() const
2319 {
2320 // means not disappeared or not out
2321   if (IsTrackDisappeared() || IsTrackOut() ) return 0;
2322   else return 1;
2323 }
2324
2325 //
2326 // secondaries
2327 //
2328
2329 Int_t TFluka::NSecondaries() const
2330 // Number of secondary particles generated in the current step
2331 // FINUC.np = number of secondaries except light and heavy ions
2332 // FHEAVY.npheav = number of secondaries for light and heavy secondary ions
2333 {
2334   Int_t caller = GetCaller();
2335   if (caller == 6)  // valid only after usdraw
2336     return FINUC.np + FHEAVY.npheav;
2337   else
2338     return 0;
2339 } // end of NSecondaries
2340
2341 void     TFluka::GetSecondary(Int_t isec, Int_t& particleId,
2342                 TLorentzVector& position, TLorentzVector& momentum)
2343 {
2344   Int_t caller = GetCaller();
2345   if (caller == 6) {  // valid only after usdraw
2346     if (isec >= 0 && isec < FINUC.np) {
2347       particleId = PDGFromId(FINUC.kpart[isec]);
2348       position.SetX(fXsco);
2349       position.SetY(fYsco);
2350       position.SetZ(fZsco);
2351       position.SetT(TRACKR.atrack);
2352 //    position.SetT(TRACKR.atrack+FINUC.agesec[isec]); //not yet implem.
2353       momentum.SetPx(FINUC.plr[isec]*FINUC.cxr[isec]);
2354       momentum.SetPy(FINUC.plr[isec]*FINUC.cyr[isec]);
2355       momentum.SetPz(FINUC.plr[isec]*FINUC.czr[isec]);
2356       momentum.SetE(FINUC.tki[isec] + PAPROP.am[FINUC.kpart[isec]+6]);
2357     }
2358     else if (isec >= FINUC.np && isec < FINUC.np + FHEAVY.npheav) {
2359       Int_t jsec = isec - FINUC.np;
2360       particleId = FHEAVY.kheavy[jsec]; // this is Fluka id !!!
2361       position.SetX(fXsco);
2362       position.SetY(fYsco);
2363       position.SetZ(fZsco);
2364       position.SetT(TRACKR.atrack);
2365 //    position.SetT(TRACKR.atrack+FHEAVY.agheav[jsec]); //not yet implem.
2366       momentum.SetPx(FHEAVY.pheavy[jsec]*FHEAVY.cxheav[jsec]);
2367       momentum.SetPy(FHEAVY.pheavy[jsec]*FHEAVY.cyheav[jsec]);
2368       momentum.SetPz(FHEAVY.pheavy[jsec]*FHEAVY.czheav[jsec]);
2369       if (FHEAVY.tkheav[jsec] >= 3 && FHEAVY.tkheav[jsec] <= 6) 
2370         momentum.SetE(FHEAVY.tkheav[jsec] + PAPROP.am[jsec+6]);
2371       else if (FHEAVY.tkheav[jsec] > 6)
2372         momentum.SetE(FHEAVY.tkheav[jsec] + FHEAVY.amnhea[jsec]); // to be checked !!!
2373     }
2374     else
2375       Warning("GetSecondary","isec out of range");
2376   }
2377   else
2378     Warning("GetSecondary","no secondaries available");
2379 } // end of GetSecondary
2380
2381 TMCProcess TFluka::ProdProcess(Int_t) const
2382 // Name of the process that has produced the secondary particles
2383 // in the current step
2384 {
2385     const TMCProcess kIpNoProc = kPNoProcess;
2386     const TMCProcess kIpPDecay = kPDecay;
2387     const TMCProcess kIpPPair = kPPair;
2388 //  const TMCProcess kIpPPairFromPhoton = kPPairFromPhoton;
2389 //  const TMCProcess kIpPPairFromVirtualPhoton = kPPairFromVirtualPhoton;
2390     const TMCProcess kIpPCompton = kPCompton;
2391     const TMCProcess kIpPPhotoelectric = kPPhotoelectric;
2392     const TMCProcess kIpPBrem = kPBrem;
2393 //  const TMCProcess kIpPBremFromHeavy = kPBremFromHeavy;
2394 //  const TMCProcess kIpPBremFromElectronOrPositron = kPBremFromElectronOrPositron;
2395     const TMCProcess kIpPDeltaRay = kPDeltaRay;
2396 //  const TMCProcess kIpPMoller = kPMoller;
2397 //  const TMCProcess kIpPBhabha = kPBhabha;
2398     const TMCProcess kIpPAnnihilation = kPAnnihilation;
2399 //  const TMCProcess kIpPAnnihilInFlight = kPAnnihilInFlight;
2400 //  const TMCProcess kIpPAnnihilAtRest = kPAnnihilAtRest;
2401     const TMCProcess kIpPHadronic = kPHadronic;
2402     const TMCProcess kIpPMuonNuclear = kPMuonNuclear;
2403     const TMCProcess kIpPPhotoFission = kPPhotoFission;
2404     const TMCProcess kIpPRayleigh = kPRayleigh;
2405 //  const TMCProcess kIpPCerenkov = kPCerenkov;
2406 //  const TMCProcess kIpPSynchrotron = kPSynchrotron;
2407
2408     Int_t mugamma = TRACKR.jtrack == 7 || TRACKR.jtrack == 10 || TRACKR.jtrack == 11;
2409     if (fIcode == 102) return kIpPDecay;
2410     else if (fIcode == 104 || fIcode == 217) return kIpPPair;
2411 //  else if (fIcode == 104) return kIpPairFromPhoton;
2412 //  else if (fIcode == 217) return kIpPPairFromVirtualPhoton;
2413     else if (fIcode == 219) return kIpPCompton;
2414     else if (fIcode == 221) return kIpPPhotoelectric;
2415     else if (fIcode == 105 || fIcode == 208) return kIpPBrem;
2416 //  else if (fIcode == 105) return kIpPBremFromHeavy;
2417 //  else if (fIcode == 208) return kPBremFromElectronOrPositron;
2418     else if (fIcode == 103 || fIcode == 400) return kIpPDeltaRay;
2419     else if (fIcode == 210 || fIcode == 212) return kIpPDeltaRay;
2420 //  else if (fIcode == 210) return kIpPMoller;
2421 //  else if (fIcode == 212) return kIpPBhabha;
2422     else if (fIcode == 214 || fIcode == 215) return kIpPAnnihilation;
2423 //  else if (fIcode == 214) return kIpPAnnihilInFlight;
2424 //  else if (fIcode == 215) return kIpPAnnihilAtRest;
2425     else if (fIcode == 101) return kIpPHadronic;
2426     else if (fIcode == 101) {
2427       if (!mugamma) return kIpPHadronic;
2428       else if (TRACKR.jtrack == 7) return kIpPPhotoFission;
2429       else return kIpPMuonNuclear;
2430     }
2431     else if (fIcode == 225) return kIpPRayleigh;
2432 // Fluka codes 100, 300 and 400 still to be investigasted
2433     else return kIpNoProc;
2434 }
2435
2436 //Int_t StepProcesses(TArrayI &proc) const
2437 // Return processes active in the current step
2438 //{
2439 //ck = total energy of the particl ???????????????? 
2440 //}
2441
2442
2443 Int_t TFluka::VolId2Mate(Int_t id) const
2444 {
2445 //
2446 // Returns the material number for a given volume ID
2447 //
2448     if (fVerbosityLevel >= 3)
2449     printf("VolId2Mate %d %d\n", id, fMediaByRegion[id-1]); 
2450     return fMediaByRegion[id-1];
2451 }
2452
2453 const char* TFluka::VolName(Int_t id) const
2454 {
2455 //
2456 // Returns the volume name for a given volume ID
2457 //
2458     FlukaVolume* vol = dynamic_cast<FlukaVolume*>((*fVolumeMediaMap)[id-1]);
2459     const char* name = vol->GetName();
2460     if (fVerbosityLevel >= 3)
2461     printf("VolName %d %s \n", id, name);
2462     return name;
2463 }
2464
2465 Int_t TFluka::VolId(const Text_t* volName) const
2466 {
2467 //
2468 // Converts from volume name to volume ID.
2469 // Time consuming. (Only used during set-up)
2470 // Could be replaced by hash-table
2471 //
2472     char tmp[5];
2473     Int_t i =0;
2474     for (i = 0; i < fNVolumes; i++)
2475   {
2476       FlukaVolume* vol = dynamic_cast<FlukaVolume*>((*fVolumeMediaMap)[i]);
2477       TString name = vol->GetName();
2478       strcpy(tmp, name.Data());
2479       tmp[4] = '\0';
2480       if (!strcmp(tmp, volName)) break;
2481   }
2482     i++;
2483
2484     return i;
2485 }
2486
2487
2488 Int_t TFluka::CurrentVolID(Int_t& copyNo) const
2489 {
2490 //
2491 // Return the logical id and copy number corresponding to the current fluka region
2492 //
2493     int ir = fCurrentFlukaRegion;
2494     int id = (FGeometryInit::GetInstance())->CurrentVolID(ir, copyNo);
2495     copyNo++;
2496     if (fVerbosityLevel >= 3)
2497     printf("CurrentVolID: %d %d %d \n", ir, id, copyNo); 
2498     return id;
2499
2500
2501 Int_t TFluka::CurrentVolOffID(Int_t off, Int_t& copyNo) const
2502 {
2503 //
2504 // Return the logical id and copy number of off'th mother 
2505 // corresponding to the current fluka region
2506 //
2507     if (off == 0) 
2508         return CurrentVolID(copyNo);
2509     
2510     int ir = fCurrentFlukaRegion;
2511     int id = (FGeometryInit::GetInstance())->CurrentVolOffID(ir, off, copyNo);
2512     copyNo++;
2513     if (fVerbosityLevel >= 3)
2514     printf("CurrentVolOffID: %d %d %d \n", ir, id, copyNo); 
2515     if (id == -1) 
2516         if (fVerbosityLevel >= 0)
2517         printf("CurrentVolOffID: Warning Mother not found !!!\n"); 
2518     return id;
2519 }
2520
2521
2522 const char* TFluka::CurrentVolName() const
2523 {
2524 //
2525 // Return the current volume name
2526 //
2527     Int_t copy;
2528     Int_t id = TFluka::CurrentVolID(copy);
2529     const char* name = TFluka::VolName(id);
2530     if (fVerbosityLevel >= 3)
2531     printf("CurrentVolumeName: %d %s \n", fCurrentFlukaRegion,  name); 
2532     return name;
2533 }
2534
2535 const char* TFluka::CurrentVolOffName(Int_t off) const
2536 {
2537 //
2538 // Return the volume name of the off'th mother of the current volume
2539 //
2540     Int_t copy;
2541     Int_t id = TFluka::CurrentVolOffID(off, copy);
2542     const char* name = TFluka::VolName(id);
2543     if (fVerbosityLevel >= 3)
2544     printf("CurrentVolumeOffName: %d %s \n", fCurrentFlukaRegion,  name); 
2545     return name;
2546 }
2547
2548 Int_t TFluka::CurrentMaterial(Float_t & /*a*/, Float_t & /*z*/, 
2549                       Float_t & /*dens*/, Float_t & /*radl*/, Float_t & /*absl*/) const
2550 {
2551 //
2552 //  Return the current medium number
2553 //
2554     Int_t copy;
2555     Int_t id  =  TFluka::CurrentVolID(copy);
2556     Int_t med =  TFluka::VolId2Mate(id);
2557     if (fVerbosityLevel >= 3)
2558     printf("CurrentMaterial: %d %d \n", fCurrentFlukaRegion,  med); 
2559     return med;
2560 }
2561
2562 void TFluka::Gmtod(Float_t* xm, Float_t* xd, Int_t iflag)
2563     {
2564 // Transforms a position from the world reference frame
2565 // to the current volume reference frame.
2566 //
2567 //  Geant3 desription:
2568 //  ==================
2569 //       Computes coordinates XD (in DRS) 
2570 //       from known coordinates XM in MRS 
2571 //       The local reference system can be initialized by
2572 //         - the tracking routines and GMTOD used in GUSTEP
2573 //         - a call to GMEDIA(XM,NUMED)
2574 //         - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER) 
2575 //             (inverse routine is GDTOM) 
2576 //
2577 //        If IFLAG=1  convert coordinates 
2578 //           IFLAG=2  convert direction cosinus
2579 //
2580 // ---
2581         Double_t xmD[3], xdD[3];        
2582         xmD[0] = xm[0]; xmD[1] = xm[1]; xmD[2] = xm[2]; 
2583         (FGeometryInit::GetInstance())->Gmtod(xmD, xdD, iflag);
2584         xd[0] = xdD[0]; xd[1] = xdD[1]; xd[2] = xdD[2]; 
2585     }
2586
2587   
2588 void TFluka::Gmtod(Double_t* xm, Double_t* xd, Int_t iflag)
2589     {
2590 // Transforms a position from the world reference frame
2591 // to the current volume reference frame.
2592 //
2593 //  Geant3 desription:
2594 //  ==================
2595 //       Computes coordinates XD (in DRS) 
2596 //       from known coordinates XM in MRS 
2597 //       The local reference system can be initialized by
2598 //         - the tracking routines and GMTOD used in GUSTEP
2599 //         - a call to GMEDIA(XM,NUMED)
2600 //         - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER) 
2601 //             (inverse routine is GDTOM) 
2602 //
2603 //        If IFLAG=1  convert coordinates 
2604 //           IFLAG=2  convert direction cosinus
2605 //
2606 // ---
2607         (FGeometryInit::GetInstance())->Gmtod(xm, xd, iflag);
2608     }
2609
2610 void TFluka::Gdtom(Float_t* xd, Float_t* xm, Int_t iflag)
2611     {
2612 // Transforms a position from the current volume reference frame
2613 // to the world reference frame.
2614 //
2615 //  Geant3 desription:
2616 //  ==================
2617 //  Computes coordinates XM (Master Reference System
2618 //  knowing the coordinates XD (Detector Ref System)
2619 //  The local reference system can be initialized by
2620 //    - the tracking routines and GDTOM used in GUSTEP
2621 //    - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2622 //        (inverse routine is GMTOD)
2623 // 
2624 //   If IFLAG=1  convert coordinates
2625 //      IFLAG=2  convert direction cosinus
2626 //
2627 // ---
2628         Double_t xmD[3], xdD[3];        
2629         xdD[0] = xd[0]; xdD[1] = xd[1]; xdD[2] = xd[2]; 
2630         (FGeometryInit::GetInstance())->Gdtom(xdD, xmD, iflag);
2631         xm[0] = xmD[0]; xm[1] = xmD[1]; xm[2] = xmD[2]; 
2632     }
2633 void TFluka::Gdtom(Double_t* xd, Double_t* xm, Int_t iflag)
2634     {
2635 // Transforms a position from the current volume reference frame
2636 // to the world reference frame.
2637 //
2638 //  Geant3 desription:
2639 //  ==================
2640 //  Computes coordinates XM (Master Reference System
2641 //  knowing the coordinates XD (Detector Ref System)
2642 //  The local reference system can be initialized by
2643 //    - the tracking routines and GDTOM used in GUSTEP
2644 //    - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2645 //        (inverse routine is GMTOD)
2646 // 
2647 //   If IFLAG=1  convert coordinates
2648 //      IFLAG=2  convert direction cosinus
2649 //
2650 // ---
2651
2652         (FGeometryInit::GetInstance())->Gdtom(xd, xm, iflag);
2653     }
2654
2655 // ===============================================================
2656 void TFluka::FutoTest() 
2657 {
2658     Int_t icode, mreg, newreg, particleId;
2659     Double_t rull, xsco, ysco, zsco;
2660     TLorentzVector position, momentum;
2661     icode = GetIcode();
2662     if (icode == 0) {
2663         if (fVerbosityLevel >=3)
2664             cout << " icode=" << icode << endl;
2665     } else if (icode > 0 && icode <= 5) {
2666 // mgdraw
2667         mreg = GetMreg();
2668         if (fVerbosityLevel >=3)
2669             cout << " icode=" << icode
2670                  << " mreg=" << mreg
2671                  << endl;
2672         TrackPosition(position);
2673         TrackMomentum(momentum);
2674         if (fVerbosityLevel >=3) {
2675             cout << "TLorentzVector positionX=" << position.X()
2676                  << "positionY=" << position.Y()
2677                  << "positionZ=" << position.Z()
2678                  << "timeT=" << position.T() << endl;
2679             cout << "TLorentzVector momentumX=" << momentum.X()
2680                  << "momentumY=" << momentum.Y()
2681                  << "momentumZ=" << momentum.Z()
2682                  << "energyE=" << momentum.E() << endl;
2683             cout << "TrackStep=" << TrackStep() << endl;
2684             cout << "TrackLength=" << TrackLength() << endl;
2685             cout << "TrackTime=" << TrackTime() << endl;
2686             cout << "Edep=" << Edep() << endl;
2687             cout << "TrackPid=" << TrackPid() << endl;
2688             cout << "TrackCharge=" << TrackCharge() << endl;
2689             cout << "TrackMass=" << TrackMass() << endl;
2690             cout << "Etot=" << Etot() << endl;
2691             cout << "IsNewTrack=" << IsNewTrack() << endl;
2692             cout << "IsTrackInside=" << IsTrackInside() << endl;
2693             cout << "IsTrackEntering=" << IsTrackEntering() << endl;
2694             cout << "IsTrackExiting=" << IsTrackExiting() << endl;
2695             cout << "IsTrackOut=" << IsTrackOut() << endl;
2696             cout << "IsTrackDisappeared=" << IsTrackDisappeared() << endl;
2697             cout << "IsTrackAlive=" << IsTrackAlive() << endl;
2698         }
2699         
2700         Float_t x = position.X();
2701         Float_t y = position.Y();
2702         Float_t z = position.Z();
2703         Float_t xm[3];
2704         Float_t xd[3];
2705         xm[0] = x; xm[1] = y; xm[2] = z;
2706         if (fVerbosityLevel >= 3)
2707             printf("Global trackPosition: %f %f %f \n", x, y, z);
2708         Gmtod(xm, xd, 1);
2709         if (fVerbosityLevel >= 3)
2710             printf("Local trackPosition: %f %f %f \n", xd[0], xd[1], xd[2]);
2711         Gdtom(xd, xm, 1);
2712         if (fVerbosityLevel >= 3)
2713             printf("New trackPosition: %f %f %f \n", xm[0], xm[1], xm[2]);
2714     } else if((icode >= 10 && icode <= 15) ||
2715               (icode >= 20 && icode <= 24) ||
2716               (icode >= 30 && icode <= 33) ||
2717               (icode >= 40 && icode <= 41) ||
2718               (icode >= 50 && icode <= 52)) {
2719 // endraw
2720         mreg = GetMreg();
2721         rull = GetRull();
2722         xsco = GetXsco();
2723         ysco = GetYsco();
2724         zsco = GetZsco();
2725         
2726         if (fVerbosityLevel >=3) {     
2727             cout << " icode=" << icode
2728                  << " mreg=" << mreg
2729                  << " rull=" << rull
2730                  << " xsco=" << xsco
2731                  << " ysco=" << ysco
2732                  << " zsco=" << zsco << endl;
2733         }
2734         TrackPosition(position);
2735         TrackMomentum(momentum);
2736         if (fVerbosityLevel >=3) {
2737             cout << "Edep=" << Edep() << endl;
2738             cout << "Etot=" << Etot() << endl;
2739             cout << "TrackPid=" << TrackPid() << endl;
2740             cout << "TrackCharge=" << TrackCharge() << endl;
2741             cout << "TrackMass=" << TrackMass() << endl;
2742             cout << "IsTrackOut=" << IsTrackOut() << endl;
2743             cout << "IsTrackDisappeared=" << IsTrackDisappeared() << endl;
2744             cout << "IsTrackStop=" << IsTrackStop() << endl;
2745             cout << "IsTrackAlive=" << IsTrackAlive() << endl;
2746         }
2747     } else if((icode >= 100 && icode <= 105) ||
2748               (icode == 208) ||
2749               (icode == 210) ||
2750               (icode == 212) ||
2751               (icode >= 214 && icode <= 215) ||
2752               (icode == 217) ||
2753               (icode == 219) ||
2754               (icode == 221) ||
2755               (icode == 225) ||
2756               (icode == 300) ||
2757               (icode == 400)) {
2758 // usdraw
2759           mreg = GetMreg();
2760           xsco = GetXsco();
2761           ysco = GetYsco();
2762           zsco = GetZsco();
2763           
2764           if (fVerbosityLevel >=3) {
2765               cout << " icode=" << icode
2766                    << " mreg=" << mreg
2767                    << " xsco=" << xsco
2768                    << " ysco=" << ysco
2769                    << " zsco=" << zsco << endl;
2770               cout << "TrackPid=" << TrackPid() << endl;
2771               cout << "NSecondaries=" << NSecondaries() << endl;
2772           }
2773           
2774           for (Int_t isec=0; isec< NSecondaries(); isec++) {
2775               TFluka::GetSecondary(isec, particleId, position, momentum);
2776               if (fVerbosityLevel >=3) {
2777                   cout << "TLorentzVector positionX=" << position.X()
2778                        << "positionY=" << position.Y()
2779                        << "positionZ=" << position.Z()
2780                        << "timeT=" << position.T() << endl;
2781                   cout << "TLorentzVector momentumX=" << momentum.X()
2782                        << "momentumY=" << momentum.Y()
2783                        << "momentumZ=" << momentum.Z()
2784                        << "energyE=" << momentum.E() << endl;
2785                   cout << "TrackPid=" << particleId << endl;
2786               }
2787           }
2788     } else if((icode == 19) ||
2789               (icode == 29) ||
2790               (icode == 39) ||
2791               (icode == 49) ||
2792               (icode == 59)) {
2793         mreg = GetMreg();
2794         newreg = GetNewreg();
2795         xsco = GetXsco();
2796         ysco = GetYsco();
2797         zsco = GetZsco();
2798         if (fVerbosityLevel >=3) {
2799             cout << " icode=" << icode
2800                  << " mreg=" << mreg
2801                  << " newreg=" << newreg
2802                  << " xsco=" << xsco
2803                  << " ysco=" << ysco
2804                  << " zsco=" << zsco << endl;
2805         }
2806     }
2807 } // end of FutoTest
2808