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