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