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