1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 // Realisation of the TVirtualMC interface for the FLUKA code
20 // (See official web side http://www.fluka.org/).
22 // This implementation makes use of the TGeo geometry modeller.
23 // User configuration is via automatic generation of FLUKA input cards.
32 #include <Riostream.h>
35 #include "TCallf77.h" //For the fortran calls
36 #include "Fdblprc.h" //(DBLPRC) fluka common
37 #include "Fepisor.h" //(EPISOR) fluka common
38 #include "Ffinuc.h" //(FINUC) fluka common
39 #include "Fiounit.h" //(IOUNIT) fluka common
40 #include "Fpaprop.h" //(PAPROP) fluka common
41 #include "Fpart.h" //(PART) fluka common
42 #include "Ftrackr.h" //(TRACKR) fluka common
43 #include "Fpaprop.h" //(PAPROP) fluka common
44 #include "Ffheavy.h" //(FHEAVY) fluka common
45 #include "Fopphst.h" //(OPPHST) fluka common
46 #include "Fstack.h" //(STACK) fluka common
47 #include "Fstepsz.h" //(STEPSZ) fluka common
48 #include "Fopphst.h" //(OPPHST) fluka common
50 #include "TVirtualMC.h"
51 #include "TMCProcess.h"
52 #include "TGeoManager.h"
53 #include "TGeoMaterial.h"
54 #include "TGeoMedium.h"
55 #include "TFlukaMCGeometry.h"
56 #include "TGeoMCGeometry.h"
57 #include "TFlukaCerenkov.h"
58 #include "TFlukaConfigOption.h"
59 #include "TFlukaScoringOption.h"
60 #include "TLorentzVector.h"
63 // Fluka methods that may be needed.
65 # define flukam flukam_
66 # define fluka_openinp fluka_openinp_
67 # define fluka_closeinp fluka_closeinp_
68 # define mcihad mcihad_
69 # define mpdgha mpdgha_
71 # define flukam FLUKAM
72 # define fluka_openinp FLUKA_OPENINP
73 # define fluka_closeinp FLUKA_CLOSEINP
74 # define mcihad MCIHAD
75 # define mpdgha MPDGHA
81 // Prototypes for FLUKA functions
83 void type_of_call flukam(const int&);
84 void type_of_call fluka_openinp(const int&, DEFCHARA);
85 void type_of_call fluka_closeinp(const int&);
86 int type_of_call mcihad(const int&);
87 int type_of_call mpdgha(const int&);
91 // Class implementation for ROOT
96 //----------------------------------------------------------------------------
97 // TFluka constructors and destructors.
98 //______________________________________________________________________________
108 // Default constructor
110 fGeneratePemf = kFALSE;
112 fCurrentFlukaRegion = -1;
124 //______________________________________________________________________________
125 TFluka::TFluka(const char *title, Int_t verbosity, Bool_t isRootGeometrySupported)
126 :TVirtualMC("TFluka",title, isRootGeometrySupported),
127 fVerbosityLevel(verbosity),
132 fProcesses(new TObjArray(100)),
133 fCuts(new TObjArray(100)),
134 fUserScore(new TObjArray(100))
136 // create geometry interface
137 if (fVerbosityLevel >=3)
138 cout << "<== TFluka::TFluka(" << title << ") constructor called." << endl;
139 SetCoreInputFileName();
141 SetGeneratePemf(kFALSE);
143 fCurrentFlukaRegion = -1;
146 fGeneratePemf = kFALSE;
147 fMCGeo = new TGeoMCGeometry("MCGeo", "TGeo Implementation of VirtualMCGeometry", kTRUE);
148 fGeom = new TFlukaMCGeometry("geom", "FLUKA VMC Geometry");
149 if (verbosity > 2) fGeom->SetDebugMode(kTRUE);
157 //______________________________________________________________________________
160 if (fVerbosityLevel >=3)
161 cout << "<== TFluka::~TFluka() destructor called." << endl;
172 fProcesses->Delete();
180 //______________________________________________________________________________
181 // TFluka control methods
182 //______________________________________________________________________________
183 void TFluka::Init() {
185 // Geometry initialisation
187 if (fVerbosityLevel >=3) cout << "==> TFluka::Init() called." << endl;
189 if (!gGeoManager) new TGeoManager("geom", "FLUKA geometry");
190 fApplication->ConstructGeometry();
191 TGeoVolume *top = (TGeoVolume*)gGeoManager->GetListOfVolumes()->First();
192 gGeoManager->SetTopVolume(top);
193 gGeoManager->CloseGeometry("di");
194 gGeoManager->DefaultColors(); // to be removed
195 fNVolumes = fGeom->NofVolumes();
196 fGeom->CreateFlukaMatFile("flukaMat.inp");
197 if (fVerbosityLevel >=3) {
198 printf("== Number of volumes: %i\n ==", fNVolumes);
199 cout << "\t* InitPhysics() - Prepare input file to be called" << endl;
201 // now we have TGeo geometry created and we have to patch FlukaVmc.inp
202 // with the material mapping file FlukaMat.inp
206 //______________________________________________________________________________
207 void TFluka::FinishGeometry() {
209 // Build-up table with region to medium correspondance
211 if (fVerbosityLevel >=3) {
212 cout << "==> TFluka::FinishGeometry() called." << endl;
213 printf("----FinishGeometry - nothing to do with TGeo\n");
214 cout << "<== TFluka::FinishGeometry() called." << endl;
218 //______________________________________________________________________________
219 void TFluka::BuildPhysics() {
221 // Prepare FLUKA input files and call FLUKA physics initialisation
224 if (fVerbosityLevel >=3)
225 cout << "==> TFluka::BuildPhysics() called." << endl;
226 // Prepare input file with the current physics settings
228 cout << "\t* InitPhysics() - Prepare input file was called" << endl;
230 if (fVerbosityLevel >=2)
231 cout << "\t* Changing lfdrtr = (" << (GLOBAL.lfdrtr?'T':'F')
232 << ") in fluka..." << endl;
233 GLOBAL.lfdrtr = true;
235 if (fVerbosityLevel >=2)
236 cout << "\t* Opening file " << fInputFileName << endl;
237 const char* fname = fInputFileName;
238 fluka_openinp(lunin, PASSCHARA(fname));
240 if (fVerbosityLevel >=2)
241 cout << "\t* Calling flukam..." << endl;
244 if (fVerbosityLevel >=2)
245 cout << "\t* Closing file " << fInputFileName << endl;
246 fluka_closeinp(lunin);
250 if (fVerbosityLevel >=3)
251 cout << "<== TFluka::Init() called." << endl;
254 if (fVerbosityLevel >=3)
255 cout << "<== TFluka::BuildPhysics() called." << endl;
258 //______________________________________________________________________________
259 void TFluka::ProcessEvent() {
264 printf("User Run Abortion: No more events handled !\n");
269 if (fVerbosityLevel >=3)
270 cout << "==> TFluka::ProcessEvent() called." << endl;
271 fApplication->GeneratePrimaries();
272 EPISOR.lsouit = true;
274 if (fVerbosityLevel >=3)
275 cout << "<== TFluka::ProcessEvent() called." << endl;
277 // Increase event number
282 //______________________________________________________________________________
283 Bool_t TFluka::ProcessRun(Int_t nevent) {
288 if (fVerbosityLevel >=3)
289 cout << "==> TFluka::ProcessRun(" << nevent << ") called."
292 if (fVerbosityLevel >=2) {
293 cout << "\t* GLOBAL.fdrtr = " << (GLOBAL.lfdrtr?'T':'F') << endl;
294 cout << "\t* Calling flukam again..." << endl;
297 fApplication->InitGeometry();
298 Int_t todo = TMath::Abs(nevent);
299 for (Int_t ev = 0; ev < todo; ev++) {
300 fApplication->BeginEvent();
302 fApplication->FinishEvent();
305 if (fVerbosityLevel >=3)
306 cout << "<== TFluka::ProcessRun(" << nevent << ") called."
311 //_____________________________________________________________________________
312 // methods for building/management of geometry
314 // functions from GCONS
315 //____________________________________________________________________________
316 void TFluka::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
317 Float_t &dens, Float_t &radl, Float_t &absl,
318 Float_t* /*ubuf*/, Int_t& /*nbuf*/) {
321 TIter next (gGeoManager->GetListOfMaterials());
322 while ((mat = (TGeoMaterial*)next())) {
323 if (mat->GetUniqueID() == (UInt_t)imat) break;
326 Error("Gfmate", "no material with index %i found", imat);
329 sprintf(name, "%s", mat->GetName());
332 dens = mat->GetDensity();
333 radl = mat->GetRadLen();
334 absl = mat->GetIntLen();
337 //______________________________________________________________________________
338 void TFluka::Gfmate(Int_t imat, char *name, Double_t &a, Double_t &z,
339 Double_t &dens, Double_t &radl, Double_t &absl,
340 Double_t* /*ubuf*/, Int_t& /*nbuf*/) {
343 TIter next (gGeoManager->GetListOfMaterials());
344 while ((mat = (TGeoMaterial*)next())) {
345 if (mat->GetUniqueID() == (UInt_t)imat) break;
348 Error("Gfmate", "no material with index %i found", imat);
351 sprintf(name, "%s", mat->GetName());
354 dens = mat->GetDensity();
355 radl = mat->GetRadLen();
356 absl = mat->GetIntLen();
359 // detector composition
360 //______________________________________________________________________________
361 void TFluka::Material(Int_t& kmat, const char* name, Double_t a,
362 Double_t z, Double_t dens, Double_t radl, Double_t absl,
363 Float_t* buf, Int_t nwbuf) {
365 Double_t* dbuf = fGeom->CreateDoubleArray(buf, nwbuf);
366 Material(kmat, name, a, z, dens, radl, absl, dbuf, nwbuf);
370 //______________________________________________________________________________
371 void TFluka::Material(Int_t& kmat, const char* name, Double_t a,
372 Double_t z, Double_t dens, Double_t radl, Double_t absl,
373 Double_t* /*buf*/, Int_t /*nwbuf*/) {
376 kmat = gGeoManager->GetListOfMaterials()->GetSize();
377 if ((z-Int_t(z)) > 1E-3) {
378 mat = fGeom->GetMakeWrongMaterial(z);
380 mat->SetRadLen(radl,absl);
381 mat->SetUniqueID(kmat);
385 gGeoManager->Material(name, a, z, dens, kmat, radl, absl);
388 //______________________________________________________________________________
389 void TFluka::Mixture(Int_t& kmat, const char *name, Float_t *a,
390 Float_t *z, Double_t dens, Int_t nlmat, Float_t *wmat) {
392 Double_t* da = fGeom->CreateDoubleArray(a, TMath::Abs(nlmat));
393 Double_t* dz = fGeom->CreateDoubleArray(z, TMath::Abs(nlmat));
394 Double_t* dwmat = fGeom->CreateDoubleArray(wmat, TMath::Abs(nlmat));
396 Mixture(kmat, name, da, dz, dens, nlmat, dwmat);
397 for (Int_t i=0; i<nlmat; i++) {
398 a[i] = da[i]; z[i] = dz[i]; wmat[i] = dwmat[i];
406 //______________________________________________________________________________
407 void TFluka::Mixture(Int_t& kmat, const char *name, Double_t *a,
408 Double_t *z, Double_t dens, Int_t nlmat, Double_t *wmat) {
410 // Defines mixture OR COMPOUND IMAT as composed by
411 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
413 // If NLMAT > 0 then wmat contains the proportion by
414 // weights of each basic material in the mixture.
416 // If nlmat < 0 then WMAT contains the number of atoms
417 // of a given kind into the molecule of the COMPOUND
418 // In this case, WMAT in output is changed to relative
425 for (i=0;i<nlmat;i++) {
426 amol += a[i]*wmat[i];
428 for (i=0;i<nlmat;i++) {
429 wmat[i] *= a[i]/amol;
432 kmat = gGeoManager->GetListOfMaterials()->GetSize();
433 // Check if we have elements with fractional Z
434 TGeoMaterial *mat = 0;
435 TGeoMixture *mix = 0;
436 Bool_t mixnew = kFALSE;
437 for (i=0; i<nlmat; i++) {
438 if (z[i]-Int_t(z[i]) < 1E-3) continue;
439 // We have found an element with fractional Z -> loop mixtures to look for it
440 for (j=0; j<kmat; j++) {
441 mat = (TGeoMaterial*)gGeoManager->GetListOfMaterials()->At(j);
443 if (!mat->IsMixture()) continue;
444 mix = (TGeoMixture*)mat;
445 if (TMath::Abs(z[i]-mix->GetZ()) >1E-3) continue;
446 // printf(" FOUND component %i as mixture %s\n", i, mat->GetName());
450 if (!mixnew) Warning("Mixture","%s : cannot find component %i with fractional Z=%f\n", name, i, z[i]);
454 Int_t nlmatnew = nlmat+mix->GetNelements()-1;
455 Double_t *anew = new Double_t[nlmatnew];
456 Double_t *znew = new Double_t[nlmatnew];
457 Double_t *wmatnew = new Double_t[nlmatnew];
459 for (j=0; j<nlmat; j++) {
463 wmatnew[ind] = wmat[j];
466 for (j=0; j<mix->GetNelements(); j++) {
467 anew[ind] = mix->GetAmixt()[j];
468 znew[ind] = mix->GetZmixt()[j];
469 wmatnew[ind] = wmat[i]*mix->GetWmixt()[j];
472 Mixture(kmat, name, anew, znew, dens, nlmatnew, wmatnew);
478 // Now we need to compact identical elements within the mixture
479 // First check if this happens
481 for (i=0; i<nlmat-1; i++) {
482 for (j=i+1; j<nlmat; j++) {
492 Double_t *anew = new Double_t[nlmat];
493 Double_t *znew = new Double_t[nlmat];
494 memset(znew, 0, nlmat*sizeof(Double_t));
495 Double_t *wmatnew = new Double_t[nlmat];
497 for (i=0; i<nlmat; i++) {
499 for (j=0; j<nlmatnew; j++) {
501 wmatnew[j] += wmat[i];
507 anew[nlmatnew] = a[i];
508 znew[nlmatnew] = z[i];
509 wmatnew[nlmatnew] = wmat[i];
512 Mixture(kmat, name, anew, znew, dens, nlmatnew, wmatnew);
518 gGeoManager->Mixture(name, a, z, dens, nlmat, wmat, kmat);
521 //______________________________________________________________________________
522 void TFluka::Medium(Int_t& kmed, const char *name, Int_t nmat,
523 Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd,
524 Double_t stemax, Double_t deemax, Double_t epsil,
525 Double_t stmin, Float_t* ubuf, Int_t nbuf) {
528 kmed = gGeoManager->GetListOfMedia()->GetSize()+1;
529 fMCGeo->Medium(kmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax,
530 epsil, stmin, ubuf, nbuf);
533 //______________________________________________________________________________
534 void TFluka::Medium(Int_t& kmed, const char *name, Int_t nmat,
535 Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd,
536 Double_t stemax, Double_t deemax, Double_t epsil,
537 Double_t stmin, Double_t* ubuf, Int_t nbuf) {
540 kmed = gGeoManager->GetListOfMedia()->GetSize()+1;
541 fMCGeo->Medium(kmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax,
542 epsil, stmin, ubuf, nbuf);
545 //______________________________________________________________________________
546 void TFluka::Matrix(Int_t& krot, Double_t thetaX, Double_t phiX,
547 Double_t thetaY, Double_t phiY, Double_t thetaZ,
550 krot = gGeoManager->GetListOfMatrices()->GetEntriesFast();
551 fMCGeo->Matrix(krot, thetaX, phiX, thetaY, phiY, thetaZ, phiZ);
554 //______________________________________________________________________________
555 void TFluka::Gstpar(Int_t itmed, const char* param, Double_t parval) {
558 // Check if material is used
559 if (fVerbosityLevel >=3)
560 printf("Gstpar called with %6d %5s %12.4e %6d\n", itmed, param, parval, fGeom->GetFlukaMaterial(itmed));
563 reglist = fGeom->GetMaterialList(fGeom->GetFlukaMaterial(itmed), nreg);
569 Bool_t process = kFALSE;
570 if (strncmp(param, "DCAY", 4) == 0 ||
571 strncmp(param, "PAIR", 4) == 0 ||
572 strncmp(param, "COMP", 4) == 0 ||
573 strncmp(param, "PHOT", 4) == 0 ||
574 strncmp(param, "PFIS", 4) == 0 ||
575 strncmp(param, "DRAY", 4) == 0 ||
576 strncmp(param, "ANNI", 4) == 0 ||
577 strncmp(param, "BREM", 4) == 0 ||
578 strncmp(param, "MUNU", 4) == 0 ||
579 strncmp(param, "CKOV", 4) == 0 ||
580 strncmp(param, "HADR", 4) == 0 ||
581 strncmp(param, "LOSS", 4) == 0 ||
582 strncmp(param, "MULS", 4) == 0 ||
583 strncmp(param, "RAYL", 4) == 0)
588 SetProcess(param, Int_t (parval), fGeom->GetFlukaMaterial(itmed));
590 SetCut(param, parval, fGeom->GetFlukaMaterial(itmed));
594 // functions from GGEOM
595 //_____________________________________________________________________________
596 void TFluka::Gsatt(const char *name, const char *att, Int_t val)
598 // Set visualisation attributes for one volume
600 fGeom->Vname(name,vname);
602 fGeom->Vname(att,vatt);
603 gGeoManager->SetVolumeAttribute(vname, vatt, val);
606 //______________________________________________________________________________
607 Int_t TFluka::Gsvolu(const char *name, const char *shape, Int_t nmed,
608 Float_t *upar, Int_t np) {
610 return fMCGeo->Gsvolu(name, shape, nmed, upar, np);
613 //______________________________________________________________________________
614 Int_t TFluka::Gsvolu(const char *name, const char *shape, Int_t nmed,
615 Double_t *upar, Int_t np) {
617 return fMCGeo->Gsvolu(name, shape, nmed, upar, np);
620 //______________________________________________________________________________
621 void TFluka::Gsdvn(const char *name, const char *mother, Int_t ndiv,
624 fMCGeo->Gsdvn(name, mother, ndiv, iaxis);
627 //______________________________________________________________________________
628 void TFluka::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
629 Int_t iaxis, Double_t c0i, Int_t numed) {
631 fMCGeo->Gsdvn2(name, mother, ndiv, iaxis, c0i, numed);
634 //______________________________________________________________________________
635 void TFluka::Gsdvt(const char *name, const char *mother, Double_t step,
636 Int_t iaxis, Int_t numed, Int_t ndvmx) {
638 fMCGeo->Gsdvt(name, mother, step, iaxis, numed, ndvmx);
641 //______________________________________________________________________________
642 void TFluka::Gsdvt2(const char *name, const char *mother, Double_t step,
643 Int_t iaxis, Double_t c0, Int_t numed, Int_t ndvmx) {
645 fMCGeo->Gsdvt2(name, mother, step, iaxis, c0, numed, ndvmx);
648 //______________________________________________________________________________
649 void TFluka::Gsord(const char * /*name*/, Int_t /*iax*/) {
651 // Nothing to do with TGeo
654 //______________________________________________________________________________
655 void TFluka::Gspos(const char *name, Int_t nr, const char *mother,
656 Double_t x, Double_t y, Double_t z, Int_t irot,
659 fMCGeo->Gspos(name, nr, mother, x, y, z, irot, konly);
662 //______________________________________________________________________________
663 void TFluka::Gsposp(const char *name, Int_t nr, const char *mother,
664 Double_t x, Double_t y, Double_t z, Int_t irot,
665 const char *konly, Float_t *upar, Int_t np) {
667 fMCGeo->Gsposp(name, nr, mother, x, y, z, irot, konly, upar, np);
670 //______________________________________________________________________________
671 void TFluka::Gsposp(const char *name, Int_t nr, const char *mother,
672 Double_t x, Double_t y, Double_t z, Int_t irot,
673 const char *konly, Double_t *upar, Int_t np) {
675 fMCGeo->Gsposp(name, nr, mother, x, y, z, irot, konly, upar, np);
678 //______________________________________________________________________________
679 void TFluka::Gsbool(const char* /*onlyVolName*/, const char* /*manyVolName*/) {
681 // Nothing to do with TGeo
684 //______________________________________________________________________________
685 void TFluka::SetCerenkov(Int_t itmed, Int_t npckov, Float_t* ppckov,
686 Float_t* absco, Float_t* effic, Float_t* rindex) {
688 // Set Cerenkov properties for medium itmed
690 // npckov: number of sampling points
691 // ppckov: energy values
692 // absco: absorption length
693 // effic: quantum efficiency
694 // rindex: refraction index
698 // Create object holding Cerenkov properties
700 TFlukaCerenkov* cerenkovProperties = new TFlukaCerenkov(npckov, ppckov, absco, effic, rindex);
702 // Pass object to medium
703 TGeoMedium* medium = gGeoManager->GetMedium(itmed);
704 medium->SetCerenkovProperties(cerenkovProperties);
707 //______________________________________________________________________________
708 void TFluka::SetCerenkov(Int_t /*itmed*/, Int_t /*npckov*/, Double_t * /*ppckov*/,
709 Double_t * /*absco*/, Double_t * /*effic*/, Double_t * /*rindex*/) {
711 // Not implemented with TGeo - what G4 did ? Any FLUKA card generated?
712 Warning("SetCerenkov", "Not implemented with TGeo");
716 //______________________________________________________________________________
717 void TFluka::WriteEuclid(const char* /*fileName*/, const char* /*topVol*/,
718 Int_t /*number*/, Int_t /*nlevel*/) {
721 Warning("WriteEuclid", "Not implemented with TGeo");
726 //_____________________________________________________________________________
727 // methods needed by the stepping
728 //____________________________________________________________________________
730 Int_t TFluka::GetMedium() const {
732 // Get the medium number for the current fluka region
734 return fGeom->GetMedium(); // this I need to check due to remapping !!!
739 //____________________________________________________________________________
740 // particle table usage
741 // ID <--> PDG transformations
742 //_____________________________________________________________________________
743 Int_t TFluka::IdFromPDG(Int_t pdg) const
746 // Return Fluka code from PDG and pseudo ENDF code
748 // Catch the feedback photons
749 if (pdg == 50000051) return (-1);
750 // MCIHAD() goes from pdg to fluka internal.
751 Int_t intfluka = mcihad(pdg);
752 // KPTOIP array goes from internal to official
753 return GetFlukaKPTOIP(intfluka);
756 //______________________________________________________________________________
757 Int_t TFluka::PDGFromId(Int_t id) const
760 // Return PDG code and pseudo ENDF code from Fluka code
761 // Alpha He3 Triton Deuteron gen. ion opt. photon
762 Int_t idSpecial[6] = {10020040, 10020030, 10010030, 10010020, 10000000, 50000050};
763 // IPTOKP array goes from official to internal
767 if (fVerbosityLevel >= 3)
768 printf("\n PDGFromId: Cerenkov Photon \n");
772 if (id == 0 || id < -6 || id > 250) {
773 if (fVerbosityLevel >= 3)
774 printf("PDGFromId: Error id = 0\n");
779 Int_t intfluka = GetFlukaIPTOKP(id);
781 if (fVerbosityLevel >= 3)
782 printf("PDGFromId: Error intfluka = 0: %d\n", id);
784 } else if (intfluka < 0) {
785 if (fVerbosityLevel >= 3)
786 printf("PDGFromId: Error intfluka < 0: %d\n", id);
789 if (fVerbosityLevel >= 3)
790 printf("mpdgha called with %d %d \n", id, intfluka);
791 // MPDGHA() goes from fluka internal to pdg.
792 return mpdgha(intfluka);
794 // ions and optical photons
795 return idSpecial[id + 6];
799 void TFluka::StopTrack()
801 // Set stopping conditions
802 // Works for photons and charged particles
806 //_____________________________________________________________________________
807 // methods for physics management
808 //____________________________________________________________________________
813 void TFluka::SetProcess(const char* flagName, Int_t flagValue, Int_t imed)
815 // Set process user flag for material imat
817 TFlukaConfigOption* proc = new TFlukaConfigOption(flagName, flagValue, imed);
818 fProcesses->Add(proc);
821 //______________________________________________________________________________
822 Bool_t TFluka::SetProcess(const char* flagName, Int_t flagValue)
824 // Set process user flag
827 // Update if already in the list
830 TIter next(fProcesses);
831 TFlukaConfigOption* proc;
832 while((proc = (TFlukaConfigOption*)next()))
834 if (strcmp(proc->GetName(), flagName) == 0) {
835 proc->SetFlag(flagValue);
841 // If not create a new process
844 proc = new TFlukaConfigOption(flagName, flagValue);
845 fProcesses->Add(proc);
850 //______________________________________________________________________________
851 void TFluka::SetCut(const char* cutName, Double_t cutValue, Int_t imed)
853 // Set user cut value for material imed
855 TFlukaConfigOption* cut = new TFlukaConfigOption(cutName, cutValue, imed);
859 //______________________________________________________________________________
860 Bool_t TFluka::SetCut(const char* cutName, Double_t cutValue)
862 // Set user cut value
865 // Update if already in the list
869 TFlukaConfigOption* cut;
870 while((cut = (TFlukaConfigOption*)next()))
872 if (strcmp(cut->GetName(), cutName) == 0) {
873 cut->SetCut(cutValue);
878 // If not create a new process
881 cut = new TFlukaConfigOption(cutName, cutValue);
887 void TFluka::SetUserScoring(const char* option, Int_t npar, Float_t what[12])
890 // Ads a user scoring option to th list
892 TFlukaScoringOption* opt = new TFlukaScoringOption(option, "User Scoring", npar, what);
893 fUserScore->Add(opt);
897 //______________________________________________________________________________
898 Double_t TFluka::Xsec(char*, Double_t, Int_t, Int_t)
900 printf("WARNING: Xsec not yet implemented !\n"); return -1.;
904 //______________________________________________________________________________
905 void TFluka::InitPhysics()
908 // Physics initialisation with preparation of FLUKA input cards
910 printf("=>InitPhysics\n");
914 FILE *pFlukaVmcCoreInp, *pFlukaVmcFlukaMat, *pFlukaVmcInp;
919 Double_t three = 3.0;
921 Float_t fLastMaterial = fGeom->GetLastMaterialIndex();
922 if (fVerbosityLevel >= 3) printf(" last FLUKA material is %g\n", fLastMaterial);
925 TObjArray *matList = GetFlukaMaterials();
926 Int_t nmaterial = matList->GetEntriesFast();
927 fMaterials = new Int_t[nmaterial+3];
929 // construct file names
931 TString sFlukaVmcCoreInp = getenv("ALICE_ROOT");
932 sFlukaVmcCoreInp +="/TFluka/input/";
933 TString sFlukaVmcTmp = "flukaMat.inp";
934 TString sFlukaVmcInp = GetInputFileName();
935 sFlukaVmcCoreInp += GetCoreInputFileName();
939 if ((pFlukaVmcCoreInp = fopen(sFlukaVmcCoreInp.Data(),"r")) == NULL) {
940 printf("\nCannot open file %s\n",sFlukaVmcCoreInp.Data());
943 if ((pFlukaVmcFlukaMat = fopen(sFlukaVmcTmp.Data(),"r")) == NULL) {
944 printf("\nCannot open file %s\n",sFlukaVmcTmp.Data());
947 if ((pFlukaVmcInp = fopen(sFlukaVmcInp.Data(),"w")) == NULL) {
948 printf("\nCannot open file %s\n",sFlukaVmcInp.Data());
952 // copy core input file
954 Float_t fEventsPerRun;
956 while ((fgets(sLine,255,pFlukaVmcCoreInp)) != NULL) {
957 if (strncmp(sLine,"GEOEND",6) != 0)
958 fprintf(pFlukaVmcInp,"%s",sLine); // copy until GEOEND card
960 fprintf(pFlukaVmcInp,"GEOEND\n"); // add GEOEND card
963 } // end of while until GEOEND card
967 while ((fgets(sLine,255,pFlukaVmcFlukaMat)) != NULL) { // copy flukaMat.inp file
968 fprintf(pFlukaVmcInp,"%s\n",sLine);
971 while ((fgets(sLine,255,pFlukaVmcCoreInp)) != NULL) {
972 if (strncmp(sLine,"START",5) != 0)
973 fprintf(pFlukaVmcInp,"%s\n",sLine);
975 sscanf(sLine+10,"%10f",&fEventsPerRun);
978 } //end of while until START card
981 // in G3 the process control values meaning can be different for
982 // different processes, but for most of them is:
983 // 0 process is not activated
984 // 1 process is activated WITH generation of secondaries
985 // 2 process is activated WITHOUT generation of secondaries
986 // if process does not generate secondaries => 1 same as 2
995 // Loop over number of SetProcess calls
996 fprintf(pFlukaVmcInp,"*----------------------------------------------------------------------------- \n");
997 fprintf(pFlukaVmcInp,"*----- The following data are generated from SetProcess and SetCut calls ----- \n");
998 fprintf(pFlukaVmcInp,"*----------------------------------------------------------------------------- \n");
1000 // Outer loop over processes
1001 TIter next(fProcesses);
1002 TFlukaConfigOption *proc;
1003 // Inner loop over processes
1004 TIter nextp(fProcesses);
1005 TFlukaConfigOption *procp;
1008 TFlukaConfigOption *cut = 0x0;
1010 while((proc = (TFlukaConfigOption*)next())) {
1011 Float_t matMin = three;
1012 Float_t matMax = fLastMaterial;
1013 Bool_t global = kTRUE;
1014 if (proc->Medium() != -1) {
1016 if ((mat = proc->Medium()) >= GetFlukaMaterials()->GetEntries()) continue;
1017 matMin = Float_t(mat);
1021 fprintf(pFlukaVmcInp,"*\n*Material specific process setting for #%8d \n", mat);
1025 // G3 default value: 1
1026 // G4 processes: G4eplusAnnihilation/G4IeplusAnnihilation
1029 // flag = 0 no annihilation
1030 // flag = 1 annihilation, decays processed
1031 // flag = 2 annihilation, no decay product stored
1032 // gMC ->SetProcess("ANNI",1); // EMFCUT -1. 0. 0. 3. lastmat 0. ANNH-THR
1033 if (strncmp(proc->GetName(),"ANNI",4) == 0) {
1034 if (proc->Flag() == 1 || proc->Flag() == 2) {
1035 fprintf(pFlukaVmcInp,"*\n*Kinetic energy threshold (GeV) for e+ annihilation - resets to default=0.\n");
1036 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('ANNI',1) or SetProcess('ANNI',2)\n");
1037 // -one = kinetic energy threshold (GeV) for e+ annihilation (resets to default=0)
1040 // matMin = lower bound of the material indices in which the respective thresholds apply
1041 // matMax = upper bound of the material indices in which the respective thresholds apply
1042 // one = step length in assigning indices
1044 fprintf(pFlukaVmcInp,"EMFCUT %10.1f%10.1f%10.1f%10.1f%10.1f%10.1fANNH-THR\n",-one,zero,zero,matMin,matMax,one);
1046 else if (proc->Flag() == 0) {
1047 fprintf(pFlukaVmcInp,"*\n*No annihilation - no FLUKA card generated\n");
1048 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('ANNI',0)\n");
1051 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('ANNI',?) call.\n");
1052 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1056 // bremsstrahlung and pair production are both activated
1057 // G3 default value: 1
1058 // G4 processes: G4eBremsstrahlung/G4IeBremsstrahlung,
1059 // G4MuBremsstrahlung/G4IMuBremsstrahlung,
1060 // G4LowEnergyBremstrahlung
1061 // Particles: e-/e+; mu+/mu-
1063 // flag = 0 no bremsstrahlung
1064 // flag = 1 bremsstrahlung, photon processed
1065 // flag = 2 bremsstrahlung, no photon stored
1066 // gMC ->SetProcess("BREM",1); // PAIRBREM 2. 0. 0. 3. lastmat
1067 // EMFCUT -1. 0. 0. 3. lastmat 0. ELPO-THR
1068 // G3 default value: 1
1069 // G4 processes: G4GammaConversion,
1070 // G4MuPairProduction/G4IMuPairProduction
1071 // G4LowEnergyGammaConversion
1072 // Particles: gamma, mu
1074 // flag = 0 no delta rays
1075 // flag = 1 delta rays, secondaries processed
1076 // flag = 2 delta rays, no secondaries stored
1078 else if ((strncmp(proc->GetName(),"PAIR",4) == 0) && (proc->Flag() == 1 || proc->Flag() == 2)) {
1082 while ((procp = (TFlukaConfigOption*)nextp())) {
1083 if ((strncmp(procp->GetName(),"BREM",4) == 0) &&
1084 (proc->Flag() == 1) &&
1085 (procp->Medium() == proc->Medium())) {
1086 fprintf(pFlukaVmcInp,"*\n*Bremsstrahlung and pair production by muons and charged hadrons both activated\n");
1087 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('BREM',1) and SetProcess('PAIR',1)\n");
1088 fprintf(pFlukaVmcInp,"*Energy threshold set by call SetCut('BCUTM',cut) or set to 0.\n");
1089 fprintf(pFlukaVmcInp,"*Energy threshold set by call SetCut('PPCUTM',cut) or set to 0.\n");
1090 // three = bremsstrahlung and pair production by muons and charged hadrons both are activated
1091 fprintf(pFlukaVmcInp,"PAIRBREM %10.1f",three);
1092 // direct pair production by muons
1093 // G4 particles: "e-", "e+"
1094 // G3 default value: 0.01 GeV
1095 //gMC ->SetCut("PPCUTM",cut); // total energy cut for direct pair prod. by muons
1098 while ((cut = (TFlukaConfigOption*)nextc())) {
1099 if (strncmp(cut->GetName(), "PPCUTM", 6) == 0 &&
1100 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1102 fprintf(pFlukaVmcInp,"%10.4g",theCut);
1103 // theCut; = e+, e- kinetic energy threshold (in GeV) for explicit pair production.
1104 // muon and hadron bremsstrahlung
1105 // G4 particles: "gamma"
1106 // G3 default value: CUTGAM=0.001 GeV
1107 //gMC ->SetCut("BCUTM",cut); // cut for muon and hadron bremsstrahlung
1110 while ((cut = (TFlukaConfigOption*)nextc())) {
1111 if (strncmp(cut->GetName(), "BCUTM", 5) == 0 &&
1112 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1114 fprintf(pFlukaVmcInp,"%10.4g%10.1f%10.1f\n",theCut,matMin,matMax);
1115 // theCut = photon energy threshold (GeV) for explicit bremsstrahlung production
1116 // matMin = lower bound of the material indices in which the respective thresholds apply
1117 // matMax = upper bound of the material indices in which the respective thresholds apply
1120 fprintf(pFlukaVmcInp,"*\n*Kinetic energy threshold (GeV) for e+/e- bremsstrahlung - resets to default=0.\n");
1121 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('BREM',1);\n");
1124 while ((cut = (TFlukaConfigOption*)nextc())) {
1125 if (strncmp(cut->GetName(), "BCUTE", 5) == 0 &&
1126 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1128 //theCut = kinetic energy threshold (GeV) for e+/e- bremsstrahlung (resets to default=0)
1131 // matMin = lower bound of the material indices in which the respective thresholds apply
1132 // matMax = upper bound of the material indices in which the respective thresholds apply
1133 // one = step length in assigning indices
1135 fprintf(pFlukaVmcInp,"EMFCUT %10.4g%10.1f%10.1f%10.1f%10.1f%10.1fELPO-THR\n",theCut,zero,zero,matMin,matMax,one);
1137 // for gamma -> e+ and e-
1138 fprintf(pFlukaVmcInp,"*\n*Pair production by photons is activated\n");
1139 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('PAIR',1);\n");
1142 while ((cut = (TFlukaConfigOption*)nextc())) {
1143 if (strncmp(cut->GetName(), "CUTELE", 6) == 0 &&
1144 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1146 // theCut = energy threshold (GeV) for gamma pair production (< 0.0 : resets to default, = 0.0 : ignored)
1147 // matMin = lower bound of the material indices in which the respective thresholds apply
1148 // matMax = upper bound of the material indices in which the respective thresholds apply
1149 // one = step length in assigning indices
1150 fprintf(pFlukaVmcInp,"EMFCUT %10.1f%10.1f%10.4g%10.1f%10.1f%10.1fPHOT-THR\n",zero,zero,theCut,matMin,matMax,one);
1152 } // end of if for BREM
1153 } // end of loop for BREM
1155 // only pair production by muons and charged hadrons is activated
1156 fprintf(pFlukaVmcInp,"*\n*Pair production by muons and charged hadrons is activated\n");
1157 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('PAIR',1) or SetProcess('PAIR',2)\n");
1158 fprintf(pFlukaVmcInp,"*Energy threshold set by call SetCut('PPCUTM',cut) or set to 0.\n");
1159 // direct pair production by muons
1160 // G4 particles: "e-", "e+"
1161 // G3 default value: 0.01 GeV
1162 //gMC ->SetCut("PPCUTM",cut); // total energy cut for direct pair prod. by muons
1163 // one = pair production by muons and charged hadrons is activated
1164 // zero = e+, e- kinetic energy threshold (in GeV) for explicit pair production.
1165 // zero = no explicit bremsstrahlung production is simulated
1166 // matMin = lower bound of the material indices in which the respective thresholds apply
1167 // matMax = upper bound of the material indices in which the respective thresholds apply
1168 fprintf(pFlukaVmcInp,"PAIRBREM %10.1f%10.1f%10.1f%10.1f%10.1f\n",one,zero,zero,matMin,matMax);
1170 // for gamma -> e+ and e-
1171 fprintf(pFlukaVmcInp,"*\n*Pair production by electrons is activated\n");
1172 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('PAIR',1) or SetProcess('PAIR',2)\n");
1175 while ((cut = (TFlukaConfigOption*)nextc())) {
1176 if (strncmp(cut->GetName(), "CUTELE", 6) == 0 &&
1177 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1179 // zero = energy threshold (GeV) for Compton scattering (= 0.0 : ignored)
1180 // zero = energy threshold (GeV) for Photoelectric (= 0.0 : ignored)
1181 // theCut = energy threshold (GeV) for gamma pair production (< 0.0 : resets to default, = 0.0 : ignored)
1182 // matMin = lower bound of the material indices in which the respective thresholds apply
1183 // matMax = upper bound of the material indices in which the respective thresholds apply
1184 // one = step length in assigning indices
1185 fprintf(pFlukaVmcInp,"EMFCUT %10.1f%10.1f%10.4g%10.1f%10.1f%10.1fPHOT-THR\n",zero,zero,theCut,matMin,matMax,one);
1189 } // end of if for PAIR
1194 // G3 default value: 1
1195 // G4 processes: G4eBremsstrahlung/G4IeBremsstrahlung,
1196 // G4MuBremsstrahlung/G4IMuBremsstrahlung,
1197 // G4LowEnergyBremstrahlung
1198 // Particles: e-/e+; mu+/mu-
1200 // flag = 0 no bremsstrahlung
1201 // flag = 1 bremsstrahlung, photon processed
1202 // flag = 2 bremsstrahlung, no photon stored
1203 // gMC ->SetProcess("BREM",1); // PAIRBREM 2. 0. 0. 3. lastmat
1204 // EMFCUT -1. 0. 0. 3. lastmat 0. ELPO-THR
1205 else if (strncmp(proc->GetName(),"BREM",4) == 0) {
1207 while((procp = (TFlukaConfigOption*)nextp())) {
1208 if ((strncmp(procp->GetName(),"PAIR",4) == 0) &&
1209 procp->Flag() == 1 &&
1210 (procp->Medium() == proc->Medium())) goto NOBREM;
1212 if (proc->Flag() == 1) {
1213 fprintf(pFlukaVmcInp,"*\n*Bremsstrahlung by muons and charged hadrons is activated\n");
1214 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('BREM',1) or SetProcess('BREM',2)\n");
1215 fprintf(pFlukaVmcInp,"*Energy threshold set by call SetCut('BCUTM',cut) or set to 0.\n");
1216 // two = bremsstrahlung by muons and charged hadrons is activated
1217 // zero = no meaning
1218 // muon and hadron bremsstrahlung
1219 // G4 particles: "gamma"
1220 // G3 default value: CUTGAM=0.001 GeV
1221 //gMC ->SetCut("BCUTM",cut); // cut for muon and hadron bremsstrahlung
1224 while ((cut = (TFlukaConfigOption*)nextc())) {
1225 if (strncmp(cut->GetName(), "BCUTM", 5) == 0 &&
1226 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1228 // theCut = photon energy threshold (GeV) for explicit bremsstrahlung production
1229 // matMin = lower bound of the material indices in which the respective thresholds apply
1230 // matMax = upper bound of the material indices in which the respective thresholds apply
1231 fprintf(pFlukaVmcInp,"PAIRBREM %10.1f%10.1f%10.4g%10.1f%10.1f\n",two,zero,theCut,matMin,matMax);
1234 fprintf(pFlukaVmcInp,"*\n*Kinetic energy threshold (GeV) for e+/e- bremsstrahlung - resets to default=0.\n");
1235 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('BREM',1);");
1236 // - one = kinetic energy threshold (GeV) for e+/e- bremsstrahlung (resets to default=0)
1239 // matMin = lower bound of the material indices in which the respective thresholds apply
1240 // matMax = upper bound of the material indices in which the respective thresholds apply
1241 // one = step length in assigning indices
1245 while ((cut = (TFlukaConfigOption*)nextc())) {
1246 if (strncmp(cut->GetName(), "CUTGAM", 6) == 0 &&
1247 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1249 fprintf(pFlukaVmcInp,"EMFCUT %10.4g%10.1f%10.1f%10.1f%10.1f%10.1fELPO-THR\n", theCut,zero,zero,matMin,matMax,one);
1251 else if (proc->Flag() == 0) {
1252 fprintf(pFlukaVmcInp,"*\n*No bremsstrahlung - no FLUKA card generated\n");
1253 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('BREM',0)\n");
1256 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('BREM',?) call.\n");
1257 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1261 } // end of else if (strncmp(proc->GetName(),"BREM",4) == 0)
1263 // Cerenkov photon generation
1264 // G3 default value: 0
1265 // G4 process: G4Cerenkov
1267 // Particles: charged
1269 // flag = 0 no Cerenkov photon generation
1270 // flag = 1 Cerenkov photon generation
1271 // flag = 2 Cerenkov photon generation with primary stopped at each step
1272 //xx gMC ->SetProcess("CKOV",1); // ??? Cerenkov photon generation
1274 else if (strncmp(proc->GetName(),"CKOV",4) == 0) {
1275 if ((proc->Flag() == 1 || proc->Flag() == 2) && global) {
1277 fprintf(pFlukaVmcInp, "* \n");
1278 fprintf(pFlukaVmcInp, "*Cerenkov photon generation\n");
1279 fprintf(pFlukaVmcInp, "*Generated from call: SetProcess('CKOV',1) or SetProcess('CKOV',2)\n");
1281 for (Int_t im = 0; im < nmaterial; im++)
1283 TGeoMaterial* material = dynamic_cast<TGeoMaterial*> (matList->At(im));
1284 Int_t idmat = material->GetIndex();
1286 if (!global && idmat != proc->Medium()) continue;
1288 fMaterials[idmat] = im;
1289 // Skip media with no Cerenkov properties
1290 TFlukaCerenkov* cerenkovProp;
1291 if (!(cerenkovProp = dynamic_cast<TFlukaCerenkov*>(material->GetCerenkovProperties()))) continue;
1293 // This medium has Cerenkov properties
1296 // Write OPT-PROD card for each medium
1297 Float_t emin = cerenkovProp->GetMinimumEnergy();
1298 Float_t emax = cerenkovProp->GetMaximumEnergy();
1299 fprintf(pFlukaVmcInp, "OPT-PROD %10.4g%10.4g%10.4g%10.4g%10.4g%10.4gCERENKOV\n", emin, emax, 0.,
1300 Float_t(idmat), Float_t(idmat), 0.);
1302 // Write OPT-PROP card for each medium
1303 // Forcing FLUKA to call user routines (queffc.cxx, rflctv.cxx, rfrndx.cxx)
1305 fprintf(pFlukaVmcInp, "OPT-PROP %10.4g%10.4g%10.4g%10.1f%10.1f%10.1fWV-LIMIT\n",
1306 cerenkovProp->GetMinimumWavelength(),
1307 cerenkovProp->GetMaximumWavelength(),
1308 cerenkovProp->GetMaximumWavelength(),
1309 Float_t(idmat), Float_t(idmat), 0.0);
1311 if (cerenkovProp->IsMetal()) {
1312 fprintf(pFlukaVmcInp, "OPT-PROP %10.1f%10.1f%10.1f%10.1f%10.1f%10.1fMETAL\n",
1313 -100., -100., -100.,
1314 Float_t(idmat), Float_t(idmat), 0.0);
1316 fprintf(pFlukaVmcInp, "OPT-PROP %10.1f%10.1f%10.1f%10.1f%10.1f%10.1f\n",
1317 -100., -100., -100.,
1318 Float_t(idmat), Float_t(idmat), 0.0);
1322 for (Int_t j = 0; j < 3; j++) {
1323 fprintf(pFlukaVmcInp, "OPT-PROP %10.1f%10.1f%10.1f%10.1f%10.1f%10.1f&\n",
1324 -100., -100., -100.,
1325 Float_t(idmat), Float_t(idmat), 0.0);
1327 // Photon detection efficiency user defined
1329 if (cerenkovProp->IsSensitive())
1330 fprintf(pFlukaVmcInp, "OPT-PROP %10.1f%10.1f%10.1f%10.1f%10.1f%10.1fSENSITIV\n",
1331 -100., -100., -100.,
1332 Float_t(idmat), Float_t(idmat), 0.0);
1335 } else if (proc->Flag() == 0) {
1336 fprintf(pFlukaVmcInp,"*\n*No Cerenkov photon generation\n");
1337 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('CKOV',0)\n");
1341 // matMin = lower bound of the material indices in which the respective thresholds apply
1342 // matMax = upper bound of the material indices in which the respective thresholds apply
1343 // one = step length in assigning indices
1345 fprintf(pFlukaVmcInp,"OPT-PROD %10.1f%10.1f%10.1f%10.1f%10.1f%10.1fCERE-OFF\n",zero,zero,zero,matMin,matMax,one);
1348 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('CKOV',?) call.\n");
1349 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1351 } // end of else if (strncmp(proc->GetName(),"CKOV",4) == 0)
1353 // Compton scattering
1354 // G3 default value: 1
1355 // G4 processes: G4ComptonScattering,
1356 // G4LowEnergyCompton,
1357 // G4PolarizedComptonScattering
1360 // flag = 0 no Compton scattering
1361 // flag = 1 Compton scattering, electron processed
1362 // flag = 2 Compton scattering, no electron stored
1363 // gMC ->SetProcess("COMP",1); // EMFCUT -1. 0. 0. 3. lastmat 0. PHOT-THR
1364 else if (strncmp(proc->GetName(),"COMP",4) == 0) {
1365 if (proc->Flag() == 1 || proc->Flag() == 2) {
1366 fprintf(pFlukaVmcInp,"*\n*Energy threshold (GeV) for Compton scattering - resets to default=0.\n");
1367 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('COMP',1);\n");
1368 // - one = energy threshold (GeV) for Compton scattering - resets to default=0.
1371 // matMin = lower bound of the material indices in which the respective thresholds apply
1372 // matMax = upper bound of the material indices in which the respective thresholds apply
1373 // one = step length in assigning indices
1377 while ((cut = (TFlukaConfigOption*)nextc())) {
1378 if (strncmp(cut->GetName(), "CUTELE", 6) == 0 &&
1379 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1381 fprintf(pFlukaVmcInp,"EMFCUT %10.4g%10.1f%10.1f%10.1f%10.1f%10.1fPHOT-THR\n",theCut,zero,zero,matMin,matMax,one);
1383 else if (proc->Flag() == 0) {
1384 fprintf(pFlukaVmcInp,"*\n*No Compton scattering - no FLUKA card generated\n");
1385 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('COMP',0)\n");
1388 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('COMP',?) call.\n");
1389 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1391 } // end of else if (strncmp(proc->GetName(),"COMP",4) == 0)
1394 // G3 default value: 1
1395 // G4 process: G4Decay
1397 // Particles: all which decay is applicable for
1399 // flag = 0 no decays
1400 // flag = 1 decays, secondaries processed
1401 // flag = 2 decays, no secondaries stored
1402 //gMC ->SetProcess("DCAY",0); // not available
1403 else if ((strncmp(proc->GetName(),"DCAY",4) == 0) && proc->Flag() == 0)
1404 cout << "SetProcess for flag =" << proc->GetName() << " value=" << proc->Flag() << " not avaliable!" << endl;
1405 else if ((strncmp(proc->GetName(),"DCAY",4) == 0) && proc->Flag() == 1) {
1406 // Nothing to do decays are switched on by default
1411 // G3 default value: 2
1412 // !! G4 treats delta rays in different way
1413 // G4 processes: G4eIonisation/G4IeIonization,
1414 // G4MuIonisation/G4IMuIonization,
1415 // G4hIonisation/G4IhIonisation
1416 // Particles: charged
1418 // flag = 0 no energy loss
1419 // flag = 1 restricted energy loss fluctuations
1420 // flag = 2 complete energy loss fluctuations
1421 // flag = 3 same as 1
1422 // flag = 4 no energy loss fluctuations
1423 // gMC ->SetProcess("DRAY",0); // DELTARAY 1.E+6 0. 0. 3. lastmat 0.
1424 else if (strncmp(proc->GetName(),"DRAY",4) == 0) {
1425 if (proc->Flag() == 0 || proc->Flag() == 4) {
1426 fprintf(pFlukaVmcInp,"*\n*Kinetic energy threshold (GeV) for delta ray production\n");
1427 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('DRAY',0) or SetProcess('DRAY',4)\n");
1428 fprintf(pFlukaVmcInp,"*No delta ray production by muons - threshold set artificially high\n");
1429 Double_t emin = 1.0e+6; // kinetic energy threshold (GeV) for delta ray production (discrete energy transfer)
1432 // matMin = lower bound of the material indices in which the respective thresholds apply
1433 // matMax = upper bound of the material indices in which the respective thresholds apply
1434 // one = step length in assigning indices
1435 fprintf(pFlukaVmcInp,"DELTARAY %10.4g%10.1f%10.1f%10.1f%10.1f%10.1f\n",emin,zero,zero,matMin,matMax,one);
1437 else if (proc->Flag() == 1 || proc->Flag() == 2 || proc->Flag() == 3) {
1438 fprintf(pFlukaVmcInp,"*\n*Kinetic energy threshold (GeV) for delta ray production\n");
1439 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('DRAY',flag), flag=1,2,3\n");
1440 fprintf(pFlukaVmcInp,"*Delta ray production by muons switched on\n");
1441 fprintf(pFlukaVmcInp,"*Energy threshold set by call SetCut('DCUTM',cut) or set to 1.0e+6.\n");
1445 // Check cut one delta-rays from electrons
1447 while ((cut = (TFlukaConfigOption*)nextc())) {
1448 if (strncmp(cut->GetName(), "DCUTM", 5) == 0 &&
1449 cut->Medium() == proc->Medium()) theCut = cut->Cut();
1451 // theCut = kinetic energy threshold (GeV) for delta ray production (discrete energy transfer)
1454 // matMin = lower bound of the material indices in which the respective thresholds apply
1455 // matMax = upper bound of the material indices in which the respective thresholds apply
1456 // one = step length in assigning indices
1457 fprintf(pFlukaVmcInp,"DELTARAY %10.4g%10.1f%10.1f%10.1f%10.1f%10.1f\n",theCut,zero,zero,matMin,matMax,one);
1460 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('DRAY',?) call.\n");
1461 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1463 } // end of else if (strncmp(proc->GetName(),"DRAY",4) == 0)
1466 // G3 default value: 1
1467 // G4 processes: all defined by TG4PhysicsConstructorHadron
1469 // Particles: hadrons
1471 // flag = 0 no multiple scattering
1472 // flag = 1 hadronic interactions, secondaries processed
1473 // flag = 2 hadronic interactions, no secondaries stored
1474 // gMC ->SetProcess("HADR",1); // ??? hadronic process
1475 //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3) ?????
1476 else if (strncmp(proc->GetName(),"HADR",4) == 0) {
1477 if (proc->Flag() == 1 || proc->Flag() == 2) {
1478 fprintf(pFlukaVmcInp,"*\n*Hadronic interaction is ON by default in FLUKA\n");
1479 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1481 else if (proc->Flag() == 0) {
1482 fprintf(pFlukaVmcInp,"*\n*Hadronic interaction is set OFF\n");
1483 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('HADR',0);\n");
1484 fprintf(pFlukaVmcInp,"*Switching off hadronic interactions not foreseen in FLUKA\n");
1485 fprintf(pFlukaVmcInp,"THRESHOL %10.1f%10.1f%10.1f%10.1e%10.1f\n",zero, zero, zero, 1.e10, zero);
1488 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('HADR',?) call.\n");
1489 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1491 } // end of else if (strncmp(proc->GetName(),"HADR",4) == 0)
1495 // G3 default value: 2
1496 // G4 processes: G4eIonisation/G4IeIonization,
1497 // G4MuIonisation/G4IMuIonization,
1498 // G4hIonisation/G4IhIonisation
1500 // Particles: charged
1502 // flag=0 no energy loss
1503 // flag=1 restricted energy loss fluctuations
1504 // flag=2 complete energy loss fluctuations
1506 // flag=4 no energy loss fluctuations
1507 // If the value ILOSS is changed, then (in G3) cross-sections and energy
1508 // loss tables must be recomputed via the command 'PHYSI'
1509 // gMC ->SetProcess("LOSS",2); // ??? IONFLUCT ? energy loss
1510 else if (strncmp(proc->GetName(),"LOSS",4) == 0) {
1511 if (proc->Flag() > 0 || proc->Flag() < 4) { // restricted energy loss fluctuations
1512 fprintf(pFlukaVmcInp,"*\n*Restricted energy loss fluctuations\n");
1513 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('LOSS',1) or SetProcess('LOSS',3)\n");
1514 // one = restricted energy loss fluctuations (for hadrons and muons) switched on
1515 // one = restricted energy loss fluctuations (for e+ and e-) switched on
1516 // one = minimal accuracy
1517 // matMin = lower bound of the material indices in which the respective thresholds apply
1518 // upper bound of the material indices in which the respective thresholds apply
1519 fprintf(pFlukaVmcInp,"IONFLUCT %10.1f%10.1f%10.1f%10.1f%10.1f\n",one, one, 4., matMin, matMax);
1521 else if (proc->Flag() == 4) { // no energy loss fluctuations
1522 fprintf(pFlukaVmcInp,"*\n*No energy loss fluctuations\n");
1523 fprintf(pFlukaVmcInp,"*\n*Generated from call: SetProcess('LOSS',4)\n");
1524 // - one = restricted energy loss fluctuations (for hadrons and muons) switched off
1525 // - one = restricted energy loss fluctuations (for e+ and e-) switched off
1526 // one = minimal accuracy
1527 // matMin = lower bound of the material indices in which the respective thresholds apply
1528 // matMax = upper bound of the material indices in which the respective thresholds apply
1529 fprintf(pFlukaVmcInp,"IONFLUCT %10.1f%10.1f%10.1f%10.1f%10.1f\n",-one,-one,one,matMin,matMax);
1532 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('LOSS',?) call.\n");
1533 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1535 } // end of else if (strncmp(proc->GetName(),"LOSS",4) == 0)
1538 // multiple scattering
1539 // G3 default value: 1
1540 // G4 process: G4MultipleScattering/G4IMultipleScattering
1542 // Particles: charged
1544 // flag = 0 no multiple scattering
1545 // flag = 1 Moliere or Coulomb scattering
1546 // flag = 2 Moliere or Coulomb scattering
1547 // flag = 3 Gaussian scattering
1548 // gMC ->SetProcess("MULS",1); // MULSOPT multiple scattering
1549 else if (strncmp(proc->GetName(),"MULS",4) == 0) {
1550 if (proc->Flag() == 1 || proc->Flag() == 2 || proc->Flag() == 3) {
1551 fprintf(pFlukaVmcInp,"*\n*Multiple scattering is ON by default for e+e- and for hadrons/muons\n");
1552 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1554 else if (proc->Flag() == 0) {
1555 fprintf(pFlukaVmcInp,"*\n*Multiple scattering is set OFF\n");
1556 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('MULS',0);\n");
1558 // three = multiple scattering for hadrons and muons is completely suppressed
1559 // three = multiple scattering for e+ and e- is completely suppressed
1560 // matMin = lower bound of the material indices in which the respective thresholds apply
1561 // matMax = upper bound of the material indices in which the respective thresholds apply
1562 fprintf(pFlukaVmcInp,"MULSOPT %10.1f%10.1f%10.1f%10.1f%10.1f\n",zero,three,three,matMin,matMax);
1565 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('MULS',?) call.\n");
1566 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1568 } // end of else if (strncmp(proc->GetName(),"MULS",4) == 0)
1571 // muon nuclear interaction
1572 // G3 default value: 0
1573 // G4 processes: G4MuNuclearInteraction,
1574 // G4MuonMinusCaptureAtRest
1578 // flag = 0 no muon-nuclear interaction
1579 // flag = 1 nuclear interaction, secondaries processed
1580 // flag = 2 nuclear interaction, secondaries not processed
1581 // gMC ->SetProcess("MUNU",1); // MUPHOTON 1. 0. 0. 3. lastmat
1582 else if (strncmp(proc->GetName(),"MUNU",4) == 0) {
1583 if (proc->Flag() == 1) {
1584 fprintf(pFlukaVmcInp,"*\n*Muon nuclear interactions with production of secondary hadrons\n");
1585 fprintf(pFlukaVmcInp,"*\n*Generated from call: SetProcess('MUNU',1);\n");
1586 // one = full simulation of muon nuclear interactions and production of secondary hadrons
1587 // zero = ratio of longitudinal to transverse virtual photon cross-section - Default = 0.25.
1588 // zero = fraction of rho-like interactions ( must be < 1) - Default = 0.75.
1589 // matMin = lower bound of the material indices in which the respective thresholds apply
1590 // matMax = upper bound of the material indices in which the respective thresholds apply
1591 fprintf(pFlukaVmcInp,"MUPHOTON %10.1f%10.1f%10.1f%10.1f%10.1f\n",one,zero,zero,matMin,matMax);
1593 else if (proc->Flag() == 2) {
1594 fprintf(pFlukaVmcInp,"*\n*Muon nuclear interactions without production of secondary hadrons\n");
1595 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('MUNU',2);\n");
1596 // two = full simulation of muon nuclear interactions and production of secondary hadrons
1597 // zero = ratio of longitudinal to transverse virtual photon cross-section - Default = 0.25.
1598 // zero = fraction of rho-like interactions ( must be < 1) - Default = 0.75.
1599 // matMin = lower bound of the material indices in which the respective thresholds apply
1600 // matMax = upper bound of the material indices in which the respective thresholds apply
1601 fprintf(pFlukaVmcInp,"MUPHOTON %10.1f%10.1f%10.1f%10.1f%10.1f\n",two,zero,zero,matMin,matMax);
1603 else if (proc->Flag() == 0) {
1604 fprintf(pFlukaVmcInp,"*\n*No muon nuclear interaction - no FLUKA card generated\n");
1605 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('MUNU',0)\n");
1608 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('MUNU',?) call.\n");
1609 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1611 } // end of else if (strncmp(proc->GetName(),"MUNU",4) == 0)
1615 // G3 default value: 0
1620 // gMC ->SetProcess("PFIS",0); // PHOTONUC -1. 0. 0. 3. lastmat 0.
1621 // flag = 0 no photon fission
1622 // flag = 1 photon fission, secondaries processed
1623 // flag = 2 photon fission, no secondaries stored
1624 else if (strncmp(proc->GetName(),"PFIS",4) == 0) {
1625 if (proc->Flag() == 0) {
1626 fprintf(pFlukaVmcInp,"*\n*No photonuclear interactions\n");
1627 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('PFIS',0);\n");
1628 // - one = no photonuclear interactions
1631 // matMin = lower bound of the material indices in which the respective thresholds apply
1632 // matMax = upper bound of the material indices in which the respective thresholds apply
1633 fprintf(pFlukaVmcInp,"PHOTONUC %10.1f%10.1f%10.1f%10.1f%10.1f\n",-one,zero,zero,matMin,matMax);
1635 else if (proc->Flag() == 1) {
1636 fprintf(pFlukaVmcInp,"*\n*Photon nuclear interactions are activated at all energies\n");
1637 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('PFIS',1);\n");
1638 // one = photonuclear interactions are activated at all energies
1641 // matMin = lower bound of the material indices in which the respective thresholds apply
1642 // matMax = upper bound of the material indices in which the respective thresholds apply
1643 fprintf(pFlukaVmcInp,"PHOTONUC %10.1f%10.1f%10.1f%10.1f%10.1f\n",one,zero,zero,matMin,matMax);
1645 else if (proc->Flag() == 0) {
1646 fprintf(pFlukaVmcInp,"*\n*No photofission - no FLUKA card generated\n");
1647 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('PFIS',0)\n");
1650 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('PFIS',?) call.\n");
1651 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1656 // photo electric effect
1657 // G3 default value: 1
1658 // G4 processes: G4PhotoElectricEffect
1659 // G4LowEnergyPhotoElectric
1662 // flag = 0 no photo electric effect
1663 // flag = 1 photo electric effect, electron processed
1664 // flag = 2 photo electric effect, no electron stored
1665 // gMC ->SetProcess("PHOT",1); // EMFCUT 0. -1. 0. 3. lastmat 0. PHOT-THR
1666 else if (strncmp(proc->GetName(),"PHOT",4) == 0) {
1667 if (proc->Flag() == 1) {
1668 fprintf(pFlukaVmcInp,"*\n*Photo electric effect is activated\n");
1669 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('PHOT',1);\n");
1671 // - one = resets to default=0.
1673 // matMin = lower bound of the material indices in which the respective thresholds apply
1674 // matMax = upper bound of the material indices in which the respective thresholds apply
1675 // one = step length in assigning indices
1679 while ((cut = (TFlukaConfigOption*)nextc())) {
1680 if (strncmp(cut->GetName(), "CUTELE", 6) == 0 &&
1681 (cut->Medium() == proc->Medium())) theCut = cut->Cut();
1683 fprintf(pFlukaVmcInp,"EMFCUT %10.1f%10.4g%10.1f%10.1f%10.1f%10.1fPHOT-THR\n",zero,theCut,zero,matMin,matMax,one);
1685 else if (proc->Flag() == 0) {
1686 fprintf(pFlukaVmcInp,"*\n*No photo electric effect - no FLUKA card generated\n");
1687 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('PHOT',0)\n");
1690 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('PHOT',?) call.\n");
1691 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1693 } // else if (strncmp(proc->GetName(),"PHOT",4) == 0)
1696 // Rayleigh scattering
1697 // G3 default value: 0
1698 // G4 process: G4OpRayleigh
1700 // Particles: optical photon
1702 // flag = 0 Rayleigh scattering off
1703 // flag = 1 Rayleigh scattering on
1704 //xx gMC ->SetProcess("RAYL",1);
1705 else if (strncmp(proc->GetName(),"RAYL",4) == 0) {
1706 if (proc->Flag() == 1) {
1707 fprintf(pFlukaVmcInp,"*\n*Rayleigh scattering is ON by default in FLUKA\n");
1708 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1710 else if (proc->Flag() == 0) {
1711 fprintf(pFlukaVmcInp,"*\n*Rayleigh scattering is set OFF\n");
1712 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('RAYL',0);\n");
1713 // - one = no Rayleigh scattering and no binding corrections for Compton
1714 // matMin = lower bound of the material indices in which the respective thresholds apply
1715 // matMax = upper bound of the material indices in which the respective thresholds apply
1716 fprintf(pFlukaVmcInp,"EMFRAY %10.1f%10.1f%10.1f%10.1f\n",-one,three,matMin,matMax);
1719 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('RAYL',?) call.\n");
1720 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1722 } // end of else if (strncmp(proc->GetName(),"RAYL",4) == 0)
1725 // synchrotron radiation in magnetic field
1726 // G3 default value: 0
1727 // G4 process: G4SynchrotronRadiation
1731 // flag = 0 no synchrotron radiation
1732 // flag = 1 synchrotron radiation
1733 //xx gMC ->SetProcess("SYNC",1); // synchrotron radiation generation
1734 else if (strncmp(proc->GetName(),"SYNC",4) == 0) {
1735 fprintf(pFlukaVmcInp,"*\n*Synchrotron radiation generation is NOT implemented in FLUKA\n");
1736 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1740 // Automatic calculation of tracking medium parameters
1741 // flag = 0 no automatic calculation
1742 // flag = 1 automatic calculation
1743 //xx gMC ->SetProcess("AUTO",1); // ??? automatic computation of the tracking medium parameters
1744 else if (strncmp(proc->GetName(),"AUTO",4) == 0) {
1745 fprintf(pFlukaVmcInp,"*\n*Automatic calculation of tracking medium parameters is always ON in FLUKA\n");
1746 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1750 // To control energy loss fluctuation model
1751 // flag = 0 Urban model
1752 // flag = 1 PAI model
1753 // flag = 2 PAI+ASHO model (not active at the moment)
1754 //xx gMC ->SetProcess("STRA",1); // ??? energy fluctuation model
1755 else if (strncmp(proc->GetName(),"STRA",4) == 0) {
1756 if (proc->Flag() == 0 || proc->Flag() == 2 || proc->Flag() == 3) {
1757 fprintf(pFlukaVmcInp,"*\n*Ionization energy losses calculation is activated\n");
1758 fprintf(pFlukaVmcInp,"*Generated from call: SetProcess('STRA',n);, n=0,1,2\n");
1759 // one = restricted energy loss fluctuations (for hadrons and muons) switched on
1760 // one = restricted energy loss fluctuations (for e+ and e-) switched on
1761 // one = minimal accuracy
1762 // matMin = lower bound of the material indices in which the respective thresholds apply
1763 // matMax = upper bound of the material indices in which the respective thresholds apply
1764 fprintf(pFlukaVmcInp,"IONFLUCT %10.1f%10.1f%10.1f%10.1f%10.1f\n",one,one,one,matMin,matMax);
1767 fprintf(pFlukaVmcInp,"*\n*Illegal flag value in SetProcess('STRA',?) call.\n");
1768 fprintf(pFlukaVmcInp,"*No FLUKA card generated\n");
1770 } // else if (strncmp(proc->GetName(),"STRA",4) == 0)
1775 else { // processes not yet treated
1777 // light photon absorption (Cerenkov photons)
1778 // it is turned on when Cerenkov process is turned on
1779 // G3 default value: 0
1780 // G4 process: G4OpAbsorption, G4OpBoundaryProcess
1782 // Particles: optical photon
1784 // flag = 0 no absorption of Cerenkov photons
1785 // flag = 1 absorption of Cerenkov photons
1786 // gMC ->SetProcess("LABS",2); // ??? Cerenkov light absorption
1790 cout << "SetProcess for flag=" << proc->GetName() << " value=" << proc->Flag() << " not yet implemented!" << endl;
1792 } //end of loop number of SetProcess calls
1795 // Loop over number of SetCut calls
1798 while ((cut = (TFlukaConfigOption*)nextc())) {
1799 Float_t matMin = three;
1800 Float_t matMax = fLastMaterial;
1801 Bool_t global = kTRUE;
1802 if (cut->Medium() != -1) {
1804 if ((mat = cut->Medium()) >= GetFlukaMaterials()->GetEntries()) continue;
1805 matMin = Float_t(mat);
1808 TGeoMaterial* material = (TGeoMaterial*) (GetFlukaMaterials())->At(GetMaterialIndex(mat));
1809 fprintf(pFlukaVmcInp,"*\n*Material specific cut setting for #%8d %s %s %13.3e\n",
1810 mat, material->GetName(), cut->GetName(), cut->Cut());
1814 // cuts handled in SetProcess calls
1815 if (strncmp(cut->GetName(),"BCUTM",5) == 0) continue;
1816 else if (strncmp(cut->GetName(),"BCUTE",5) == 0) continue;
1817 else if (strncmp(cut->GetName(),"DCUTM",5) == 0) continue;
1818 else if (strncmp(cut->GetName(),"PPCUTM",6) == 0) continue;
1821 // G4 particles: "gamma"
1822 // G3 default value: 0.001 GeV
1823 // gMC ->SetCut("CUTGAM",cut); // cut for gammas
1825 else if (strncmp(cut->GetName(),"CUTGAM",6) == 0 && global) {
1826 fprintf(pFlukaVmcInp,"*\n*Cut for gamma\n");
1827 fprintf(pFlukaVmcInp,"*Generated from call: SetCut('CUTGAM',cut);\n");
1828 fprintf(pFlukaVmcInp,"EMFCUT %10.4g%10.4g%10.1f%10.1f%10.1f%10.1f\n",
1829 zero, cut->Cut(), zero, zero, Float_t(fGeom->NofVolumes()), one);
1831 else if (strncmp(cut->GetName(),"CUTGAM",6) == 0 && !global) {
1832 // loop over materials for EMFCUT FLUKA cards
1833 for (j=0; j < matMax-matMin+1; j++) {
1834 Int_t nreg, imat, *reglist;
1836 imat = (Int_t) matMin + j;
1837 reglist = fGeom->GetMaterialList(imat, nreg);
1838 // loop over regions of a given material
1839 for (Int_t k = 0; k < nreg; k++) {
1841 fprintf(pFlukaVmcInp,"EMFCUT %10.4g%10.4g%10.1f%10.1f%10.1f%10.1f\n", zero, cut->Cut(), zero, ireg, ireg, one);
1844 } // end of else if for gamma
1848 // G4 particles: "e-"
1850 // G3 default value: 0.001 GeV
1851 //gMC ->SetCut("CUTELE",cut); // cut for e+,e-
1852 else if (strncmp(cut->GetName(),"CUTELE",6) == 0 && global) {
1853 fprintf(pFlukaVmcInp,"*\n*Cut for electrons\n");
1854 fprintf(pFlukaVmcInp,"*Generated from call: SetCut('CUTELE',cut);\n");
1855 fprintf(pFlukaVmcInp,"EMFCUT %10.4g%10.4g%10.1f%10.1f%10.1f%10.1f\n",
1856 -cut->Cut(), zero, zero, zero, Float_t(fGeom->NofVolumes()), one);
1858 else if (strncmp(cut->GetName(),"CUTELE",6) == 0 && !global) {
1859 // loop over materials for EMFCUT FLUKA cards
1860 for (j=0; j < matMax-matMin+1; j++) {
1861 Int_t nreg, imat, *reglist;
1863 imat = (Int_t) matMin + j;
1864 reglist = fGeom->GetMaterialList(imat, nreg);
1865 // loop over regions of a given material
1866 for (k=0; k<nreg; k++) {
1868 fprintf(pFlukaVmcInp,"EMFCUT %10.4g%10.4g%10.1f%10.1f%10.1f%10.1f\n", -cut->Cut(), zero, zero, ireg, ireg, one);
1871 } // end of else if for electrons
1875 // G4 particles: of type "baryon", "meson", "nucleus" with zero charge
1876 // G3 default value: 0.01 GeV
1877 //gMC ->SetCut("CUTNEU",cut); // cut for neutral hadrons
1878 else if (strncmp(cut->GetName(),"CUTNEU",6) == 0 && global) {
1879 fprintf(pFlukaVmcInp,"*\n*Cut for neutral hadrons\n");
1880 fprintf(pFlukaVmcInp,"*Generated from call: SetCut('CUTNEU',cut);\n");
1883 // 9.0 = Antineutron
1884 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),8.0,9.0);
1886 // 12.0 = Kaon zero long
1887 // 12.0 = Kaon zero long
1888 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),12.0,12.0);
1890 // 17.0 = Lambda, 18.0 = Antilambda
1891 // 19.0 = Kaon zero short
1892 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),17.0,19.0);
1894 // 22.0 = Sigma zero, Pion zero, Kaon zero
1895 // 25.0 = Antikaon zero
1896 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),22.0,25.0);
1898 // 32.0 = Antisigma zero
1899 // 32.0 = Antisigma zero
1900 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),32.0,32.0);
1903 // 35.0 = AntiXi zero
1904 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),34.0,35.0);
1907 // 48.0 = AntiD zero
1908 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),47.0,48.0);
1912 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),53.0,53.0);
1914 // 55.0 = Xi'_c zero
1915 // 56.0 = Omega_c zero
1916 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),55.0,56.0);
1918 // 59.0 = AntiXi_c zero
1919 // 59.0 = AntiXi_c zero
1920 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),59.0,59.0);
1922 // 61.0 = AntiXi'_c zero
1923 // 62.0 = AntiOmega_c zero
1924 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),61.0,62.0);
1928 // G4 particles: of type "baryon", "meson", "nucleus" with non-zero charge
1929 // G3 default value: 0.01 GeV
1930 //gMC ->SetCut("CUTHAD",cut); // cut for charged hadrons
1931 else if (strncmp(cut->GetName(),"CUTHAD",6) == 0 && global) {
1932 fprintf(pFlukaVmcInp,"*\n*Cut for charged hadrons\n");
1933 fprintf(pFlukaVmcInp,"*Generated from call: SetCut('CUTHAD',cut);\n");
1937 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),1.0,2.0);
1939 // 13.0 = Positive Pion, Negative Pion, Positive Kaon
1940 // 16.0 = Negative Kaon
1941 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),13.0,16.0);
1943 // 20.0 = Negative Sigma
1944 // 21.0 = Positive Sigma
1945 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),20.0,21.0);
1947 // 31.0 = Antisigma minus
1948 // 33.0 = Antisigma plus
1949 // 2.0 = step length
1950 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f%10.1f\n",-cut->Cut(),31.0,33.0,2.0);
1952 // 36.0 = Negative Xi, Positive Xi, Omega minus
1954 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),36.0,39.0);
1958 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),45.0,46.0);
1960 // 49.0 = D_s plus, D_s minus, Lambda_c plus
1962 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),49.0,52.0);
1964 // 54.0 = Xi'_c plus
1965 // 60.0 = AntiXi'_c minus
1966 // 6.0 = step length
1967 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f%10.1f\n",-cut->Cut(),54.0,60.0,6.0);
1969 // 57.0 = Antilambda_c minus
1970 // 58.0 = AntiXi_c minus
1971 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),57.0,58.0);
1975 // G4 particles: "mu+", "mu-"
1976 // G3 default value: 0.01 GeV
1977 //gMC ->SetCut("CUTMUO",cut); // cut for mu+, mu-
1978 else if (strncmp(cut->GetName(),"CUTMUO",6)== 0 && global) {
1979 fprintf(pFlukaVmcInp,"*\n*Cut for muons\n");
1980 fprintf(pFlukaVmcInp,"*Generated from call: SetCut('CUTMUO',cut);\n");
1983 fprintf(pFlukaVmcInp,"PART-THR %10.4g%10.1f%10.1f\n",-cut->Cut(),10.0,11.0);
1987 // time of flight cut in seconds
1988 // G4 particles: all
1989 // G3 default value: 0.01 GeV
1990 //gMC ->SetCut("TOFMAX",tofmax); // time of flight cuts in seconds
1991 else if (strncmp(cut->GetName(),"TOFMAX",6) == 0) {
1992 fprintf(pFlukaVmcInp,"*\n*Time of flight cuts in seconds\n");
1993 fprintf(pFlukaVmcInp,"*Generated from call: SetCut('TOFMAX',tofmax);\n");
1996 // -6.0 = lower bound of the particle numbers for which the transport time cut-off and/or the start signal is to be applied
1997 // 64.0 = upper bound of the particle numbers for which the transport time cut-off and/or the start signal is to be applied
1998 fprintf(pFlukaVmcInp,"TIME-CUT %10.4g%10.1f%10.1f%10.1f%10.1f\n",cut->Cut()*1.e9,zero,zero,-6.0,64.0);
2002 cout << "SetCut for flag=" << cut->GetName() << " value=" << cut->Cut() << " not yet implemented!" << endl;
2005 cout << "SetCut for flag=" << cut->GetName() << " value=" << cut->Cut() << " (material specific) not yet implemented!" << endl;
2008 } //end of loop over SetCut calls
2010 // Add START and STOP card
2011 fprintf(pFlukaVmcInp,"START %10.1f\n",fEventsPerRun);
2012 fprintf(pFlukaVmcInp,"STOP \n");
2017 fclose(pFlukaVmcCoreInp);
2018 fclose(pFlukaVmcFlukaMat);
2019 fclose(pFlukaVmcInp);
2021 } // end of InitPhysics
2024 //______________________________________________________________________________
2025 void TFluka::SetMaxStep(Double_t step)
2027 // Set the maximum step size
2028 if (step > 1.e4) return;
2031 fGeom->GetCurrentRegion(mreg, latt);
2032 STEPSZ.stepmx[mreg - 1] = step;
2036 Double_t TFluka::MaxStep() const
2038 // Return the maximum for current medium
2040 fGeom->GetCurrentRegion(mreg, latt);
2041 return (STEPSZ.stepmx[mreg - 1]);
2044 //______________________________________________________________________________
2045 void TFluka::SetMaxNStep(Int_t)
2047 // SetMaxNStep is dummy procedure in TFluka !
2048 if (fVerbosityLevel >=3)
2049 cout << "SetMaxNStep is dummy procedure in TFluka !" << endl;
2052 //______________________________________________________________________________
2053 void TFluka::SetUserDecay(Int_t)
2055 // SetUserDecay is dummy procedure in TFluka !
2056 if (fVerbosityLevel >=3)
2057 cout << "SetUserDecay is dummy procedure in TFluka !" << endl;
2061 // dynamic properties
2063 //______________________________________________________________________________
2064 void TFluka::TrackPosition(TLorentzVector& position) const
2066 // Return the current position in the master reference frame of the
2067 // track being transported
2068 // TRACKR.atrack = age of the particle
2069 // TRACKR.xtrack = x-position of the last point
2070 // TRACKR.ytrack = y-position of the last point
2071 // TRACKR.ztrack = z-position of the last point
2072 Int_t caller = GetCaller();
2073 if (caller == 3 || caller == 6 || caller == 11 || caller == 12) { //bxdraw,endraw,usdraw
2074 position.SetX(GetXsco());
2075 position.SetY(GetYsco());
2076 position.SetZ(GetZsco());
2077 position.SetT(TRACKR.atrack);
2079 else if (caller == 4) { // mgdraw
2080 position.SetX(TRACKR.xtrack[TRACKR.ntrack]);
2081 position.SetY(TRACKR.ytrack[TRACKR.ntrack]);
2082 position.SetZ(TRACKR.ztrack[TRACKR.ntrack]);
2083 position.SetT(TRACKR.atrack);
2085 else if (caller == 5) { // sodraw
2086 position.SetX(TRACKR.xtrack[TRACKR.ntrack]);
2087 position.SetY(TRACKR.ytrack[TRACKR.ntrack]);
2088 position.SetZ(TRACKR.ztrack[TRACKR.ntrack]);
2092 Warning("TrackPosition","position not available");
2095 //______________________________________________________________________________
2096 void TFluka::TrackPosition(Double_t& x, Double_t& y, Double_t& z) const
2098 // Return the current position in the master reference frame of the
2099 // track being transported
2100 // TRACKR.atrack = age of the particle
2101 // TRACKR.xtrack = x-position of the last point
2102 // TRACKR.ytrack = y-position of the last point
2103 // TRACKR.ztrack = z-position of the last point
2104 Int_t caller = GetCaller();
2105 if (caller == 3 || caller == 6 || caller == 11 || caller == 12) { //bxdraw,endraw,usdraw
2110 else if (caller == 4 || caller == 5) { // mgdraw, sodraw
2111 x = TRACKR.xtrack[TRACKR.ntrack];
2112 y = TRACKR.ytrack[TRACKR.ntrack];
2113 z = TRACKR.ztrack[TRACKR.ntrack];
2116 Warning("TrackPosition","position not available");
2119 //______________________________________________________________________________
2120 void TFluka::TrackMomentum(TLorentzVector& momentum) const
2122 // Return the direction and the momentum (GeV/c) of the track
2123 // currently being transported
2124 // TRACKR.ptrack = momentum of the particle (not always defined, if
2125 // < 0 must be obtained from etrack)
2126 // TRACKR.cx,y,ztrck = direction cosines of the current particle
2127 // TRACKR.etrack = total energy of the particle
2128 // TRACKR.jtrack = identity number of the particle
2129 // PAPROP.am[TRACKR.jtrack] = particle mass in gev
2130 Int_t caller = GetCaller();
2131 if (caller != 2) { // not eedraw
2132 if (TRACKR.ptrack >= 0) {
2133 momentum.SetPx(TRACKR.ptrack*TRACKR.cxtrck);
2134 momentum.SetPy(TRACKR.ptrack*TRACKR.cytrck);
2135 momentum.SetPz(TRACKR.ptrack*TRACKR.cztrck);
2136 momentum.SetE(TRACKR.etrack);
2140 Double_t p = sqrt(TRACKR.etrack*TRACKR.etrack - PAPROP.am[TRACKR.jtrack+6]*PAPROP.am[TRACKR.jtrack+6]);
2141 momentum.SetPx(p*TRACKR.cxtrck);
2142 momentum.SetPy(p*TRACKR.cytrck);
2143 momentum.SetPz(p*TRACKR.cztrck);
2144 momentum.SetE(TRACKR.etrack);
2149 Warning("TrackMomentum","momentum not available");
2152 //______________________________________________________________________________
2153 void TFluka::TrackMomentum(Double_t& px, Double_t& py, Double_t& pz, Double_t& e) const
2155 // Return the direction and the momentum (GeV/c) of the track
2156 // currently being transported
2157 // TRACKR.ptrack = momentum of the particle (not always defined, if
2158 // < 0 must be obtained from etrack)
2159 // TRACKR.cx,y,ztrck = direction cosines of the current particle
2160 // TRACKR.etrack = total energy of the particle
2161 // TRACKR.jtrack = identity number of the particle
2162 // PAPROP.am[TRACKR.jtrack] = particle mass in gev
2163 Int_t caller = GetCaller();
2164 if (caller != 2) { // not eedraw
2165 if (TRACKR.ptrack >= 0) {
2166 px = TRACKR.ptrack*TRACKR.cxtrck;
2167 py = TRACKR.ptrack*TRACKR.cytrck;
2168 pz = TRACKR.ptrack*TRACKR.cztrck;
2173 Double_t p = sqrt(TRACKR.etrack*TRACKR.etrack - PAPROP.am[TRACKR.jtrack+6]*PAPROP.am[TRACKR.jtrack+6]);
2174 px = p*TRACKR.cxtrck;
2175 py = p*TRACKR.cytrck;
2176 pz = p*TRACKR.cztrck;
2182 Warning("TrackMomentum","momentum not available");
2185 //______________________________________________________________________________
2186 Double_t TFluka::TrackStep() const
2188 // Return the length in centimeters of the current step
2189 // TRACKR.ctrack = total curved path
2190 Int_t caller = GetCaller();
2191 if (caller == 11 || caller==12 || caller == 3 || caller == 6) //bxdraw,endraw,usdraw
2193 else if (caller == 4) //mgdraw
2194 return TRACKR.ctrack;
2199 //______________________________________________________________________________
2200 Double_t TFluka::TrackLength() const
2202 // TRACKR.cmtrck = cumulative curved path since particle birth
2203 Int_t caller = GetCaller();
2204 if (caller == 11 || caller==12 || caller == 3 || caller == 4 || caller == 6) //bxdraw,endraw,mgdraw,usdraw
2205 return TRACKR.cmtrck;
2210 //______________________________________________________________________________
2211 Double_t TFluka::TrackTime() const
2213 // Return the current time of flight of the track being transported
2214 // TRACKR.atrack = age of the particle
2215 Int_t caller = GetCaller();
2216 if (caller == 11 || caller==12 || caller == 3 || caller == 4 || caller == 6) //bxdraw,endraw,mgdraw,usdraw
2217 return TRACKR.atrack;
2222 //______________________________________________________________________________
2223 Double_t TFluka::Edep() const
2225 // Energy deposition
2226 // if TRACKR.ntrack = 0, TRACKR.mtrack = 0:
2227 // -->local energy deposition (the value and the point are not recorded in TRACKR)
2228 // but in the variable "rull" of the procedure "endraw.cxx"
2229 // if TRACKR.ntrack > 0, TRACKR.mtrack = 0:
2230 // -->no energy loss along the track
2231 // if TRACKR.ntrack > 0, TRACKR.mtrack > 0:
2232 // -->energy loss distributed along the track
2233 // TRACKR.dtrack = energy deposition of the jth deposition event
2235 // If coming from bxdraw we have 2 steps of 0 length and 0 edep
2236 Int_t caller = GetCaller();
2237 if (caller == 11 || caller==12) return 0.0;
2239 for ( Int_t j=0;j<TRACKR.mtrack;j++) {
2240 sum +=TRACKR.dtrack[j];
2242 if (TRACKR.ntrack == 0 && TRACKR.mtrack == 0)
2249 //______________________________________________________________________________
2250 Int_t TFluka::TrackPid() const
2252 // Return the id of the particle transported
2253 // TRACKR.jtrack = identity number of the particle
2254 Int_t caller = GetCaller();
2255 if (caller != 2) { // not eedraw
2256 return PDGFromId(TRACKR.jtrack);
2262 //______________________________________________________________________________
2263 Double_t TFluka::TrackCharge() const
2265 // Return charge of the track currently transported
2266 // PAPROP.ichrge = electric charge of the particle
2267 // TRACKR.jtrack = identity number of the particle
2268 Int_t caller = GetCaller();
2269 if (caller != 2) // not eedraw
2270 return PAPROP.ichrge[TRACKR.jtrack+6];
2275 //______________________________________________________________________________
2276 Double_t TFluka::TrackMass() const
2278 // PAPROP.am = particle mass in GeV
2279 // TRACKR.jtrack = identity number of the particle
2280 Int_t caller = GetCaller();
2281 if (caller != 2) // not eedraw
2282 return PAPROP.am[TRACKR.jtrack+6];
2287 //______________________________________________________________________________
2288 Double_t TFluka::Etot() const
2290 // TRACKR.etrack = total energy of the particle
2291 Int_t caller = GetCaller();
2292 if (caller != 2) // not eedraw
2293 return TRACKR.etrack;
2301 //______________________________________________________________________________
2302 Bool_t TFluka::IsNewTrack() const
2304 // Return true for the first call of Stepping()
2308 void TFluka::SetTrackIsNew(Bool_t flag)
2310 // Return true for the first call of Stepping()
2316 //______________________________________________________________________________
2317 Bool_t TFluka::IsTrackInside() const
2319 // True if the track is not at the boundary of the current volume
2320 // In Fluka a step is always inside one kind of material
2321 // If the step would go behind the region of one material,
2322 // it will be shortened to reach only the boundary.
2323 // Therefore IsTrackInside() is always true.
2324 Int_t caller = GetCaller();
2325 if (caller == 11 || caller==12) // bxdraw
2331 //______________________________________________________________________________
2332 Bool_t TFluka::IsTrackEntering() const
2334 // True if this is the first step of the track in the current volume
2336 Int_t caller = GetCaller();
2337 if (caller == 11) // bxdraw entering
2342 //______________________________________________________________________________
2343 Bool_t TFluka::IsTrackExiting() const
2345 // True if track is exiting volume
2347 Int_t caller = GetCaller();
2348 if (caller == 12) // bxdraw exiting
2353 //______________________________________________________________________________
2354 Bool_t TFluka::IsTrackOut() const
2356 // True if the track is out of the setup
2358 // Icode = 14: escape - call from Kaskad
2359 // Icode = 23: escape - call from Emfsco
2360 // Icode = 32: escape - call from Kasneu
2361 // Icode = 40: escape - call from Kashea
2362 // Icode = 51: escape - call from Kasoph
2367 fIcode == 51) return 1;
2371 //______________________________________________________________________________
2372 Bool_t TFluka::IsTrackDisappeared() const
2374 // means all inelastic interactions and decays
2375 // fIcode from usdraw
2376 if (fIcode == 101 || // inelastic interaction
2377 fIcode == 102 || // particle decay
2378 fIcode == 103 || // delta ray generation by hadron
2379 fIcode == 104 || // direct pair production
2380 fIcode == 105 || // bremsstrahlung (muon)
2381 fIcode == 208 || // bremsstrahlung (electron)
2382 fIcode == 214 || // in-flight annihilation
2383 fIcode == 215 || // annihilation at rest
2384 fIcode == 217 || // pair production
2385 fIcode == 219 || // Compton scattering
2386 fIcode == 221 || // Photoelectric effect
2387 fIcode == 300 || // hadronic interaction
2388 fIcode == 400 // delta-ray
2393 //______________________________________________________________________________
2394 Bool_t TFluka::IsTrackStop() const
2396 // True if the track energy has fallen below the threshold
2397 // means stopped by signal or below energy threshold
2398 // Icode = 12: stopping particle - call from Kaskad
2399 // Icode = 15: time kill - call from Kaskad
2400 // Icode = 21: below threshold, iarg=1 - call from Emfsco
2401 // Icode = 22: below threshold, iarg=2 - call from Emfsco
2402 // Icode = 24: time kill - call from Emfsco
2403 // Icode = 31: below threshold - call from Kasneu
2404 // Icode = 33: time kill - call from Kasneu
2405 // Icode = 41: time kill - call from Kashea
2406 // Icode = 52: time kill - call from Kasoph
2415 fIcode == 52) return 1;
2419 //______________________________________________________________________________
2420 Bool_t TFluka::IsTrackAlive() const
2422 // means not disappeared or not out
2423 if (IsTrackDisappeared() || IsTrackOut() ) return 0;
2431 //______________________________________________________________________________
2432 Int_t TFluka::NSecondaries() const
2435 // Number of secondary particles generated in the current step
2436 // FINUC.np = number of secondaries except light and heavy ions
2437 // FHEAVY.npheav = number of secondaries for light and heavy secondary ions
2438 Int_t caller = GetCaller();
2439 if (caller == 6) // valid only after usdraw
2440 return FINUC.np + FHEAVY.npheav;
2441 else if (caller == 50) {
2442 // Cerenkov Photon production
2446 } // end of NSecondaries
2448 //______________________________________________________________________________
2449 void TFluka::GetSecondary(Int_t isec, Int_t& particleId,
2450 TLorentzVector& position, TLorentzVector& momentum)
2452 // Copy particles from secondary stack to vmc stack
2455 Int_t caller = GetCaller();
2456 if (caller == 6) { // valid only after usdraw
2458 // Hadronic interaction
2459 if (isec >= 0 && isec < FINUC.np) {
2460 particleId = PDGFromId(FINUC.kpart[isec]);
2461 position.SetX(fXsco);
2462 position.SetY(fYsco);
2463 position.SetZ(fZsco);
2464 position.SetT(TRACKR.atrack);
2465 momentum.SetPx(FINUC.plr[isec]*FINUC.cxr[isec]);
2466 momentum.SetPy(FINUC.plr[isec]*FINUC.cyr[isec]);
2467 momentum.SetPz(FINUC.plr[isec]*FINUC.czr[isec]);
2468 momentum.SetE(FINUC.tki[isec] + PAPROP.am[FINUC.kpart[isec]+6]);
2470 else if (isec >= FINUC.np && isec < FINUC.np + FHEAVY.npheav) {
2471 Int_t jsec = isec - FINUC.np;
2472 particleId = FHEAVY.kheavy[jsec]; // this is Fluka id !!!
2473 position.SetX(fXsco);
2474 position.SetY(fYsco);
2475 position.SetZ(fZsco);
2476 position.SetT(TRACKR.atrack);
2477 momentum.SetPx(FHEAVY.pheavy[jsec]*FHEAVY.cxheav[jsec]);
2478 momentum.SetPy(FHEAVY.pheavy[jsec]*FHEAVY.cyheav[jsec]);
2479 momentum.SetPz(FHEAVY.pheavy[jsec]*FHEAVY.czheav[jsec]);
2480 if (FHEAVY.tkheav[jsec] >= 3 && FHEAVY.tkheav[jsec] <= 6)
2481 momentum.SetE(FHEAVY.tkheav[jsec] + PAPROP.am[jsec+6]);
2482 else if (FHEAVY.tkheav[jsec] > 6)
2483 momentum.SetE(FHEAVY.tkheav[jsec] + FHEAVY.amnhea[jsec]); // to be checked !!!
2486 Warning("GetSecondary","isec out of range");
2488 } else if (caller == 50) {
2489 Int_t index = OPPHST.lstopp - isec;
2490 position.SetX(OPPHST.xoptph[index]);
2491 position.SetY(OPPHST.yoptph[index]);
2492 position.SetZ(OPPHST.zoptph[index]);
2493 position.SetT(OPPHST.agopph[index]);
2494 Double_t p = OPPHST.poptph[index];
2496 momentum.SetPx(p * OPPHST.txopph[index]);
2497 momentum.SetPy(p * OPPHST.tyopph[index]);
2498 momentum.SetPz(p * OPPHST.tzopph[index]);
2502 Warning("GetSecondary","no secondaries available");
2504 } // end of GetSecondary
2507 //______________________________________________________________________________
2508 TMCProcess TFluka::ProdProcess(Int_t) const
2511 // Name of the process that has produced the secondary particles
2512 // in the current step
2514 Int_t mugamma = (TRACKR.jtrack == 7 || TRACKR.jtrack == 10 || TRACKR.jtrack == 11);
2516 if (fIcode == 102) return kPDecay;
2517 else if (fIcode == 104 || fIcode == 217) return kPPair;
2518 else if (fIcode == 219) return kPCompton;
2519 else if (fIcode == 221) return kPPhotoelectric;
2520 else if (fIcode == 105 || fIcode == 208) return kPBrem;
2521 else if (fIcode == 103 || fIcode == 400) return kPDeltaRay;
2522 else if (fIcode == 210 || fIcode == 212) return kPDeltaRay;
2523 else if (fIcode == 214 || fIcode == 215) return kPAnnihilation;
2524 else if (fIcode == 101) return kPHadronic;
2525 else if (fIcode == 101) {
2526 if (!mugamma) return kPHadronic;
2527 else if (TRACKR.jtrack == 7) return kPPhotoFission;
2528 else return kPMuonNuclear;
2530 else if (fIcode == 225) return kPRayleigh;
2531 // Fluka codes 100, 300 and 400 still to be investigasted
2532 else return kPNoProcess;
2536 Int_t TFluka::StepProcesses(TArrayI &proc) const
2539 // Return processes active in the current step
2563 iproc = kPLightAbsorption;
2566 iproc = kPLightRefraction;
2568 iproc = kPPhotoelectric;
2571 iproc = ProdProcess(0);
2576 //______________________________________________________________________________
2577 Int_t TFluka::VolId2Mate(Int_t id) const
2580 // Returns the material number for a given volume ID
2582 return fMCGeo->VolId2Mate(id);
2585 //______________________________________________________________________________
2586 const char* TFluka::VolName(Int_t id) const
2589 // Returns the volume name for a given volume ID
2591 return fMCGeo->VolName(id);
2594 //______________________________________________________________________________
2595 Int_t TFluka::VolId(const Text_t* volName) const
2598 // Converts from volume name to volume ID.
2599 // Time consuming. (Only used during set-up)
2600 // Could be replaced by hash-table
2602 return fMCGeo->VolId(volName);
2605 //______________________________________________________________________________
2606 Int_t TFluka::CurrentVolID(Int_t& copyNo) const
2609 // Return the logical id and copy number corresponding to the current fluka region
2611 if (gGeoManager->IsOutside()) return 0;
2612 TGeoNode *node = gGeoManager->GetCurrentNode();
2613 copyNo = node->GetNumber();
2614 Int_t id = node->GetVolume()->GetNumber();
2618 //______________________________________________________________________________
2619 Int_t TFluka::CurrentVolOffID(Int_t off, Int_t& copyNo) const
2622 // Return the logical id and copy number of off'th mother
2623 // corresponding to the current fluka region
2625 if (off<0 || off>gGeoManager->GetLevel()) return 0;
2626 if (off==0) return CurrentVolID(copyNo);
2627 TGeoNode *node = gGeoManager->GetMother(off);
2628 if (!node) return 0;
2629 copyNo = node->GetNumber();
2630 return node->GetVolume()->GetNumber();
2633 //______________________________________________________________________________
2634 const char* TFluka::CurrentVolName() const
2637 // Return the current volume name
2639 if (gGeoManager->IsOutside()) return 0;
2640 return gGeoManager->GetCurrentVolume()->GetName();
2643 //______________________________________________________________________________
2644 const char* TFluka::CurrentVolOffName(Int_t off) const
2647 // Return the volume name of the off'th mother of the current volume
2649 if (off<0 || off>gGeoManager->GetLevel()) return 0;
2650 if (off==0) return CurrentVolName();
2651 TGeoNode *node = gGeoManager->GetMother(off);
2652 if (!node) return 0;
2653 return node->GetVolume()->GetName();
2656 //______________________________________________________________________________
2657 Int_t TFluka::CurrentMaterial(Float_t & /*a*/, Float_t & /*z*/,
2658 Float_t & /*dens*/, Float_t & /*radl*/, Float_t & /*absl*/) const
2661 // Return the current medium number ??? what about material properties
2664 Int_t id = TFluka::CurrentVolID(copy);
2665 Int_t med = TFluka::VolId2Mate(id);
2669 //______________________________________________________________________________
2670 void TFluka::Gmtod(Float_t* xm, Float_t* xd, Int_t iflag)
2672 // Transforms a position from the world reference frame
2673 // to the current volume reference frame.
2675 // Geant3 desription:
2676 // ==================
2677 // Computes coordinates XD (in DRS)
2678 // from known coordinates XM in MRS
2679 // The local reference system can be initialized by
2680 // - the tracking routines and GMTOD used in GUSTEP
2681 // - a call to GMEDIA(XM,NUMED)
2682 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2683 // (inverse routine is GDTOM)
2685 // If IFLAG=1 convert coordinates
2686 // IFLAG=2 convert direction cosinus
2689 Double_t xmL[3], xdL[3];
2691 for (i=0;i<3;i++) xmL[i]=xm[i];
2692 if (iflag == 1) gGeoManager->MasterToLocal(xmL,xdL);
2693 else gGeoManager->MasterToLocalVect(xmL,xdL);
2694 for (i=0;i<3;i++) xd[i] = xdL[i];
2697 //______________________________________________________________________________
2698 void TFluka::Gmtod(Double_t* xm, Double_t* xd, Int_t iflag)
2700 if (iflag == 1) gGeoManager->MasterToLocal(xm,xd);
2701 else gGeoManager->MasterToLocalVect(xm,xd);
2704 //______________________________________________________________________________
2705 void TFluka::Gdtom(Float_t* xd, Float_t* xm, Int_t iflag)
2707 // Transforms a position from the current volume reference frame
2708 // to the world reference frame.
2710 // Geant3 desription:
2711 // ==================
2712 // Computes coordinates XM (Master Reference System
2713 // knowing the coordinates XD (Detector Ref System)
2714 // The local reference system can be initialized by
2715 // - the tracking routines and GDTOM used in GUSTEP
2716 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2717 // (inverse routine is GMTOD)
2719 // If IFLAG=1 convert coordinates
2720 // IFLAG=2 convert direction cosinus
2723 Double_t xmL[3], xdL[3];
2725 for (i=0;i<3;i++) xdL[i] = xd[i];
2726 if (iflag == 1) gGeoManager->LocalToMaster(xdL,xmL);
2727 else gGeoManager->LocalToMasterVect(xdL,xmL);
2728 for (i=0;i<3;i++) xm[i]=xmL[i];
2731 //______________________________________________________________________________
2732 void TFluka::Gdtom(Double_t* xd, Double_t* xm, Int_t iflag)
2734 if (iflag == 1) gGeoManager->LocalToMaster(xd,xm);
2735 else gGeoManager->LocalToMasterVect(xd,xm);
2738 //______________________________________________________________________________
2739 TObjArray *TFluka::GetFlukaMaterials()
2741 return fGeom->GetMatList();
2744 //______________________________________________________________________________
2745 void TFluka::SetMreg(Int_t l)
2747 // Set current fluka region
2748 fCurrentFlukaRegion = l;
2755 TString TFluka::ParticleName(Int_t pdg) const
2757 // Return particle name for particle with pdg code pdg.
2758 Int_t ifluka = IdFromPDG(pdg);
2759 return TString((CHPPRP.btype[ifluka+6]), 8);
2763 Double_t TFluka::ParticleMass(Int_t pdg) const
2765 // Return particle mass for particle with pdg code pdg.
2766 Int_t ifluka = IdFromPDG(pdg);
2767 return (PAPROP.am[ifluka+6]);
2770 Double_t TFluka::ParticleCharge(Int_t pdg) const
2772 // Return particle charge for particle with pdg code pdg.
2773 Int_t ifluka = IdFromPDG(pdg);
2774 return Double_t(PAPROP.ichrge[ifluka+6]);
2777 Double_t TFluka::ParticleLifeTime(Int_t pdg) const
2779 // Return particle lifetime for particle with pdg code pdg.
2780 Int_t ifluka = IdFromPDG(pdg);
2781 return (PAPROP.thalf[ifluka+6]);
2784 void TFluka::Gfpart(Int_t pdg, char* name, Int_t& type, Float_t& mass, Float_t& charge, Float_t& tlife)
2786 // Retrieve particle properties for particle with pdg code pdg.
2788 strcpy(name, ParticleName(pdg).Data());
2789 type = ParticleMCType(pdg);
2790 mass = ParticleMass(pdg);
2791 charge = ParticleCharge(pdg);
2792 tlife = ParticleLifeTime(pdg);
2797 #define pushcerenkovphoton pushcerenkovphoton_
2798 #define usersteppingckv usersteppingckv_
2802 void pushcerenkovphoton(Double_t & px, Double_t & py, Double_t & pz, Double_t & e,
2803 Double_t & vx, Double_t & vy, Double_t & vz, Double_t & tof,
2804 Double_t & polx, Double_t & poly, Double_t & polz, Double_t & wgt, Int_t& ntr)
2807 // Pushes one cerenkov photon to the stack
2810 TFluka* fluka = (TFluka*) gMC;
2811 TVirtualMCStack* cppstack = fluka->GetStack();
2812 Int_t parent = TRACKR.ispusr[mkbmx2-1];
2813 cppstack->PushTrack(0, parent, 50000050,
2817 kPCerenkov, ntr, wgt, 0);
2820 void usersteppingckv(Int_t & nphot, Int_t & mreg, Double_t & x, Double_t & y, Double_t & z)
2823 // Calls stepping in order to signal cerenkov production
2825 TFluka *fluka = (TFluka*)gMC;
2826 fluka->SetMreg(mreg);
2830 fluka->SetNCerenkov(nphot);
2831 fluka->SetCaller(50);
2832 printf("userstepping ckv: %10d %10d %13.3f %13.3f %13.2f\n", nphot, mreg, x, y, z);
2833 (TVirtualMCApplication::Instance())->Stepping();