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