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