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