3 // Configuration for the Geant4 production 2010
7 #if !defined(__CINT__) || defined(__MAKECINT__)
12 #include <TVirtualMC.h>
13 #include <TGeant3TGeo.h>
14 #include "STEER/AliRunLoader.h"
15 #include "STEER/AliRun.h"
16 #include "STEER/AliConfig.h"
17 #include "PYTHIA6/AliDecayerPythia.h"
18 #include "PYTHIA6/AliGenPythia.h"
19 #include "TDPMjet/AliGenDPMjet.h"
20 #include "STEER/AliMagFCheb.h"
21 #include "STRUCT/AliBODY.h"
22 #include "STRUCT/AliMAG.h"
23 #include "STRUCT/AliABSOv3.h"
24 #include "STRUCT/AliDIPOv3.h"
25 #include "STRUCT/AliHALLv3.h"
26 #include "STRUCT/AliFRAMEv2.h"
27 #include "STRUCT/AliSHILv3.h"
28 #include "STRUCT/AliPIPEv3.h"
29 #include "ITS/AliITSv11.h"
30 #include "TPC/AliTPCv2.h"
31 #include "TOF/AliTOFv6T0.h"
32 #include "HMPID/AliHMPIDv3.h"
33 #include "ZDC/AliZDCv3.h"
34 #include "TRD/AliTRDv1.h"
35 #include "TRD/AliTRDgeometry.h"
36 #include "FMD/AliFMDv1.h"
37 #include "MUON/AliMUONv1.h"
38 #include "PHOS/AliPHOSv1.h"
39 #include "PHOS/AliPHOSSimParam.h"
40 #include "PMD/AliPMDv1.h"
41 #include "T0/AliT0v1.h"
42 #include "EMCAL/AliEMCALv2.h"
43 #include "ACORDE/AliACORDEv1.h"
44 #include "VZERO/AliVZEROv7.h"
50 kPythia6, kPythia6D6T, kPythia6ATLAS, kPythia6ATLAS_Flat, kPythiaPerugia0, kPhojet, kRunMax
53 const char * pprRunName[] = {
54 "kPythia6", "kPythia6D6T", "kPythia6ATLAS", "kPythia6ATLAS_Flat", "kPythiaPerugia0", "kPhojet"
59 kNoField, k5kG, kFieldMax
62 const char * pprField[] = {
68 QGSP_BERT_EMV, CHIPS, QGSP_BERT_CHIPS, QGSP_FTFP_BERT, FTFP_BERT, FTFP_BERT_EMV,
69 QGSP_BERT_EMV_OPTICAL, CHIPS_OPTICAL, QGSP_BERT_CHIPS_OPTICAL, QGSP_FTFP_BERT_OPTICAL, FTFP_BERT_OPTICAL, FTFP_BERT_EMV_OPTICAL, kListMax
72 const char * physicsListName[] = {
73 "QGSP_BERT_EMV", "CHIPS", "QGSP_BERT_CHIPS", "QGSP_FTFP_BERT", "FTFP_BERT", "FTFP_BERT_EMV",
74 "QGSP_BERT_EMV_OPTICAL", "CHIPS_OPTICAL", "QGSP_BERT_CHIPS_OPTICAL", "QGSP_FTFP_BERT_OPTICAL", "FTFP_BERT_OPTICAL", "FTFP_BERT_EMV_OPTICAL",
79 kDefaultPPTrig, kDefaultPbPbTrig
82 const char * pprTrigConfName[] = {
88 AliGenerator *MbPythia();
89 AliGenerator *MbPythiaTuneD6T();
90 AliGenerator *MbPhojet();
91 void ProcessEnvironmentVars();
93 // Geterator, field, beam energy
94 static PDC06Proc_t proc = kPhojet;
95 static Mag_t mag = k5kG;
96 static Float_t energy = 10000; // energy in CMS
97 static Int_t runNumber = 0;
98 static PprTrigConf_t strig = kDefaultPPTrig; // default pp trigger configuration
99 static PhysicsList_t physicslist = QGSP_BERT_EMV;
101 //========================//
102 // Set Random Number seed //
103 //========================//
105 static UInt_t seed = dt.Get();
108 static TString comment;
114 // Get settings from environment variables
115 ProcessEnvironmentVars();
117 gRandom->SetSeed(seed);
118 cerr<<"Seed for random number generation= "<<seed<<endl;
120 // Libraries required by geant321
121 #if defined(__CINT__)
122 gSystem->Load("liblhapdf"); // Parton density functions
123 gSystem->Load("libEGPythia6"); // TGenerator interface
124 if (proc == kPythia6 || proc == kPhojet) {
125 gSystem->Load("libpythia6"); // Pythia 6.2
127 gSystem->Load("libpythia6.4.21"); // Pythia 6.4
129 gSystem->Load("libAliPythia6"); // ALICE specific implementations
130 // gSystem->Load("libgeant321");
133 // new TGeant3TGeo("C++ Interface to Geant3");
135 //=======================================================================
136 // Create the output file
139 AliRunLoader* rl=0x0;
141 cout<<"Config.C: Creating Run Loader ..."<<endl;
142 rl = AliRunLoader::Open("galice.root",
143 AliConfig::GetDefaultEventFolderName(),
147 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
150 rl->SetCompressionLevel(2);
151 rl->SetNumberOfEventsPerFile(1000);
152 gAlice->SetRunLoader(rl);
153 // gAlice->SetGeometryFromFile("geometry.root");
154 // gAlice->SetGeometryFromCDB();
156 // Set the trigger configuration: proton-proton
158 AliSimulation::Instance()->SetTriggerConfig(pprTrigConfName[strig]);
159 cout <<"Trigger configuration is set to "<<pprTrigConfName[strig]<<endl;
187 //=================== Alice BODY parameters =============================
188 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
193 //=================== MAG parameters ============================
194 // --- Start with Magnet since detector layouts may be depending ---
195 // --- on the selected Magnet dimensions ---
196 AliMAG *MAG = new AliMAG("MAG", "Magnet");
202 //=================== ABSO parameters ============================
203 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
208 //=================== DIPO parameters ============================
210 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
215 //=================== HALL parameters ============================
217 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
223 //=================== FRAME parameters ============================
225 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
231 //=================== SHIL parameters ============================
233 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
239 //=================== PIPE parameters ============================
241 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
246 //=================== ITS parameters ============================
248 AliITS *ITS = new AliITSv11("ITS","ITS v11");
253 //============================ TPC parameters =====================
255 AliTPC *TPC = new AliTPCv2("TPC", "Default");
256 TPC->SetPrimaryIonisation();// not used with Geant3
261 //=================== TOF parameters ============================
263 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
269 //=================== HMPID parameters ===========================
271 AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
278 //=================== ZDC parameters ============================
280 AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
281 //Collimators aperture
282 ZDC->SetVCollSideCAperture(0.85);
283 ZDC->SetVCollSideCCentre(0.);
284 ZDC->SetVCollSideAAperture(0.75);
285 ZDC->SetVCollSideACentre(0.);
295 //=================== TRD parameters ============================
297 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
298 AliTRDgeometry *geoTRD = TRD->GetGeometry();
299 // Partial geometry: modules at 0,1,7,8,9,16,17
300 // starting at 3h in positive direction
301 geoTRD->SetSMstatus(2,0);
302 geoTRD->SetSMstatus(3,0);
303 geoTRD->SetSMstatus(4,0);
304 geoTRD->SetSMstatus(5,0);
305 geoTRD->SetSMstatus(6,0);
306 geoTRD->SetSMstatus(11,0);
307 geoTRD->SetSMstatus(12,0);
308 geoTRD->SetSMstatus(13,0);
309 geoTRD->SetSMstatus(14,0);
310 geoTRD->SetSMstatus(15,0);
311 geoTRD->SetSMstatus(16,0);
316 //=================== FMD parameters ============================
318 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
323 //=================== MUON parameters ===========================
324 // New MUONv1 version (geometry defined via builders)
325 AliMUON *MUON = new AliMUONv1("MUON", "default");
326 // activate trigger efficiency by cells
327 MUON->SetTriggerEffCells(1);
332 //=================== PHOS parameters ===========================
334 AliPHOS *PHOS = new AliPHOSv1("PHOS", "noCPV_Modules123");
341 //=================== PMD parameters ============================
343 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
348 //=================== T0 parameters ============================
349 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
354 //=================== EMCAL parameters ============================
356 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_FIRSTYEARV1");
361 //=================== ACORDE parameters ============================
363 AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
368 //=================== ACORDE parameters ============================
370 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
374 // Load Geant4 + Geant4 VMC libraries
376 if (gClassTable->GetID("TGeant4") == -1) {
377 TString g4libsMacro = "$G4INSTALL/macro/g4libs.C";
378 //TString g4libsMacro = "$ALICE/geant4_vmc/examples/macro/g4libs.C";
379 //Load Geant4 libraries
380 if (!gInterpreter->IsLoaded(g4libsMacro.Data())) {
381 gROOT->LoadMacro(g4libsMacro.Data());
382 gInterpreter->ProcessLine("g4libs()");
391 TG4RunConfiguration* runConfiguration=0x0;
392 for (Int_t iList = 0; iList < kListMax; iList++) {
393 if(iList<kListMax/2){
394 if(physicslist == iList){
396 new TG4RunConfiguration("geomRoot",
397 physicsListName[iList],
398 "specialCuts+stackPopper+stepLimiter",
402 else if(iList>=kListMax/2){//add "optical" PL to HadronPhysicsList
403 if(physicslist == iList){
405 new TG4RunConfiguration("geomRoot",
406 Form("%s+optical",physicsListName[iList-kListMax/2]),
407 "specialCuts+stackPopper+stepLimiter",
412 geant4 = new TGeant4("TGeant4", "The Geant4 Monte Carlo", runConfiguration);
413 cout << "Geant4 has been created." << endl;
416 cout << "Monte Carlo has been already created." << endl;
421 // Customization of Geant4 VMC
423 geant4->ProcessGeantCommand("/mcVerbose/all 1");
424 geant4->ProcessGeantCommand("/mcVerbose/geometryManager 1");
425 geant4->ProcessGeantCommand("/mcVerbose/opGeometryManager 1");
426 geant4->ProcessGeantCommand("/mcTracking/loopVerbose 1");
427 geant4->ProcessGeantCommand("/mcPhysics/rangeCuts 0.01 mm");
428 // for Geant4 <= 9.4.p03
429 //geant4->ProcessGeantCommand("/mcPhysics/selectOpProcess Scintillation");
430 //geant4->ProcessGeantCommand("/mcPhysics/setOpProcessActivation false");
432 geant4->ProcessGeantCommand("/optics_engine/selectOpProcess Scintillation");
433 geant4->ProcessGeantCommand("/optics_engine/setOpProcessUse false");
434 geant4->ProcessGeantCommand("/optics_engine/selectOpProcess OpWLS");
435 geant4->ProcessGeantCommand("/optics_engine/setOpProcessUse false");
436 geant4->ProcessGeantCommand("/optics_engine/selectOpProcess OpMieHG");
437 geant4->ProcessGeantCommand("/optics_engine/setOpProcessUse false");
439 geant4->ProcessGeantCommand("/mcVerbose/composedPhysicsList 2");
440 geant4->ProcessGeantCommand("/mcTracking/skipNeutrino true");
441 // geant4->ProcessGeantCommand("/mcDet/setMaxStepInLowDensityMaterials 1 cm");
445 //=======================================================================
446 // ************* STEERING parameters FOR ALICE SIMULATION **************
447 // --- Specify event type to be tracked through the ALICE setup
448 // --- All positions are in cm, angles in degrees, and P and E in GeV
451 gMC->SetProcess("DCAY",1);
452 gMC->SetProcess("PAIR",1);
453 gMC->SetProcess("COMP",1);
454 gMC->SetProcess("PHOT",1);
455 gMC->SetProcess("PFIS",0);
456 gMC->SetProcess("DRAY",0);
457 gMC->SetProcess("ANNI",1);
458 gMC->SetProcess("BREM",1);
459 gMC->SetProcess("MUNU",1);
460 gMC->SetProcess("CKOV",1);
461 gMC->SetProcess("HADR",1);
462 gMC->SetProcess("LOSS",2);
463 gMC->SetProcess("MULS",1);
464 gMC->SetProcess("RAYL",1);
466 Float_t cut = 1.e-3; // 1MeV cut by default
467 Float_t tofmax = 1.e10;
469 gMC->SetCut("CUTGAM", cut);
470 gMC->SetCut("CUTELE", cut);
471 gMC->SetCut("CUTNEU", cut);
472 gMC->SetCut("CUTHAD", cut);
473 gMC->SetCut("CUTMUO", cut);
474 gMC->SetCut("BCUTE", cut);
475 gMC->SetCut("BCUTM", cut);
476 gMC->SetCut("DCUTE", cut);
477 gMC->SetCut("DCUTM", cut);
478 gMC->SetCut("PPCUTM", cut);
479 gMC->SetCut("TOFMAX", tofmax);
484 //======================//
485 // Set External decayer //
486 //======================//
487 TVirtualMCDecayer* decayer = new AliDecayerPythia();
488 decayer->SetForceDecay(kAll);
490 gMC->SetExternalDecayer(decayer);
492 //=========================//
493 // Generator Configuration //
494 //=========================//
495 AliGenerator* gener = 0x0;
497 if (proc == kPythia6) {
499 } else if (proc == kPythia6D6T) {
500 gener = MbPythiaTuneD6T();
501 } else if (proc == kPythia6ATLAS) {
502 gener = MbPythiaTuneATLAS();
503 } else if (proc == kPythiaPerugia0) {
504 gener = MbPythiaTunePerugia0();
505 } else if (proc == kPythia6ATLAS_Flat) {
506 gener = MbPythiaTuneATLAS_Flat();
507 } else if (proc == kPhojet) {
514 // Size of the interaction diamond
516 Float_t sigmaz = 5.4 / TMath::Sqrt(2.); // [cm]
518 //sigmaz = 10.5 / TMath::Sqrt(2.); // [cm]
521 sigmaz = 6.3 / TMath::Sqrt(2.); // [cm]
527 Float_t betast = 10.0; // beta* [m]
528 if (runNumber >= 117048) betast = 2.0;
529 if (runNumber > 122375) betast = 3.5; // starting with fill 1179
532 Float_t eps = 5.0e-6; // emittance [m]
533 Float_t gamma = energy / 2.0 / 0.938272; // relativistic gamma [1]
534 Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.; // [cm]
535 printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
537 gener->SetSigma(sigmaxy, sigmaxy, sigmaz); // Sigma in (X,Y,Z) (cm) on IP position
538 gener->SetVertexSmear(kPerEvent);
541 printf("\n \n Comment: %s \n \n", comment.Data());
549 AliGenerator* MbPythia()
551 comment = comment.Append(" pp: Pythia low-pt");
554 AliGenPythia* pythia = new AliGenPythia(-1);
555 pythia->SetMomentumRange(0, 999999.);
556 pythia->SetThetaRange(0., 180.);
557 pythia->SetYRange(-12.,12.);
558 pythia->SetPtRange(0,1000.);
559 pythia->SetProcess(kPyMb);
560 pythia->SetEnergyCMS(energy);
565 AliGenerator* MbPythiaTuneD6T()
567 comment = comment.Append(" pp: Pythia low-pt");
570 AliGenPythia* pythia = new AliGenPythia(-1);
571 pythia->SetMomentumRange(0, 999999.);
572 pythia->SetThetaRange(0., 180.);
573 pythia->SetYRange(-12.,12.);
574 pythia->SetPtRange(0,1000.);
575 pythia->SetProcess(kPyMb);
576 pythia->SetEnergyCMS(energy);
578 // 109 D6T : Rick Field's CDF Tune D6T (NB: needs CTEQ6L pdfs externally)
579 pythia->SetTune(109); // F I X
580 pythia->SetStrucFunc(kCTEQ6l);
585 AliGenerator* MbPythiaTunePerugia0()
587 comment = comment.Append(" pp: Pythia low-pt (Perugia0)");
590 AliGenPythia* pythia = new AliGenPythia(-1);
591 pythia->SetMomentumRange(0, 999999.);
592 pythia->SetThetaRange(0., 180.);
593 pythia->SetYRange(-12.,12.);
594 pythia->SetPtRange(0,1000.);
595 pythia->SetProcess(kPyMb);
596 pythia->SetEnergyCMS(energy);
599 pythia->SetTune(320);
600 pythia->UseNewMultipleInteractionsScenario();
601 pythia->SetCrossingAngle(0,0.000280);
608 AliGenerator* MbPythiaTuneATLAS()
610 comment = comment.Append(" pp: Pythia low-pt");
613 AliGenPythia* pythia = new AliGenPythia(-1);
614 pythia->SetMomentumRange(0, 999999.);
615 pythia->SetThetaRange(0., 180.);
616 pythia->SetYRange(-12.,12.);
617 pythia->SetPtRange(0,1000.);
618 pythia->SetProcess(kPyMb);
619 pythia->SetEnergyCMS(energy);
621 // C 306 ATLAS-CSC: Arthur Moraes' (new) ATLAS tune (needs CTEQ6L externally)
622 pythia->SetTune(306);
623 pythia->SetStrucFunc(kCTEQ6l);
628 AliGenerator* MbPythiaTuneATLAS_Flat()
630 AliGenPythia* pythia = MbPythiaTuneATLAS();
632 comment = comment.Append("; flat multiplicity distribution");
634 // set high multiplicity trigger
635 // this weight achieves a flat multiplicity distribution
638 0.234836, 0.103484, 0.00984802, 0.0199906, 0.0260018, 0.0208481, 0.0101477, 0.00146998, -0.00681513, -0.0114928,
639 -0.0113352, -0.0102012, -0.00895238, -0.00534961, -0.00261648, -0.000819048, 0.00130816, 0.00177978, 0.00373838, 0.00566255,
640 0.00628156, 0.00687458, 0.00668017, 0.00702917, 0.00810848, 0.00876167, 0.00832783, 0.00848518, 0.0107709, 0.0106849,
641 0.00933702, 0.00912525, 0.0106553, 0.0102785, 0.0101756, 0.010962, 0.00957103, 0.00970448, 0.0117133, 0.012271,
642 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113,
643 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113,
644 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113,
645 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113,
646 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113,
647 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113, 0.0113,
648 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
649 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
654 0.00168216, 0.000743548, 0.00103125, 0.00108605, 0.00117101, 0.00124577, 0.00129119, 0.00128341, 0.00121421, 0.00112431,
655 0.00100588, 0.000895078, 0.000790314, 0.000711673, 0.000634547, 0.000589133, 0.000542763, 0.000509552, 0.000487375, 0.000468906,
656 0.000460196, 0.000455806, 0.00044843, 0.000449317, 0.00045007, 0.000458016, 0.000461167, 0.000474742, 0.00050161, 0.00051637,
657 0.000542336, 0.000558854, 0.000599169, 0.000611982, 0.000663855, 0.000696322, 0.000722825, 0.000771686, 0.000838023, 0.000908317,
658 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003,
659 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003,
660 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003,
661 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003,
662 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003,
663 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003, 0.0003,
664 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
665 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
668 TH1F *weight = new TH1F("newweight","newweight",120,-0.5,119.5);
670 weight->SetContent(cont);
671 weight->SetError(err);
673 Int_t limit = weight->GetRandom();
674 pythia->SetTriggerChargedMultiplicity(limit, 1.4);
676 comment = comment.Append(Form("; multiplicity threshold set to %d in |eta| < 1.4", limit));
681 AliGenerator* MbPhojet()
683 comment = comment.Append(" pp: Pythia low-pt");
686 #if defined(__CINT__)
687 gSystem->Load("libdpmjet"); // Parton density functions
688 gSystem->Load("libTDPMjet"); // Parton density functions
690 AliGenDPMjet* dpmjet = new AliGenDPMjet(-1);
691 dpmjet->SetMomentumRange(0, 999999.);
692 dpmjet->SetThetaRange(0., 180.);
693 dpmjet->SetYRange(-12.,12.);
694 dpmjet->SetPtRange(0,1000.);
695 dpmjet->SetProcess(kDpmMb);
696 dpmjet->SetEnergyCMS(energy);
697 dpmjet->SetCrossingAngle(0,0.000280);
701 void ProcessEnvironmentVars()
704 if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
705 for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
706 if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
707 proc = (PDC06Proc_t)iRun;
708 cout<<"Run type set to "<<pprRunName[iRun]<<endl;
714 if (gSystem->Getenv("CONFIG_FIELD")) {
715 for (Int_t iField = 0; iField < kFieldMax; iField++) {
716 if (strcmp(gSystem->Getenv("CONFIG_FIELD"), pprField[iField])==0) {
718 cout<<"Field set to "<<pprField[iField]<<endl;
724 if (gSystem->Getenv("CONFIG_ENERGY")) {
725 energy = atoi(gSystem->Getenv("CONFIG_ENERGY"));
726 cout<<"Energy set to "<<energy<<" GeV"<<endl;
729 // Random Number seed
730 if (gSystem->Getenv("CONFIG_SEED")) {
731 seed = atoi(gSystem->Getenv("CONFIG_SEED"));
735 if (gSystem->Getenv("DC_RUN")) {
736 runNumber = atoi(gSystem->Getenv("DC_RUN"));
740 if (gSystem->Getenv("CONFIG_PHYSICSLIST")) {
741 for (Int_t iList = 0; iList < kListMax; iList++) {
742 if (strcmp(gSystem->Getenv("CONFIG_PHYSICSLIST"), physicsListName[iList])==0){
743 physicslist = (PhysicsList_t)iList;
744 cout<<"Physics list set to "<< physicsListName[iList]<<endl;