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