[988] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
| 3 | // * License and Disclaimer * |
---|
| 4 | // * * |
---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
| 7 | // * conditions of the Geant4 Software License, included in the file * |
---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
| 9 | // * include a list of copyright holders. * |
---|
| 10 | // * * |
---|
| 11 | // * Neither the authors of this software system, nor their employing * |
---|
| 12 | // * institutes,nor the agencies providing financial support for this * |
---|
| 13 | // * work make any representation or warranty, express or implied, * |
---|
| 14 | // * regarding this software system or assume any liability for its * |
---|
| 15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
| 16 | // * for the full disclaimer and the limitation of liability. * |
---|
| 17 | // * * |
---|
| 18 | // * This code implementation is the result of the scientific and * |
---|
| 19 | // * technical work of the GEANT4 collaboration. * |
---|
| 20 | // * By using, copying, modifying or distributing the software (or * |
---|
| 21 | // * any work based on the software) you agree to acknowledge its * |
---|
| 22 | // * use in resulting scientific publications, and indicate your * |
---|
| 23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 24 | // ******************************************************************** |
---|
| 25 | // |
---|
| 26 | // $Id: G4PhysListFactory.cc,v 1.6 2008/11/25 15:36:19 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-02-ref-02 $ |
---|
| 28 | // |
---|
| 29 | //--------------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // ClassName: G4PhysListFactory |
---|
| 32 | // |
---|
| 33 | // Author: 21 April 2008 V. Ivanchenko |
---|
| 34 | // |
---|
| 35 | // Modified: |
---|
| 36 | // |
---|
| 37 | //---------------------------------------------------------------------------- |
---|
| 38 | // |
---|
| 39 | |
---|
| 40 | #include "G4PhysListFactory.hh" |
---|
| 41 | #include "FTFC.hh" |
---|
| 42 | #include "FTFP.hh" |
---|
| 43 | #include "FTFP_BERT.hh" |
---|
| 44 | #include "FTFP_EMV.hh" |
---|
| 45 | #include "FTF_BIC.hh" |
---|
| 46 | #include "LBE.hh" |
---|
| 47 | #include "LHEP.hh" |
---|
| 48 | #include "LHEP_BERT.hh" |
---|
| 49 | #include "LHEP_BERT_HP.hh" |
---|
| 50 | #include "LHEP_EMV.hh" |
---|
| 51 | #include "LHEP_PRECO_HP.hh" |
---|
| 52 | #include "QBBC.hh" |
---|
| 53 | #include "QGSC.hh" |
---|
| 54 | #include "QGSC_BERT.hh" |
---|
| 55 | #include "QGSC_EFLOW.hh" |
---|
| 56 | #include "QGSC_EMV.hh" |
---|
| 57 | #include "QGSP.hh" |
---|
| 58 | #include "QGSP_BERT.hh" |
---|
| 59 | #include "QGSP_BERT_DIF.hh" |
---|
| 60 | #include "QGSP_BERT_EMV.hh" |
---|
| 61 | #include "QGSP_BERT_HP.hh" |
---|
| 62 | #include "QGSP_BERT_NQE.hh" |
---|
| 63 | #include "QGSP_BERT_TRV.hh" |
---|
| 64 | #include "QGSP_BIC.hh" |
---|
| 65 | #include "QGSP_BIC_HP.hh" |
---|
| 66 | #include "QGSP_DIF.hh" |
---|
| 67 | #include "QGSP_EMV.hh" |
---|
| 68 | #include "QGSP_EMV_NQE.hh" |
---|
| 69 | #include "QGSP_EMX.hh" |
---|
| 70 | #include "QGSP_NQE.hh" |
---|
| 71 | #include "QGSP_QEL.hh" |
---|
| 72 | #include "QGS_BIC.hh" |
---|
| 73 | |
---|
| 74 | G4PhysListFactory::G4PhysListFactory() |
---|
| 75 | { |
---|
| 76 | defName = "QGSP_BERT"; |
---|
| 77 | nlists = 35; |
---|
| 78 | G4String s[35] = { |
---|
| 79 | "FTFC","FTFP","FTFP_BERY","FTFP_EMV","FTF_BIC", |
---|
| 80 | "LBE","LHEP","LHEP_BERT","LHEP_EMV","LHEP_PRECO_HP" |
---|
| 81 | "QBBB","QBBC","QBBCG","QBBCF","QBBC_HP","QGSC", |
---|
| 82 | "QGSC_BERT","QGSC_EFLOW","QGSC_EMV","QGSP","QGSP_BERT", |
---|
| 83 | "QGSP_BERT_DIF","QGSP_BERT_EMV","QGSP_BERT_HP","QGSP_BERT_NQE","QGSP_BERT_TRV", |
---|
| 84 | "QGSP_BIC","QGSP_BIC_HP","QGSP_DIF","QGSP_EMV","QGSP_EMV_NQE", |
---|
| 85 | "QGSP_EMX","QGSP_NQE","QGSP_QEL","QGS_BIC"}; |
---|
| 86 | |
---|
| 87 | for(size_t i=0; i<nlists; i++) { |
---|
| 88 | listnames.push_back(s[i]); |
---|
| 89 | } |
---|
| 90 | } |
---|
| 91 | |
---|
| 92 | G4PhysListFactory::~G4PhysListFactory() |
---|
| 93 | {} |
---|
| 94 | |
---|
| 95 | G4VModularPhysicsList* G4PhysListFactory::ReferencePhysList() |
---|
| 96 | { |
---|
| 97 | // instantiate PhysList by environment variable "PHYSLIST" |
---|
| 98 | G4String name = ""; |
---|
| 99 | char* path = getenv("PHYSLIST"); |
---|
| 100 | if (path) { |
---|
| 101 | name = G4String(path); |
---|
| 102 | } else { |
---|
| 103 | name = defName; |
---|
| 104 | G4cout << "### G4PhysListFactory WARNING: " |
---|
| 105 | << " environment variable PHYSLIST is not defined" |
---|
| 106 | << G4endl |
---|
| 107 | << " Default Physics Lists " << name |
---|
| 108 | << " is instantiated" |
---|
| 109 | << G4endl; |
---|
| 110 | } |
---|
| 111 | return GetReferencePhysList(name); |
---|
| 112 | } |
---|
| 113 | |
---|
| 114 | G4VModularPhysicsList* G4PhysListFactory::GetReferencePhysList( |
---|
| 115 | const G4String& name) |
---|
| 116 | { |
---|
| 117 | G4VModularPhysicsList* p = 0; |
---|
| 118 | if (name == "FTFC") {p = new FTFC();} |
---|
| 119 | else if(name == "FTFP") {p = new FTFP();} |
---|
| 120 | else if(name == "FTFP_BERT") {p = new FTFP_BERT();} |
---|
| 121 | else if(name == "FTFP_EMV") {p = new FTFP_EMV();} |
---|
| 122 | else if(name == "FTF_BIC") {p = new FTF_BIC();} |
---|
| 123 | else if(name == "LBE") {p = new LBE();} |
---|
| 124 | else if(name == "LHEP") {p = new LHEP();} |
---|
| 125 | else if(name == "LHEP_BERT") {p = new LHEP_BERT();} |
---|
| 126 | else if(name == "LHEP_EMV") {p = new LHEP_EMV();} |
---|
| 127 | else if(name == "LHEP_PRECO_HP") {p = new LHEP_PRECO_HP();} |
---|
| 128 | else if(name == "QBBBG") {p = new QBBC(1, "QBBBG");} |
---|
| 129 | else if(name == "QBBC") {p = new QBBC();} |
---|
| 130 | else if(name == "QBBCG") {p = new QBBC(1, "QBBCG");} |
---|
| 131 | else if(name == "QBBCF") {p = new QBBC(1, "QBBCF");} |
---|
| 132 | else if(name == "QBBC_HP") {p = new QBBC(1, "QBBC_HP");} |
---|
| 133 | else if(name == "QGSC") {p = new QGSC();} |
---|
| 134 | else if(name == "QGSC_BERT") {p = new QGSC_BERT();} |
---|
| 135 | else if(name == "QGSC_EFLOW") {p = new QGSC_EFLOW();} |
---|
| 136 | else if(name == "QGSC_EMV") {p = new QGSC_EMV();} |
---|
| 137 | else if(name == "QGSP") {p = new QGSP();} |
---|
| 138 | else if(name == "QGSP_BERT") {p = new QGSP_BERT();} |
---|
| 139 | else if(name == "QGSP_BERT_DIF") {p = new QGSP_BERT_DIF();} |
---|
| 140 | else if(name == "QGSP_BERT_EMV") {p = new QGSP_BERT_EMV();} |
---|
| 141 | else if(name == "QGSP_BERT_HP") {p = new QGSP_BERT_HP();} |
---|
| 142 | else if(name == "QGSP_BERT_NQE") {p = new QGSP_BERT_NQE();} |
---|
| 143 | else if(name == "QGSP_BERT_TRV") {p = new QGSP_BERT_TRV();} |
---|
| 144 | else if(name == "QGSP_BIC") {p = new QGSP_BIC();} |
---|
| 145 | else if(name == "QGSP_BIC_HP") {p = new QGSP_BIC_HP();} |
---|
| 146 | else if(name == "QGSP_DIF") {p = new QGSP_DIF();} |
---|
| 147 | else if(name == "QGSP_EMV") {p = new QGSP_EMV();} |
---|
| 148 | else if(name == "QGSP_EMV_NQE") {p = new QGSP_EMV_NQE();} |
---|
| 149 | else if(name == "QGSP_EMX") {p = new QGSP_EMX();} |
---|
| 150 | else if(name == "QGSP_NQE") {p = new QGSP_NQE();} |
---|
| 151 | else if(name == "QGSP_QEL") {p = new QGSP_QEL();} |
---|
| 152 | else if(name == "QGS_BIC") {p = new QGS_BIC();} |
---|
| 153 | else { |
---|
| 154 | G4cout << "### G4PhysListFactory WARNING: " |
---|
| 155 | << "PhysicsList " << name << " is not known" |
---|
| 156 | << G4endl |
---|
| 157 | << " Default Physics Lists " << defName |
---|
| 158 | << " is instantiated" |
---|
| 159 | << G4endl; |
---|
| 160 | p = new QGSP_BERT(); |
---|
| 161 | } |
---|
| 162 | return p; |
---|
| 163 | } |
---|
| 164 | |
---|
| 165 | G4bool G4PhysListFactory::IsReferencePhysList(const G4String& name) |
---|
| 166 | { |
---|
| 167 | G4bool res = false; |
---|
| 168 | for(size_t i=0; i<nlists; i++) { |
---|
| 169 | if(name == listnames[i]) { |
---|
| 170 | res = true; |
---|
| 171 | break; |
---|
| 172 | } |
---|
| 173 | } |
---|
| 174 | return res; |
---|
| 175 | } |
---|
| 176 | |
---|
| 177 | const std::vector<G4String>& |
---|
| 178 | G4PhysListFactory::AvailablePhysLists() const |
---|
| 179 | { |
---|
| 180 | return listnames; |
---|
| 181 | } |
---|
| 182 | |
---|