Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/gorad/src/GRPhysicsList.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //  Gorad (Geant4 Open-source Radiation Analysis and Design)
 27 //
 28 //  Author : Makoto Asai (SLAC National Accelerator Laboratory)
 29 //
 30 //  Development of Gorad is funded by NASA Johnson Space Center (JSC)
 31 //  under the contract NNJ15HK11B.
 32 //
 33 // ********************************************************************
 34 //
 35 // GRPhysicsList.cc
 36 //   Gorad Physics List
 37 //
 38 // History
 39 //   September 8th, 2020 : first implementation
 40 //
 41 // ********************************************************************
 42 
 43 #include "GRPhysicsList.hh"
 44 #include "GRPhysicsListMessenger.hh"
 45 #include "G4SystemOfUnits.hh"
 46 #include "G4PhysListFactory.hh"
 47 
 48 GRPhysicsList::GRPhysicsList()
 49 : PLName("FTFP_BERT"), physList(nullptr),
 50   EM_opt("Op_0"), Had_opt("FTFP_BERT"),
 51   addHP(false), addRDM(false), addRMC(false), addOptical(false),
 52   stepLimit_opt(-1)
 53 {
 54   factory = nullptr;
 55   messenger = new GRPhysicsListMessenger(this);
 56 
 57   globalCuts[0] = 0.7*mm; // e-
 58   globalCuts[1] = 0.7*mm; // e+
 59   globalCuts[2] = 0.7*mm; // gamma
 60   globalCuts[3] = 0.7*mm; // proton
 61 }
 62 
 63 GRPhysicsList::~GRPhysicsList()
 64 {
 65   delete factory;
 66   if(physList) delete physList;
 67   delete messenger;
 68 }
 69 
 70 void GRPhysicsList::ConstructParticle()
 71 {
 72   if(!physList) GeneratePL();
 73   physList->ConstructParticle();
 74 }
 75 
 76 void GRPhysicsList::ConstructProcess()
 77 {
 78   if(!physList) GeneratePL();
 79   physList->ConstructProcess();
 80 }
 81 
 82 #include "G4Region.hh"
 83 #include "G4ProductionCuts.hh"
 84 
 85 void GRPhysicsList::SetCuts()
 86 {
 87   if(!physList) GeneratePL();
 88   physList->SetCutValue(globalCuts[2],"gamma"); // gamma should be defined first!
 89   physList->SetCutValue(globalCuts[0],"e-");
 90   physList->SetCutValue(globalCuts[1],"e+");
 91   physList->SetCutValue(globalCuts[3],"proton");
 92 }
 93   
 94 void GRPhysicsList::SetGlobalCuts(G4double val)
 95 {
 96   for(G4int i=0; i<4; i++)
 97   { SetGlobalCut(i,val); }
 98   if(physList) SetCuts();
 99 }
