Skip to content

Commit 63f91c4

Browse files
authored
Merge pull request #49185 from SeverinDiederichs/adept_integration
integrate AdePT via new physics list FTFP_BERT_EMA
2 parents e00e418 + 1f0701d commit 63f91c4

File tree

7 files changed

+290
-0
lines changed

7 files changed

+290
-0
lines changed

SimG4Core/PhysicsLists/BuildFile.xml

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -6,6 +6,7 @@
66
<use name="FWCore/MessageLogger"/>
77
<use name="SimG4Core/MagneticField"/>
88
<use name="SimG4Core/Physics"/>
9+
<use name="SimG4Core/Notification"/>
910
<export>
1011
<lib name="1"/>
1112
</export>

SimG4Core/PhysicsLists/plugins/BuildFile.xml

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -9,3 +9,10 @@
99
<use name="SimG4Core/PhysicsLists"/>
1010
<flags EDM_PLUGIN="1"/>
1111
</library>
12+
<iftool name="adept">
13+
<library file="adept/*.cc" name="SimG4CorePhysicsListsAdeptPlugins">
14+
<use name="SimG4Core/PhysicsLists"/>
15+
<use name="adept"/>
16+
<flags EDM_PLUGIN="1"/>
17+
</library>
18+
</iftool>
Lines changed: 163 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,163 @@
1+
#include "SimG4Core/PhysicsLists/plugins/adept/CMSEmStandardPhysicsA.h"
2+
#include "SimG4Core/Physics/interface/CMSG4TrackInterface.h"
3+
#include "SimG4Core/Notification/interface/CurrentG4Track.h"
4+
#include "FWCore/MessageLogger/interface/MessageLogger.h"
5+
6+
#include <CLHEP/Units/SystemOfUnits.h>
7+
#include "G4ParticleDefinition.hh"
8+
#include "G4LossTableManager.hh"
9+
#include "G4EmParameters.hh"
10+
#include "G4EmBuilder.hh"
11+
12+
#include "G4ComptonScattering.hh"
13+
#include "G4GammaConversion.hh"
14+
#include "G4PhotoElectricEffect.hh"
15+
16+
#include "G4MscStepLimitType.hh"
17+
18+
#include "G4eMultipleScattering.hh"
19+
#include "G4hMultipleScattering.hh"
20+
#include "G4eCoulombScatteringModel.hh"
21+
#include "G4CoulombScattering.hh"
22+
#include "G4WentzelVIModel.hh"
23+
#include "G4UrbanMscModel.hh"
24+
25+
#include "G4eIonisation.hh"
26+
#include "G4eBremsstrahlung.hh"
27+
#include "G4eplusAnnihilation.hh"
28+
29+
#include "G4hIonisation.hh"
30+
#include "G4ionIonisation.hh"
31+
32+
#include "G4ParticleTable.hh"
33+
#include "G4Gamma.hh"
34+
#include "G4Electron.hh"
35+
#include "G4Positron.hh"
36+
#include "G4GenericIon.hh"
37+
38+
#include "G4PhysicsListHelper.hh"
39+
#include "G4BuilderType.hh"
40+
#include "G4GammaGeneralProcess.hh"
41+
42+
#include "G4ProcessManager.hh"
43+
#include "G4TransportationWithMsc.hh"
44+
45+
#include "G4RegionStore.hh"
46+
#include "G4Region.hh"
47+
48+
#include <AdePT/core/AdePTConfiguration.hh>
49+
#include <AdePT/integration/AdePTTrackingManager.hh>
50+
51+
#include "G4HepEmConfig.hh"
52+
53+
CMSEmStandardPhysicsA::CMSEmStandardPhysicsA(G4int ver, const edm::ParameterSet& p)
54+
: G4VPhysicsConstructor("CMSEmStandard_ema") {
55+
fAdePTConfiguration = new AdePTConfiguration();
56+
SetVerboseLevel(ver);
57+
// EM parameters specific for this EM physics configuration
58+
G4EmParameters* param = G4EmParameters::Instance();
59+
param->SetDefaults();
60+
param->SetVerbose(ver);
61+
param->SetApplyCuts(true);
62+
param->SetStepFunction(0.8, 1 * CLHEP::mm);
63+
param->SetMscRangeFactor(0.2);
64+
param->SetMscStepLimitType(fMinimal);
65+
param->SetFluo(false);
66+
SetPhysicsType(bElectromagnetic);
67+
fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
68+
fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
69+
fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
70+
fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
71+
std::string msc = p.getParameter<std::string>("G4MscStepLimit");
72+
fStepLimitType = fUseSafety;
73+
if (msc == "UseSafetyPlus") {
74+
fStepLimitType = fUseSafetyPlus;
75+
}
76+
if (msc == "Minimal") {
77+
fStepLimitType = fMinimal;
78+
}
79+
double tcut = p.getParameter<double>("G4TrackingCut") * CLHEP::MeV;
80+
param->SetLowestElectronEnergy(tcut);
81+
param->SetLowestMuHadEnergy(tcut);
82+
}
83+
84+
void CMSEmStandardPhysicsA::ConstructParticle() {
85+
// minimal set of particles for EM physics
86+
G4EmBuilder::ConstructMinimalEmSet();
87+
}
88+
89+
void CMSEmStandardPhysicsA::ConstructProcess() {
90+
if (verboseLevel > 0) {
91+
int id = CMSG4TrackInterface::instance()->getThreadID();
92+
edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct EM Processes; EMA threadID=" << id;
93+
}
94+
95+
// This EM builder takes default models of Geant4 10 EMV.
96+
// Multiple scattering by WentzelVI for all particles except:
97+
// a) e+e- below 100 MeV for which the Urban model is used
98+
// b) ions for which Urban model is used
99+
G4EmBuilder::PrepareEMPhysics();
100+
101+
G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
102+
// processes used by several particles
103+
G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
104+
G4NuclearStopping* pnuc(nullptr);
105+
106+
const G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
107+
const G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
108+
109+
// G4HepEm is Active
110+
if (verboseLevel > 0) {
111+
edm::LogVerbatim("PhysicsList") << "AdePT is active, registering AdePTTrackingManager";
112+
}
113+
114+
// number of worker threads must be passed to AdePT
115+
fAdePTConfiguration->SetNumThreads(CurrentG4Track::numberOfWorkers());
116+
117+
// Construct the AdePT tracking manager
118+
auto* hepEmTM = new AdePTTrackingManager(fAdePTConfiguration, verboseLevel);
119+
120+
// now configure the G4HepEm config in the AdePT TM, as it will be used to define the physics
121+
G4HepEmConfig* config = hepEmTM->GetG4HepEmConfig();
122+
// First set global configuration parameters:
123+
// ------------------------------------------
124+
// The default MSC `RangeFactor`, `SafetyFactor`, `StepLimitType` parameters
125+
// as well as `SetApplyCuts` and `SetLowestElectronEnergy` are taken from
126+
// `G4EmParameters` so we set only the step function parameters here.
127+
128+
config->SetEnergyLossStepLimitFunctionParameters(0.8, 1.0 * CLHEP::mm);
129+
130+
// Then set special configuration for some regions:
131+
// ------------------------------------------------
132+
if (nullptr != aRegion) {
133+
// HCal region
134+
const G4String& rname = aRegion->GetName();
135+
config->SetMinimalMSCStepLimit(fStepLimitType == fMinimal, rname);
136+
config->SetMSCRangeFactor(fRangeFactor, rname);
137+
config->SetMSCSafetyFactor(fSafetyFactor, rname);
138+
}
139+
140+
if (nullptr != bRegion) {
141+
// HGCal region
142+
const G4String& rname = bRegion->GetName();
143+
config->SetMinimalMSCStepLimit(fStepLimitType == fMinimal, rname);
144+
config->SetMSCRangeFactor(fRangeFactor, rname);
145+
config->SetMSCSafetyFactor(fSafetyFactor, rname);
146+
147+
config->SetWoodcockTrackingRegion(rname);
148+
config->SetWDTEnergyLimit(0.5 * CLHEP::MeV);
149+
}
150+
151+
G4Electron::Electron()->SetTrackingManager(hepEmTM);
152+
G4Positron::Positron()->SetTrackingManager(hepEmTM);
153+
G4Gamma::Gamma()->SetTrackingManager(hepEmTM);
154+
155+
// generic ion
156+
G4ParticleDefinition* particle = G4GenericIon::GenericIon();
157+
G4ionIonisation* ionIoni = new G4ionIonisation();
158+
ph->RegisterProcess(hmsc, particle);
159+
ph->RegisterProcess(ionIoni, particle);
160+
161+
// muons, hadrons ions
162+
G4EmBuilder::ConstructCharged(hmsc, pnuc);
163+
}
Lines changed: 36 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,36 @@
1+
//--------------------------------------------------------------------
2+
//
3+
// 18.10.2025 S.Diederichs EM physics using AdePT based on
4+
// CMSEmStandardPhysicsof
5+
//
6+
//--------------------------------------------------------------------
7+
8+
#ifndef SimG4Core_PhysicsLists_CMSEmStandardPhysicsA_h
9+
#define SimG4Core_PhysicsLists_CMSEmStandardPhysicsA_h
10+
11+
#include "G4VPhysicsConstructor.hh"
12+
#include "globals.hh"
13+
#include "G4MscStepLimitType.hh"
14+
15+
#include "FWCore/ParameterSet/interface/ParameterSet.h"
16+
17+
class AdePTConfiguration;
18+
19+
class CMSEmStandardPhysicsA : public G4VPhysicsConstructor {
20+
public:
21+
CMSEmStandardPhysicsA(G4int ver, const edm::ParameterSet& p);
22+
~CMSEmStandardPhysicsA() override = default;
23+
24+
void ConstructParticle() override;
25+
void ConstructProcess() override;
26+
27+
private:
28+
G4double fRangeFactor;
29+
G4double fGeomFactor;
30+
G4double fSafetyFactor;
31+
G4double fLambdaLimit;
32+
G4MscStepLimitType fStepLimitType;
33+
AdePTConfiguration* fAdePTConfiguration;
34+
};
35+
36+
#endif
Lines changed: 65 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,65 @@
1+
#include "FTFPCMS_BERT_EMA.h"
2+
#include "SimG4Core/PhysicsLists/plugins/adept/CMSEmStandardPhysicsA.h"
3+
#include "SimG4Core/PhysicsLists/interface/CMSHadronPhysicsFTFP_BERT.h"
4+
#include "FWCore/MessageLogger/interface/MessageLogger.h"
5+
6+
#include "G4DecayPhysics.hh"
7+
#include "G4EmExtraPhysics.hh"
8+
#include "G4IonPhysics.hh"
9+
#include "G4StoppingPhysics.hh"
10+
#include "G4HadronElasticPhysics.hh"
11+
#include "G4HadronicParameters.hh"
12+
13+
FTFPCMS_BERT_EMA::FTFPCMS_BERT_EMA(const edm::ParameterSet& p) : PhysicsList(p) {
14+
int ver = p.getUntrackedParameter<int>("Verbosity", 0);
15+
bool emPhys = p.getUntrackedParameter<bool>("EMPhysics", true);
16+
bool hadPhys = p.getUntrackedParameter<bool>("HadPhysics", true);
17+
double minFTFP = p.getParameter<double>("EminFTFP") * CLHEP::GeV;
18+
double maxBERT = p.getParameter<double>("EmaxBERT") * CLHEP::GeV;
19+
double maxBERTpi = p.getParameter<double>("EmaxBERTpi") * CLHEP::GeV;
20+
std::string type = p.getParameter<std::string>("type");
21+
edm::LogVerbatim("PhysicsList") << "CMS Physics List " << type << "\n Flags for EM Physics: " << emPhys
22+
<< "; Hadronic Physics: " << hadPhys << "\n Transition energy Bertini/FTFP from "
23+
<< minFTFP / CLHEP::GeV << " to " << maxBERT / CLHEP::GeV << "; for pions to "
24+
<< maxBERTpi / CLHEP::GeV << " GeV";
25+
26+
if (emPhys) {
27+
// EM Physics
28+
RegisterPhysics(new CMSEmStandardPhysicsA(ver, p));
29+
30+
// Synchroton Radiation & GN Physics
31+
G4EmExtraPhysics* gn = new G4EmExtraPhysics(ver);
32+
RegisterPhysics(gn);
33+
bool mu = p.getParameter<bool>("G4MuonPairProductionByMuon");
34+
gn->MuonToMuMu(mu);
35+
edm::LogVerbatim("PhysicsList") << " Muon pair production by muons: " << mu;
36+
}
37+
38+
// Decays
39+
this->RegisterPhysics(new G4DecayPhysics(ver));
40+
41+
if (hadPhys) {
42+
bool ngen = p.getParameter<bool>("G4NeutronGeneralProcess");
43+
bool bc = p.getParameter<bool>("G4BCHadronicProcess");
44+
bool hn = p.getParameter<bool>("G4LightHyperNucleiTracking");
45+
auto param = G4HadronicParameters::Instance();
46+
param->SetEnableNeutronGeneralProcess(ngen);
47+
param->SetEnableBCParticles(bc);
48+
param->SetEnableHyperNuclei(hn);
49+
edm::LogVerbatim("PhysicsList") << " Eneble neutron general process: " << ngen
50+
<< "\n Enable b- and c- hadron physics: " << bc
51+
<< "\n Enable light hyper-nuclei physics: " << hn;
52+
53+
// Hadron Elastic scattering
54+
RegisterPhysics(new G4HadronElasticPhysics(ver));
55+
56+
// Hadron Physics
57+
RegisterPhysics(new CMSHadronPhysicsFTFP_BERT(minFTFP, maxBERT, maxBERTpi, minFTFP, maxBERT));
58+
59+
// Stopping Physics
60+
RegisterPhysics(new G4StoppingPhysics(ver));
61+
62+
// Ion Physics
63+
RegisterPhysics(new G4IonPhysics(ver));
64+
}
65+
}
Lines changed: 12 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,12 @@
1+
#ifndef SimG4Core_PhysicsLists_FTFPCMS_BERT_EMA_H
2+
#define SimG4Core_PhysicsLists_FTFPCMS_BERT_EMA_H
3+
4+
#include "SimG4Core/Physics/interface/PhysicsList.h"
5+
#include "FWCore/ParameterSet/interface/ParameterSet.h"
6+
7+
class FTFPCMS_BERT_EMA : public PhysicsList {
8+
public:
9+
FTFPCMS_BERT_EMA(const edm::ParameterSet& p);
10+
};
11+
12+
#endif
Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,6 @@
1+
#include "SimG4Core/Physics/interface/PhysicsListFactory.h"
2+
3+
#include "FTFPCMS_BERT_EMA.h"
4+
5+
typedef FTFPCMS_BERT_EMA FTFP_BERT_EMA;
6+
DEFINE_PHYSICSLIST(FTFP_BERT_EMA);

0 commit comments

Comments
 (0)