7 #include <G4SteppingVerbose.hh> 8 #include <G4SteppingManager.hh> 12 class Output :
public G4SteppingVerbose
20 void StepInfo() { G4SteppingVerbose::StepInfo();
22 void Reset() { trk.clear(); stp.clear(); vlm.clear(); pro.clear();
23 pdg.clear(); pid.clear(); xx.clear(); yy.clear(); zz.clear(); dt.clear();
24 de.clear(); dl.clear(); l.clear(); x.clear(); y.clear(); z.clear();
25 t.clear(); k.clear(); p.clear(); px.clear(); py.clear(); pz.clear();
26 q.clear(); et.clear(); }
57 #include <G4AnalysisManager.hh> 60 auto manager = G4AnalysisManager::Instance();
61 manager->CreateNtuple(
"t",
"Geant4 step points");
62 manager->CreateNtupleIColumn(
"n");
63 manager->CreateNtupleIColumn(
"m");
64 manager->CreateNtupleIColumn(
"trk",
trk);
65 manager->CreateNtupleIColumn(
"stp",
stp);
66 manager->CreateNtupleIColumn(
"vlm",
vlm);
67 manager->CreateNtupleIColumn(
"pro",
pro);
68 manager->CreateNtupleIColumn(
"pdg",
pdg);
69 manager->CreateNtupleIColumn(
"pid",
pid);
70 manager->CreateNtupleDColumn(
"xx",
xx);
71 manager->CreateNtupleDColumn(
"yy",
yy);
72 manager->CreateNtupleDColumn(
"zz",
zz);
73 manager->CreateNtupleDColumn(
"dt",
dt);
74 manager->CreateNtupleDColumn(
"de",
de);
75 manager->CreateNtupleDColumn(
"dl",
dl);
76 manager->CreateNtupleDColumn(
"l",
l);
77 manager->CreateNtupleDColumn(
"x",
x);
78 manager->CreateNtupleDColumn(
"y",
y);
79 manager->CreateNtupleDColumn(
"z",
z);
80 manager->CreateNtupleDColumn(
"t",
t);
81 manager->CreateNtupleDColumn(
"k",
k);
82 manager->CreateNtupleDColumn(
"p",
p);
83 manager->CreateNtupleDColumn(
"px",
px);
84 manager->CreateNtupleDColumn(
"py",
py);
85 manager->CreateNtupleDColumn(
"pz",
pz);
86 manager->CreateNtupleDColumn(
"q",
q);
87 manager->CreateNtupleDColumn(
"et",
et);
88 manager->FinishNtuple();
92 #include <G4NavigationHistory.hh> 98 G4TouchableHandle handle = fStep->GetPreStepPoint()->GetTouchableHandle();
99 int copyNo=handle->GetReplicaNumber();
100 if (copyNo<=0)
return;
101 if (
trk.size()>=10000) {
102 G4cout<<
"GEARS: # of step points >=10000. Recording stopped."<<G4endl;
103 fTrack->SetTrackStatus(fKillTrackAndSecondaries);
107 trk.push_back(fTrack->GetTrackID());
108 stp.push_back(fTrack->GetCurrentStepNumber());
109 vlm.push_back(copyNo);
110 pdg.push_back(fTrack->GetDefinition()->GetPDGEncoding());
111 pid.push_back(fTrack->GetParentID());
114 pro.push_back(fTrack->GetCreatorProcess()->GetProcessType()*1000
115 + fTrack->GetCreatorProcess()->GetProcessSubType());
116 else pro.push_back(1000);
118 const G4VProcess *pr = fStep->GetPostStepPoint()->GetProcessDefinedStep();
119 if (pr)
pro.push_back(pr->GetProcessType()*1000 + pr->GetProcessSubType());
120 else pro.push_back(900);
123 k.push_back(fTrack->GetKineticEnergy()/CLHEP::keV);
124 p.push_back(fTrack->GetMomentum().mag()/CLHEP::keV);
125 q.push_back(fStep->GetPostStepPoint()->GetCharge());
126 l.push_back(fTrack->GetTrackLength()/CLHEP::mm);
128 px.push_back(fTrack->GetMomentumDirection().x());
129 py.push_back(fTrack->GetMomentumDirection().y());
130 pz.push_back(fTrack->GetMomentumDirection().z());
131 de.push_back(fStep->GetTotalEnergyDeposit()/CLHEP::keV);
132 dl.push_back(fTrack->GetStepLength()/CLHEP::mm);
134 t.push_back(fTrack->GetGlobalTime()/CLHEP::ns);
135 x.push_back(fTrack->GetPosition().x()/CLHEP::mm);
136 y.push_back(fTrack->GetPosition().y()/CLHEP::mm);
137 z.push_back(fTrack->GetPosition().z()/CLHEP::mm);
139 G4ThreeVector pos = handle->GetHistory()->GetTopTransform()
140 .TransformPoint(fStep->GetPostStepPoint()->GetPosition());
141 xx.push_back(pos.x()/CLHEP::mm);
142 yy.push_back(pos.y()/CLHEP::mm);
143 zz.push_back(pos.z()/CLHEP::mm);
144 dt.push_back(fTrack->GetLocalTime()/CLHEP::ns);
146 if (
de.back()>0 && G4StrUtil::contains(handle->GetVolume()->GetName(),
"(S)")) {
147 if (
et.size()<(
unsigned int)copyNo+1)
et.resize((
unsigned int)copyNo+1);
148 et[copyNo]+=
de.back();
et[0]+=
de.back();
153 #include <G4OpticalSurface.hh> 169 #include <G4tgrLineProcessor.hh> 176 G4MaterialPropertiesTable* CreateMaterialPropertiesTable(
177 const vector<G4String> &words,
size_t idxOfWords);
198 G4bool ProcessLine(
const vector<G4String> &words);
203 #include <G4NistManager.hh> 204 #include <G4tgbMaterialMgr.hh> 205 #include <G4UImessenger.hh> 209 G4bool processed = G4tgrLineProcessor::ProcessLine(words);
210 if (processed)
return true;
213 G4String tag = words[0]; G4StrUtil::to_lower(tag);
214 if (tag.substr(0,5)==
":prop") {
215 G4NistManager *mgr = G4NistManager::Instance(); mgr->SetVerbose(2);
216 G4Material *m = mgr->FindOrBuildMaterial(words[1]);
218 m=G4tgbMaterialMgr::GetInstance()->FindOrBuildG4Material(words[1]);
219 G4cout<<
"GEARS: Set optical properties of "<<words[1]<<
":"<<G4endl;
220 m->SetMaterialPropertiesTable(CreateMaterialPropertiesTable(words,2));
222 }
else if (tag.substr(0,5)==
":surf") {
229 bdr->optic =
new G4OpticalSurface(words[1]);
232 while (i<words.size()) {
233 G4String setting = words[i], value = words[i+1];
234 G4StrUtil::to_lower(setting); G4StrUtil::to_lower(value);
235 if (setting==
"property") {
237 }
else if (setting==
"type") {
238 if (value==
"dielectric_metal") bdr->optic->SetType(dielectric_metal);
239 else if (value==
"dielectric_dielectric")
240 bdr->optic->SetType(dielectric_dielectric);
241 else if (value==
"firsov") bdr->optic->SetType(firsov);
242 else if (value==
"x_ray") bdr->optic->SetType(x_ray);
243 else G4cout<<
"GERAS: Ignore unknown surface type "<<value<<G4endl;
244 }
else if (setting==
"model") {
245 if (value==
"glisur") bdr->optic->SetModel(glisur);
246 else if (value==
"unified") bdr->optic->SetModel(unified);
247 else G4cout<<
"GERAS: Ignore unknown surface model "<<value<<G4endl;
248 }
else if (setting==
"finish") {
249 G4cout<<
"GERAS: Set surface finish to be "<<value<<G4endl;
250 if (value==
"polished") bdr->optic->SetFinish(polished);
251 else if (value==
"polishedfrontpainted")
252 bdr->optic->SetFinish(polishedfrontpainted);
253 else if (value==
"polishedbackpainted")
254 bdr->optic->SetFinish(polishedbackpainted);
255 else if (value==
"ground") bdr->optic->SetFinish(ground);
256 else if (value==
"groundfrontpainted")
257 bdr->optic->SetFinish(groundfrontpainted);
258 else if (value==
"groundbackpainted")
259 bdr->optic->SetFinish(groundbackpainted);
261 G4cout<<
"GERAS: Unknown surface finish: "<<value<<
"!"<<G4endl;
264 }
else if (setting==
"sigma_alpha") {
265 bdr->optic->SetSigmaAlpha(G4UIcommand::ConvertToInt(value));
267 G4cout<<
"GERAS: Ignore unknown surface setting "<<value<<G4endl;
270 if (i<words.size()) {
271 G4cout<<
"GEARS: Set optical properties of "<<bdr->name<<
":"<<G4endl;
272 bdr->optic->SetMaterialPropertiesTable(
273 CreateMaterialPropertiesTable(words,i));
281 #include <G4tgrUtils.hh> 282 G4MaterialPropertiesTable* LineProcessor::CreateMaterialPropertiesTable(
283 const vector<G4String> &words,
size_t idxOfWords)
285 bool photonEnergyUnDefined=
true;
287 double *energies=NULL;
288 auto table =
new G4MaterialPropertiesTable();
289 for (
size_t i=idxOfWords; i<words.size(); i++) {
290 G4String
property = words[i]; G4StrUtil::to_upper(property);
291 if (G4StrUtil::contains(property,
"TIMECONSTANT") ||
292 G4StrUtil::contains(property,
"SCINTILLATIONYIELD") ||
293 property==
"RESOLUTIONSCALE" || property==
"YIELDRATIO") {
294 table->AddConstProperty(property, G4tgrUtils::GetDouble(words[i+1]));
295 G4cout<<
"GEARS: "<<
property<<
"="<<words[i+1]<<G4endl;
297 }
else if (property.substr(0,12)==
"PHOTON_ENERG") {
298 photonEnergyUnDefined=
false;
299 cnt = G4UIcommand::ConvertToInt(words[i+1]);
300 energies =
new double[cnt];
301 for (
int j=0; j<cnt; j++)
302 energies[j]=G4tgrUtils::GetDouble(words[i+2+j]);
305 if (photonEnergyUnDefined) {
306 G4cout<<
"GEARS: photon energies undefined, " 307 <<
"ignore all wavelength-dependent properties!"<<G4endl;
310 double *values =
new double[cnt];
311 for (
int j=0; j<cnt; j++) values[j]=G4tgrUtils::GetDouble(words[i+1+j]);
312 G4cout<<
"GEARS: "<<
property<<
"="<<values[0]<<
", " 313 <<values[1]<<
"..."<<G4endl;
314 table->AddProperty(property, energies, values, cnt);
324 #include <G4tgbDetectorBuilder.hh> 332 G4tgbDetectorBuilder() { fLineProcessor =
new LineProcessor(); }
334 const G4tgrVolume* ReadDetector();
338 G4VPhysicalVolume* ConstructDetector(
const G4tgrVolume* topVol);
345 #include <G4tgrVolumeMgr.hh> 346 #include <G4tgrFileReader.hh> 347 #include <G4tgbVolumeMgr.hh> 350 G4tgrFileReader* reader = G4tgrFileReader::GetInstance();
351 reader->SetLineProcessor(fLineProcessor);
353 G4tgrVolumeMgr* mgr = G4tgrVolumeMgr::GetInstance();
354 const G4tgrVolume* world = mgr->GetTopVolume();
359 #include <G4LogicalBorderSurface.hh> 361 const G4tgrVolume* topVol)
363 G4VPhysicalVolume *world = G4tgbDetectorBuilder::ConstructDetector(topVol);
365 G4tgbVolumeMgr* tgbVolmgr = G4tgbVolumeMgr::GetInstance();
368 G4String physV1 = border->
v1.substr(0,border->
v1.find(
":"));
369 G4String physV2 = border->
v2.substr(0,border->
v2.find(
":"));
370 int copyNo1 = atoi(border->
v1.substr(border->
v1.find(
":")+1).data());
371 int copyNo2 = atoi(border->
v2.substr(border->
v2.find(
":")+1).data());
372 G4LogicalVolume *m1=tgbVolmgr->FindG4PhysVol(physV1)->GetMotherLogical();
373 G4LogicalVolume *m2=tgbVolmgr->FindG4PhysVol(physV2)->GetMotherLogical();
375 G4VPhysicalVolume *v1=0, *v2=0;
376 for (
int i=0; i<(int)m1->GetNoDaughters(); i++) {
377 v1 = m1->GetDaughter(i);
378 if (v1->GetCopyNo()==copyNo1)
break;
380 for (
int i=0; i<(int)m2->GetNoDaughters(); i++) {
381 v2 = m2->GetDaughter(i);
382 if (v2->GetCopyNo()==copyNo2)
break;
385 new G4LogicalBorderSurface(border->
name,v1,v2,border->
optic);
386 G4cout<<
"Border surface "<<border->
name<<
" in between " 387 <<physV1<<
":"<<copyNo1<<
" and "<<physV2<<
":"<<copyNo2
396 #include <G4UImessenger.hh> 397 #include <G4UIdirectory.hh> 398 #include <G4UIcmdWithAString.hh> 399 #include <G4UIcmdWith3VectorAndUnit.hh> 400 #include <G4VUserDetectorConstruction.hh> 413 class Detector :
public G4VUserDetectorConstruction,
public G4UImessenger
417 ~Detector() {
delete fCmdSetB;
delete fCmdSrc;
delete fCmdOut; }
418 G4VPhysicalVolume* Construct();
419 void SetNewValue(G4UIcommand* cmd, G4String value);
422 G4UIcmdWith3VectorAndUnit* fCmdSetB;
423 G4UIcmdWithAString* fCmdSrc;
424 G4UIcmdWithAString* fCmdOut;
425 G4VPhysicalVolume * fWorld;
432 fCmdOut =
new G4UIcmdWithAString(
"/geometry/export",
this);
433 fCmdOut->SetGuidance(
"Export geometry gdml file name");
434 fCmdOut->SetParameterName(
"gdml geometry output",
false);
435 fCmdOut->AvailableForStates(G4State_Idle);
440 fCmdSrc =
new G4UIcmdWithAString(
"/geometry/source",
this);
441 fCmdSrc->SetGuidance(
"Set geometry source file name");
442 fCmdSrc->SetParameterName(
"text geometry input",
false);
443 fCmdSrc->AvailableForStates(G4State_PreInit);
445 fCmdSetB =
new G4UIcmdWith3VectorAndUnit(
"/geometry/SetB",
this);
446 fCmdSetB->SetGuidance(
"Set uniform magnetic field value.");
447 fCmdSetB->SetParameterName(
"Bx",
"By",
"Bz",
false);
448 fCmdSetB->SetUnitCategory(
"Magnetic flux density");
452 #include <G4FieldManager.hh> 453 #include <G4UniformMagField.hh> 454 #include <G4TransportationManager.hh> 456 #include "G4GDMLParser.hh" 461 auto field =
new G4UniformMagField(0,0,0);
462 field->SetFieldValue(fCmdSetB->GetNew3VectorValue(value));
463 G4FieldManager* mgr =
464 G4TransportationManager::GetTransportationManager()->GetFieldManager();
465 mgr->SetDetectorField(field);
466 mgr->CreateChordFinder(field);
467 G4cout<<
"GEARS: Magnetic field is set to "<<value<<G4endl;
469 }
else if(cmd==fCmdOut) {
471 paser.Write(value,fWorld);
474 if (value.substr(value.length()-4)!=
"gdml") {
475 G4tgbVolumeMgr* mgr = G4tgbVolumeMgr::GetInstance();
476 mgr->AddTextFile(value);
478 mgr->SetDetectorBuilder(tgb);
479 fWorld = mgr->ReadAndConstructDetector();
484 fWorld=parser.GetWorldVolume();
492 #include "G4PVPlacement.hh" 496 G4cout<<
"GEARS: no detector specified, set to a 10x10x10 m^3 box."<<G4endl;
497 auto box =
new G4Box(
"hall", 5*CLHEP::m, 5*CLHEP::m, 5*CLHEP::m);
498 G4NistManager *nist = G4NistManager::Instance();
499 G4Material *vacuum = nist->FindOrBuildMaterial(
"G4_Galactic");
500 auto v =
new G4LogicalVolume(box, vacuum,
"hall");
501 fWorld =
new G4PVPlacement(0, G4ThreeVector(), v,
"hall", 0, 0, 0);
507 #include <G4VUserPrimaryGeneratorAction.hh> 508 #include <G4GeneralParticleSource.hh> 515 G4GeneralParticleSource* fSource;
518 : G4VUserPrimaryGeneratorAction(), fSource(0)
519 { fSource =
new G4GeneralParticleSource; }
522 { fSource->GeneratePrimaryVertex(evt); }
526 #include <G4RunManagerFactory.hh> 527 #include <G4UserRunAction.hh> 536 auto a = G4AnalysisManager::Instance();
if (a->GetFileName()==
"")
return;
538 Output* o = ((
Output*) G4VSteppingVerbose::GetInstance());
545 auto a = G4AnalysisManager::Instance();
546 if (a->GetFileName()!=
"") { a->Write(); a->CloseFile(); }
553 auto a = G4AnalysisManager::Instance();
554 Output* o = ((
Output*) G4VSteppingVerbose::GetInstance());
555 if (a->GetFileName()!=
"" && o->
stp.size()!=0) {
556 a->FillNtupleIColumn(0,o->
stp.size());
557 a->FillNtupleIColumn(1,o->
et.size()-1);
564 #include <G4UserEventAction.hh> 572 #include <G4UserStackingAction.hh> 573 #include <G4UIcmdWithADoubleAndUnit.hh> 582 G4UIcmdWithADoubleAndUnit *fCmdT;
585 fT0(0), fTimeWindow(0), fCmdT(0) {
586 fCmdT =
new G4UIcmdWithADoubleAndUnit(
"/grdm/setTimeWindow",
this);
587 fCmdT->SetGuidance(
"Time window to split a radioactive decay chain.");
588 fCmdT->SetGuidance(
"If a daughter nucleus appears after the window,");
589 fCmdT->SetGuidance(
"it is saved in a new entry in the output ntuple.");
590 fCmdT->SetGuidance(
"---Set it to <=0 to disable the splitting---");
591 fCmdT->SetParameterName(
"time window",
false,
true);
592 fCmdT->SetDefaultUnit(
"s");
593 fCmdT->AvailableForStates(G4State_PreInit, G4State_Idle);
597 if (fTimeWindow<=0)
return fUrgent;
598 if (trk->GetGlobalTime()>fT0+fTimeWindow)
return fWaiting;
602 if (fTimeWindow<=0)
return;
603 Output* o = ((
Output*) G4VSteppingVerbose::GetInstance());
608 {
if (cmd!=fCmdT)
return; fTimeWindow = fCmdT->GetNewDoubleValue(value); }
612 #include <G4VUserActionInitialization.hh> 613 class Action :
public G4VUserActionInitialization
624 #include <G4PhysListFactory.hh> 625 #include <G4ScoringManager.hh> 626 #include <G4VisExecutive.hh> 627 #include <G4UIExecutive.hh> 628 #include <G4UImanager.hh> 632 int main(
int argc,
char **argv)
635 G4VSteppingVerbose::SetInstance(
new Output);
636 auto run=G4RunManagerFactory::CreateRunManager(G4RunManagerType::SerialOnly);
637 G4PhysListFactory factory;
638 run->SetUserInitialization(factory.ReferencePhysList());
639 run->SetUserInitialization(
new Detector);
640 run->SetUserInitialization(
new Action);
641 G4ScoringManager::GetScoringManager();
642 G4UIExecutive* ui =
nullptr;
643 if (argc==1) { ui =
new G4UIExecutive(argc, argv); }
644 auto vis =
new G4VisExecutive(
"quiet");
650 G4String cmd =
"/control/execute ";
651 G4UImanager::GetUIpointer()->ApplyCommand(cmd+argv[1]);
653 delete vis;
delete run;
vector< double > l
length of track till this point [mm]
void EndOfEventAction(const G4Event *)
vector< double > p
momentum [keV]
vector< double > zz
z [mm] (origin: center of local volume)
void EndOfRunAction(const G4Run *)
Close output file.
vector< int > pid
parent particle's PDG encoding
BorderSurface * Border
pointer to a BorderSurface object
G4String v1
name of volume 1
vector< double > yy
y [mm] (origin: center of local volume)
G4bool ProcessLine(const vector< G4String > &words)
Overwrite G4tgrLineProcessor::ProcessLine to add new tags.
vector< int > pdg
PDG encoding.
vector< double > q
charge [elementary charge]
vector< double > px
x component of momentum direction
void Record()
Record simulated data.
G4OpticalSurface * optic
point to G4OpticalSurface object
vector< double > k
kinetic energy [keV]
vector< double > dt
time elapsed from previous step point [ns]
void BeginOfRunAction(const G4Run *)
enable output if output file name is not empty
vector< double > z
z [mm] (origin: center of the world)
Book keeping before and after a run.
void SetSteppingVerbose(int level)
vector< double > py
y component of momentum direction
Construct detector geometry.
vector< double > dl
step length [mm]
Dump simulation results to screen or a file.
vector< double > x
x [mm] (origin: center of the world)
vector< double > et
Total energy deposited in a volume [keV].
G4String v2
name of volume 2
G4VPhysicalVolume * ConstructDetector(const G4tgrVolume *topVol)
Construct detector based on text geometry description.
A link list of G4LogicalBorderSurface.
vector< double > y
y [mm] (origin: center of the world)
Book keeping before and after an event.
vector< double > de
energy deposited [keV]
const G4tgrVolume * ReadDetector()
Read text geometry input.
vector< double > xx
x [mm] (origin: center of local volume)
G4String name
name of the surface
Output()
use analysis manager to handle output
BorderSurface * next
link to next border surface
Split a radioactive decay chain to different events based on a time window.
vector< double > t
time elapsed from the beginning of an event [ns]
vector< int > vlm
volume copy number
Call Geant4 General Particle Source to generate particles.
void StepInfo()
Information of steps>0.
G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *trk)
send a daughter particle to waiting stack if it appears too late
vector< int > stp
step number
vector< int > trk
track ID
int main(int argc, char **argv)
The main function that calls individual components.
StackingAction()
created macro /grdm/setTimeWindow
vector< int > pro
process ID * 100 + sub-process ID
void NewStage()
save and reset output before processing waiting tracks
vector< double > pz
z component of momentum direction
void SetNewValue(G4UIcommand *cmd, G4String value)
virtual void GeneratePrimaries(G4Event *evt)
add sources to an event
void TrackingStarted()
Information of step 0 (initStep)
Extension to default text geometry file line processor.
Construct detector based on text geometry description.
G4VPhysicalVolume * Construct()
called at /run/initialize
void SetNewValue(G4UIcommand *cmd, G4String value)
for G4UI
void SaveAndResetEvent()
save and then reset an event