100 
101 void GRPhysicsList::SetGlobalCut(G4int i, G4double val) 
102 {
103   globalCuts[i] = val; 
104   if(physList) SetCuts();
105 }
106 
107 void GRPhysicsList::GeneratePLName()
108 {
109   G4String plname = Had_opt;
110   if(addHP && Had_opt != "Shielding") plname += "_HP";
111 
112   G4String EMopt = "";
113   if(EM_opt=="Op_1") EMopt = "_EMV"; 
114   else if(EM_opt=="Op_3") EMopt = "_EMY"; 
115   else if(EM_opt=="Op_4") EMopt = "_EMZ"; 
116   else if(EM_opt=="LIV")  EMopt = "_LIV"; 
117   else if(EM_opt=="LIV_Pol") G4cout << "EM option <LIV_Pol> is under development." << G4endl;
118   plname += EMopt;
119 
120   auto valid = factory->IsReferencePhysList(plname);
121   if(valid) 
122   { PLName = plname; }
123   else
124   {
125     G4ExceptionDescription ed;
126     ed << "Physics List <" << plname << "> is not a valid reference physics list.";
127     G4Exception("GRPhysicsList::GeneratePLName()","GRPHYS0001",
128                 FatalException,ed);
129   }
130 }
131 
132 #include "G4RadioactiveDecayPhysics.hh"
133 #include "G4OpticalPhysics.hh"
134 #include "G4StepLimiterPhysics.hh"
135 #include "G4ParallelWorldPhysics.hh"
136 #include "G4GenericBiasingPhysics.hh"
137 #include "GRParallelWorldPhysics.hh"
138 #include "G4ProcessTable.hh"
139 #include "G4EmParameters.hh"
140 #include "G4HadronicParameters.hh"
141 
142 void GRPhysicsList::GeneratePL()
143 {
144   if(physList) return;
145 
146   factory = new G4PhysListFactory();
147   GeneratePLName();
148   physList = factory->GetReferencePhysList(PLName);
149   G4cout << "Creating " << PLName << " physics list ################ " << applyGeomImpBias << G4endl; 
150 
151   if(addRDM && Had_opt != "Shielding")
152   { physList->RegisterPhysics(new G4RadioactiveDecayPhysics()); 
153     G4cout << "Adding G4RadioactiveDecayPhysics ################ " << G4endl; }
154 
155   if(addRMC)
156   { G4cout << "Reverse Monte Calro option is under development." << G4endl; }
157 
158   if(stepLimit_opt>=0)
159   { physList->RegisterPhysics(new G4StepLimiterPhysics()); 
160     G4cout << "Adding G4StepLimiterPhysics ################ " << G4endl; }
161 
162   if(addOptical) // Optical processes should come last!
163   { physList->RegisterPhysics(new G4OpticalPhysics()); 
164     G4cout << "Adding G4OpticalPhysics ################ " << G4endl; }
165 
166   if(applyGeomImpBias) // Geometry Importance Biasing with parallel world
167   {
168     physList->RegisterPhysics(new GRParallelWorldPhysics("GeomBias",false));
169     G4cout << "Adding G4GenericBiasingPhysics for GeomBias ################ " << G4endl;
170   }
171 
172   G4int verbose = G4ProcessTable::GetProcessTable()->GetVerboseLevel();
173   physList->SetVerboseLevel(verbose);
174   G4EmParameters::Instance()->SetVerbose(verbose);
175   G4HadronicParameters::Instance()->SetVerboseLevel(verbose);
176 }
177 
178 #include "G4RegionStore.hh"
179 
180 G4Region* GRPhysicsList::FindRegion(const G4String& reg) const
181 {
182   auto store = G4RegionStore::GetInstance();
183   return store->GetRegion(reg);
184 }
185 
186 G4Region* GRPhysicsList::SetLocalCut(const G4String& reg,G4int i,G4double val)
187 {
188   auto regPtr = FindRegion(reg);
189   if(!regPtr) return regPtr;
190 
191   auto cuts = regPtr->GetProductionCuts();
192   if(!cuts) 
193   {
194     cuts = new G4ProductionCuts();
195     regPtr->SetProductionCuts(cuts);
196   }
197 
198   cuts->SetProductionCut(val,i);
199   return regPtr;
200 }
201 
202 G4double GRPhysicsList::GetLocalCut(const G4String& reg,G4int i) const
203 {
204   auto regPtr = FindRegion(reg);
205   G4double val = -1.0;
206   if(regPtr) 
207   { 
208     auto cuts = regPtr->GetProductionCuts();
209     if(cuts) val = cuts->GetProductionCut(i);
210   }
211   return val;
212 }
213 
214 #include "G4UserLimits.hh"
215 
216 G4Region* GRPhysicsList::SetLocalStepLimit(const G4String& reg,G4double val)
217 {
218   auto regPtr = FindRegion(reg);
219   if(!regPtr) return regPtr;
220 
221   auto uLim = regPtr->GetUserLimits();
222   if(!uLim)
223   {
224     uLim = new G4UserLimits(val);
225     regPtr->SetUserLimits(uLim);
226   }
227   else
228   { uLim->SetMaxAllowedStep(val); }
229   return regPtr;
230 }
231 
232 #include "G4Track.hh"
233 G4double GRPhysicsList::GetLocalStepLimit(const G4String& reg) const
234 {
235   static G4Track dummyTrack;
236   auto regPtr = FindRegion(reg);
237   G4double val = -1.0;
238   if(regPtr)
239   {
240     auto uLim = regPtr->GetUserLimits();
241     if(uLim) val = uLim->GetMaxAllowedStep(dummyTrack);
242   }
243   return val;
244 }
245 
246 void GRPhysicsList::SetGlobalStepLimit(G4double val)
247 { SetLocalStepLimit("DefaultRegionForTheWorld",val); }
248 
249 G4double GRPhysicsList::GetGlobalStepLimit() const
250 { return GetLocalStepLimit("DefaultRegionForTheWorld"); }
251 
252 
253