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