GEARS
Geant4 Example Application with Rich features and Small footprints
gears.cc
Go to the documentation of this file.
1 
5 #include <vector>
6 using namespace std;
7 #include <G4SteppingVerbose.hh>
8 #include <G4SteppingManager.hh>
12 class Output : public G4SteppingVerbose
13 {
14  protected:
15  void Record();
16  public:
17  Output();
18  void TrackingStarted() { G4SteppingVerbose::TrackingStarted();
19  Record(); }
20  void StepInfo() { G4SteppingVerbose::StepInfo();
21  Record(); }
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(); }
27  void SetSteppingVerbose(int level) { fManager->SetVerboseLevel(level); }
28  int GetSteppingVerbose() { return fManager->GetverboseLevel(); }
29 
30  vector<int> trk;
31  vector<int> stp;
32  vector<int> vlm;
33  vector<int> pro;
34  vector<int> pdg;
35  vector<int> pid;
36  vector<double> xx;
37  vector<double> yy;
38  vector<double> zz;
39  vector<double> dt;
40  vector<double> de;
41  vector<double> dl;
42  vector<double> l;
43  vector<double> x;
44  vector<double> y;
45  vector<double> z;
46  vector<double> t;
47  vector<double> k;
48  vector<double> p;
49  vector<double> px;
50  vector<double> py;
51  vector<double> pz;
52  vector<double> q;
53  vector<double> et;
54 };
55 //______________________________________________________________________________
56 //
57 #include <G4AnalysisManager.hh>
58 Output::Output(): G4SteppingVerbose()
59 {
60  auto manager = G4AnalysisManager::Instance();
61  manager->CreateNtuple("t", "Geant4 step points");
62  manager->CreateNtupleIColumn("n"); // total number of recorded hits
63  manager->CreateNtupleIColumn("m"); // max copy number of sensitive volume
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();
89 }
90 //______________________________________________________________________________
91 //
92 #include <G4NavigationHistory.hh>
94 {
95  if (GetSilent()==1) // CopyState() won't be called in G4SteppingVerbose
96  CopyState(); // point fTrack, fStep, etc. to right places
97 
98  G4TouchableHandle handle = fStep->GetPreStepPoint()->GetTouchableHandle();
99  int copyNo=handle->GetReplicaNumber();
100  if (copyNo<=0) return; //skip uninteresting volumes (copy No. of world == 0)
101  if (trk.size()>=10000) {
102  G4cout<<"GEARS: # of step points >=10000. Recording stopped."<<G4endl;
103  fTrack->SetTrackStatus(fKillTrackAndSecondaries);
104  return;
105  }
106 
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());
112  if (stp.back()==0) { // step zero
113  if (pid.back()!=0) // not primary particle
114  pro.push_back(fTrack->GetCreatorProcess()->GetProcessType()*1000
115  + fTrack->GetCreatorProcess()->GetProcessSubType());
116  else pro.push_back(1000); // primary particle
117  } else {
118  const G4VProcess *pr = fStep->GetPostStepPoint()->GetProcessDefinedStep();
119  if (pr) pro.push_back(pr->GetProcessType()*1000 + pr->GetProcessSubType());
120  else pro.push_back(900); // not sure why pr can be zero
121  }
122 
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);
127 
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);
133 
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);
138 
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);
145 
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();
149  }
150 }
151 //______________________________________________________________________________
152 //
153 #include <G4OpticalSurface.hh>
160 {
161  G4String name;
162  G4String v1;
163  G4String v2;
164  G4OpticalSurface* optic;
166 };
167 //______________________________________________________________________________
168 //
169 #include <G4tgrLineProcessor.hh>
173 class LineProcessor: public G4tgrLineProcessor
174 {
175  private:
176  G4MaterialPropertiesTable* CreateMaterialPropertiesTable(
177  const vector<G4String> &words, size_t idxOfWords);
178  public:
179  LineProcessor(): G4tgrLineProcessor(), Border(0) {};
181  while (Border) { // deleting G4OpticalSurface is done in Geant4
182  BorderSurface *next=Border->next;
183  delete Border;
184  Border=next;
185  }
186  }
198  G4bool ProcessLine(const vector<G4String> &words);
200 };
201 //______________________________________________________________________________
202 //
203 #include <G4NistManager.hh>
204 #include <G4tgbMaterialMgr.hh>
205 #include <G4UImessenger.hh>
206 G4bool LineProcessor::ProcessLine(const vector<G4String> &words)
207 {
208  // process default text geometry tags
209  G4bool processed = G4tgrLineProcessor::ProcessLine(words);
210  if (processed) return true; // no new tags involved
211 
212  // process added tags: prop & surf
213  G4String tag = words[0]; G4StrUtil::to_lower(tag); // to lower cases
214  if (tag.substr(0,5)==":prop") { // set optical properties of a material
215  G4NistManager *mgr = G4NistManager::Instance(); mgr->SetVerbose(2);
216  G4Material *m = mgr->FindOrBuildMaterial(words[1]);
217  if (m==NULL) // if not in NIST, then build in tgb
218  m=G4tgbMaterialMgr::GetInstance()->FindOrBuildG4Material(words[1]);
219  G4cout<<"GEARS: Set optical properties of "<<words[1]<<":"<<G4endl;
220  m->SetMaterialPropertiesTable(CreateMaterialPropertiesTable(words,2));
221  return true;
222  } else if (tag.substr(0,5)==":surf") { // define an optical surface
223  auto bdr = new BorderSurface;
224  bdr->next=Border; // save current border pointer
225  Border=bdr; // overwrite current border pointer
226  bdr->name=words[1];
227  bdr->v1=words[2];
228  bdr->v2=words[3];
229  bdr->optic = new G4OpticalSurface(words[1]);
230  size_t i=4;
231  // loop over optical surface setup lines
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") {
236  i++; break;
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);
260  else {
261  G4cout<<"GERAS: Unknown surface finish: "<<value<<"!"<<G4endl;
262  abort();
263  }
264  } else if (setting=="sigma_alpha") {
265  bdr->optic->SetSigmaAlpha(G4UIcommand::ConvertToInt(value));
266  } else
267  G4cout<<"GERAS: Ignore unknown surface setting "<<value<<G4endl;
268  i+=2;
269  }
270  if (i<words.size()) { // break while loop because of "property"
271  G4cout<<"GEARS: Set optical properties of "<<bdr->name<<":"<<G4endl;
272  bdr->optic->SetMaterialPropertiesTable(
273  CreateMaterialPropertiesTable(words,i));
274  }
275  return true;
276  } else
277  return false;
278 }
279 //______________________________________________________________________________
280 //
281 #include <G4tgrUtils.hh>
282 G4MaterialPropertiesTable* LineProcessor::CreateMaterialPropertiesTable(
283  const vector<G4String> &words, size_t idxOfWords)
284 {
285  bool photonEnergyUnDefined=true;
286  int cnt=0; // number of photon energy values
287  double *energies=NULL; // photon energy values
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;
296  i++; // property value has been used
297  } else if (property.substr(0,12)=="PHOTON_ENERG") {
298  photonEnergyUnDefined=false;
299  cnt = G4UIcommand::ConvertToInt(words[i+1]); // get array size
300  energies = new double[cnt]; // create energy array
301  for (int j=0; j<cnt; j++)
302  energies[j]=G4tgrUtils::GetDouble(words[i+2+j]);
303  i=i+1+cnt; // array has been used
304  } else { // wavelength-dependent properties
305  if (photonEnergyUnDefined) {
306  G4cout<<"GEARS: photon energies undefined, "
307  <<"ignore all wavelength-dependent properties!"<<G4endl;
308  break;
309  }
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);
315  delete[] values;
316  i=i+cnt;
317  }
318  }
319  delete[] energies;
320  return table;
321 }
322 //______________________________________________________________________________
323 //
324 #include <G4tgbDetectorBuilder.hh>
328 class TextDetectorBuilder : public G4tgbDetectorBuilder
329 {
330  public :
332  G4tgbDetectorBuilder() { fLineProcessor = new LineProcessor(); }
333  ~TextDetectorBuilder() { delete fLineProcessor; }
334  const G4tgrVolume* ReadDetector();
335 
338  G4VPhysicalVolume* ConstructDetector(const G4tgrVolume* topVol);
339 
340  private :
341  LineProcessor* fLineProcessor;
342 };
343 //______________________________________________________________________________
344 //
345 #include <G4tgrVolumeMgr.hh>
346 #include <G4tgrFileReader.hh>
347 #include <G4tgbVolumeMgr.hh>
349 {
350  G4tgrFileReader* reader = G4tgrFileReader::GetInstance();
351  reader->SetLineProcessor(fLineProcessor);
352  reader->ReadFiles();
353  G4tgrVolumeMgr* mgr = G4tgrVolumeMgr::GetInstance();
354  const G4tgrVolume* world = mgr->GetTopVolume();
355  return world;
356 }
357 //______________________________________________________________________________
358 //
359 #include <G4LogicalBorderSurface.hh>
361  const G4tgrVolume* topVol)
362 {
363  G4VPhysicalVolume *world = G4tgbDetectorBuilder::ConstructDetector(topVol);
364 
365  G4tgbVolumeMgr* tgbVolmgr = G4tgbVolumeMgr::GetInstance();
366  BorderSurface* border = fLineProcessor->Border;
367  while (border) {
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();
374  // search for physics volumes on the sides of the border
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;
379  }
380  for (int i=0; i<(int)m2->GetNoDaughters(); i++) {
381  v2 = m2->GetDaughter(i);
382  if (v2->GetCopyNo()==copyNo2) break;
383  }
384  if (v1 && v2) {
385  new G4LogicalBorderSurface(border->name,v1,v2,border->optic);
386  G4cout<<"Border surface "<<border->name<<" in between "
387  <<physV1<<":"<<copyNo1<<" and "<<physV2<<":"<<copyNo2
388  <<" added"<<G4endl;
389  }
390  border=border->next;
391  }
392  return world;
393 }
394 //______________________________________________________________________________
395 //
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
414 {
415  public:
416  Detector();
417  ~Detector() { delete fCmdSetB; delete fCmdSrc; delete fCmdOut; }
418  G4VPhysicalVolume* Construct();
419  void SetNewValue(G4UIcommand* cmd, G4String value);
420 
421  private:
422  G4UIcmdWith3VectorAndUnit* fCmdSetB;
423  G4UIcmdWithAString *fCmdSrc;
424  G4UIcmdWithAString *fCmdOut;
425  G4UIcommand *fCmdLmt;
426  G4VPhysicalVolume *fWorld;
427 };
428 //______________________________________________________________________________
429 //
430 Detector::Detector(): G4VUserDetectorConstruction(), G4UImessenger(), fWorld(0)
431 {
432 #ifdef hasGDML
433  fCmdOut = new G4UIcmdWithAString("/geometry/export",this);
434  fCmdOut->SetGuidance("Export geometry gdml file name");
435  fCmdOut->SetParameterName("gdml geometry output",false);
436  fCmdOut->AvailableForStates(G4State_Idle);
437 #else
438  fCmdOut=0;
439 #endif
440 
441  fCmdSrc = new G4UIcmdWithAString("/geometry/source",this);
442  fCmdSrc->SetGuidance("Set geometry source file name");
443  fCmdSrc->SetParameterName("text geometry input",false);
444  fCmdSrc->AvailableForStates(G4State_PreInit);
445 
446  fCmdLmt = new G4UIcommand("/tracking/setStepLimit",this);
447  fCmdLmt->SetGuidance("set max step length for a logical volume");
448  fCmdLmt->AvailableForStates(G4State_Idle);
449  G4UIparameter *p0 = new G4UIparameter("logical volume name", 's', false);
450  G4UIparameter *p1 = new G4UIparameter("step length in mm", 'd', false);
451  fCmdLmt->SetParameter(p0); fCmdLmt->SetParameter(p1);
452 
453  fCmdSetB = new G4UIcmdWith3VectorAndUnit("/geometry/SetB",this);
454  fCmdSetB->SetGuidance("Set uniform magnetic field value.");
455  fCmdSetB->SetParameterName("Bx", "By", "Bz", false);
456  fCmdSetB->SetUnitCategory("Magnetic flux density");
457 }
458 //______________________________________________________________________________
459 //
460 #include <G4UserLimits.hh>
461 #include <G4FieldManager.hh>
462 #include <G4UniformMagField.hh>
463 #include <G4LogicalVolumeStore.hh>
464 #include <G4TransportationManager.hh>
465 #ifdef hasGDML
466 #include "G4GDMLParser.hh"
467 #endif
468 void Detector::SetNewValue(G4UIcommand* cmd, G4String value)
469 {
470  if (cmd==fCmdSetB) {
471  auto field = new G4UniformMagField(0,0,0);
472  field->SetFieldValue(fCmdSetB->GetNew3VectorValue(value));
473  G4FieldManager* mgr =
474  G4TransportationManager::GetTransportationManager()->GetFieldManager();
475  mgr->SetDetectorField(field);
476  mgr->CreateChordFinder(field);
477  G4cout<<"GEARS: Magnetic field is set to "<<value<<G4endl;
478 #ifdef hasGDML
479  } else if (cmd==fCmdOut) {
480  G4GDMLParser paser;
481  paser.Write(value, fWorld);
482 #endif
483  } else if (cmd==fCmdLmt) {
484  istringstream iss(value); G4String vlm, limit; iss>>vlm>>limit;
485  auto v = G4LogicalVolumeStore::GetInstance()->GetVolume(vlm);
486  if (v) {
487  v->SetUserLimits(new G4UserLimits(stof(limit)));
488  G4cout<<"GEARS: max step length in "<<vlm<<": "<<limit<<" mm"<<G4endl;
489  }
490  } else { // cmd==fCmdSrc
491  if (value.substr(value.length()-4)!="gdml") { // text geometry input
492  G4tgbVolumeMgr* mgr = G4tgbVolumeMgr::GetInstance();
493  mgr->AddTextFile(value);
494  auto tgb = new TextDetectorBuilder;
495  mgr->SetDetectorBuilder(tgb);
496  fWorld = mgr->ReadAndConstructDetector();
497 #ifdef hasGDML
498  } else { // GDML input
499  G4GDMLParser parser;
500  parser.Read(value);
501  fWorld=parser.GetWorldVolume();
502 #endif
503  }
504  }
505 }
506 //______________________________________________________________________________
507 //
508 #include "G4Box.hh"
509 #include "G4PVPlacement.hh"
510 G4VPhysicalVolume* Detector::Construct()
511 {
512  if (fWorld==NULL) {
513  G4cout<<"GEARS: no detector specified, set to a 10x10x10 m^3 box."<<G4endl;
514  auto box = new G4Box("hall", 5*CLHEP::m, 5*CLHEP::m, 5*CLHEP::m);
515  G4NistManager *nist = G4NistManager::Instance();
516  G4Material *vacuum = nist->FindOrBuildMaterial("G4_Galactic");
517  auto v = new G4LogicalVolume(box, vacuum, "hall");
518  fWorld = new G4PVPlacement(0, G4ThreeVector(), v, "hall", 0, 0, 0);
519  }
520  return fWorld;
521 }
522 //______________________________________________________________________________
523 //
524 #include <G4VUserPrimaryGeneratorAction.hh>
525 #include <G4GeneralParticleSource.hh>
529 class Generator : public G4VUserPrimaryGeneratorAction
530 {
531  private:
532  G4GeneralParticleSource* fSource;
533  public:
535  : G4VUserPrimaryGeneratorAction(), fSource(0)
536  { fSource = new G4GeneralParticleSource; }
537  virtual ~Generator() { delete fSource; }
538  virtual void GeneratePrimaries(G4Event* evt)
539  { fSource->GeneratePrimaryVertex(evt); }
540 };
541 //______________________________________________________________________________
542 //
543 #include <G4RunManagerFactory.hh>
544 #include <G4UserRunAction.hh>
545 #include <G4Run.hh>
549 class RunAction : public G4UserRunAction
550 {
551  public:
552  void BeginOfRunAction (const G4Run*) {
553  auto a = G4AnalysisManager::Instance(); if (a->GetFileName()=="") return;
554  a->OpenFile();
555  Output* o = ((Output*) G4VSteppingVerbose::GetInstance());
556  if (o->GetSteppingVerbose()==0) { // in case of /tracking/verbose 0
557  o->SetSilent(1); // avoid screen dump
558  o->SetSteppingVerbose(1);//enable calling StepInfo() in G4SteppingManager
559  }
560  }
561  void EndOfRunAction (const G4Run*) {
562  auto a = G4AnalysisManager::Instance();
563  if (a->GetFileName()!="") { a->Write(); a->CloseFile(); }
564  }
565 };
566 //______________________________________________________________________________
567 //
569 {
570  auto a = G4AnalysisManager::Instance();
571  Output* o = ((Output*) G4VSteppingVerbose::GetInstance());
572  if (a->GetFileName()!="" && o->stp.size()!=0) {
573  a->FillNtupleIColumn(0,o->stp.size());
574  a->FillNtupleIColumn(1,o->et.size()-1);
575  a->AddNtupleRow();
576  } // save n-tuple if it is not empty and output file name is specified
577  o->Reset(); // reset Output member variables for new record
578 }
579 //______________________________________________________________________________
580 //
581 #include <G4UserEventAction.hh>
585 class EventAction : public G4UserEventAction
586 { public: void EndOfEventAction(const G4Event*) { SaveAndResetEvent(); } };
587 //______________________________________________________________________________
588 //
589 #include <G4UserStackingAction.hh>
590 #include <G4UIcmdWithADoubleAndUnit.hh>
594 class StackingAction : public G4UserStackingAction, public G4UImessenger
595 {
596  private:
597  double fT0;
598  double fTimeWindow;
599  G4UIcmdWithADoubleAndUnit *fCmdT;
600  public:
601  StackingAction() : G4UserStackingAction(), G4UImessenger(),
602  fT0(0), fTimeWindow(0), fCmdT(0) {
603  fCmdT = new G4UIcmdWithADoubleAndUnit("/grdm/setTimeWindow", this);
604  fCmdT->SetGuidance("Time window to split a radioactive decay chain.");
605  fCmdT->SetGuidance("If a daughter nucleus appears after the window,");
606  fCmdT->SetGuidance("it is saved in a new entry in the output ntuple.");
607  fCmdT->SetGuidance("---Set it to <=0 to disable the splitting---");
608  fCmdT->SetParameterName("time window",false,true);
609  fCmdT->SetDefaultUnit("s");
610  fCmdT->AvailableForStates(G4State_PreInit, G4State_Idle);
611  }
612  ~StackingAction() { delete fCmdT; }
613  G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *trk) {
614  if (fTimeWindow<=0) return fUrgent; // no need to split
615  if (trk->GetGlobalTime()>fT0+fTimeWindow) return fWaiting; // split
616  else return fUrgent; // too fast to be split
617  }
618  void NewStage() { // called after processing urgent trk, before waiting trk
619  if (fTimeWindow<=0) return; // do nothing if no time window is specified
620  Output* o = ((Output*) G4VSteppingVerbose::GetInstance());
621  fT0 = o->t.back(); // update the reference time to the latest decay time
622  SaveAndResetEvent(); // end an event manually
623  }
624  void SetNewValue(G4UIcommand* cmd, G4String value)
625  { if (cmd!=fCmdT) return; fTimeWindow = fCmdT->GetNewDoubleValue(value); }
626 };
627 //______________________________________________________________________________
628 //
629 #include <G4VUserActionInitialization.hh>
630 class Action : public G4VUserActionInitialization
631 {
632  void Build() const {
633  SetUserAction(new RunAction);
634  SetUserAction(new Generator);
635  SetUserAction(new EventAction);
636  SetUserAction(new StackingAction);
637  }
638 };
639 //______________________________________________________________________________
640 //
641 #include <G4StepLimiterPhysics.hh>
642 #include <G4PhysListFactory.hh>
643 #include <G4ScoringManager.hh>
644 #include <G4VisExecutive.hh>
645 #include <G4UIExecutive.hh>
646 #include <G4UImanager.hh> // needed for g4.10 and above
650 int main(int argc, char **argv)
651 {
652  // inherit G4SteppingVerbose instead of G4UserSteppingAction to record data
653  G4VSteppingVerbose::SetInstance(new Output); // must be before run manager
654  auto run=G4RunManagerFactory::CreateRunManager(G4RunManagerType::SerialOnly);
655  G4PhysListFactory factory; auto physics=factory.ReferencePhysList();
656  physics->RegisterPhysics(new G4StepLimiterPhysics());
657  run->SetUserInitialization(physics); // initialize physics
658  run->SetUserInitialization(new Detector); // initialize detector
659  run->SetUserInitialization(new Action); // initialize user actions
660  G4ScoringManager::GetScoringManager(); // enable built-in scoring cmds
661  G4UIExecutive* ui = nullptr; // assume batch mode
662  if (argc==1) { ui = new G4UIExecutive(argc, argv); } // interactive mode
663  auto vis = new G4VisExecutive("quiet"); // visualization
664  vis->Initialize(); // do this after ui mode is decided
665  if (ui) { // interactive mode
666  ui->SessionStart(); // do this after vis
667  delete ui;
668  } else { // batch mode
669  G4String cmd = "/control/execute ";
670  G4UImanager::GetUIpointer()->ApplyCommand(cmd+argv[1]);
671  }
672  delete vis; delete run;
673  return 0;
674 }
675 // -*- C++; indent-tabs-mode:nil; tab-width:2 -*-
676 // vim: ft=cpp:ts=2:sts=2:sw=2:et
vector< double > l
length of track till this point [mm]
Definition: gears.cc:42
virtual ~Generator()
Definition: gears.cc:537
void EndOfEventAction(const G4Event *)
Definition: gears.cc:586
vector< double > p
momentum [keV]
Definition: gears.cc:48
vector< double > zz
z [mm] (origin: center of local volume)
Definition: gears.cc:38
void EndOfRunAction(const G4Run *)
Close output file.
Definition: gears.cc:561
vector< int > pid
parent particle&#39;s PDG encoding
Definition: gears.cc:35
BorderSurface * Border
pointer to a BorderSurface object
Definition: gears.cc:199
~LineProcessor()
Definition: gears.cc:180
G4String v1
name of volume 1
Definition: gears.cc:162
vector< double > yy
y [mm] (origin: center of local volume)
Definition: gears.cc:37
Generator()
Definition: gears.cc:534
G4bool ProcessLine(const vector< G4String > &words)
Overwrite G4tgrLineProcessor::ProcessLine to add new tags.
Definition: gears.cc:206
vector< int > pdg
PDG encoding.
Definition: gears.cc:34
LineProcessor()
Definition: gears.cc:179
vector< double > q
charge [elementary charge]
Definition: gears.cc:52
vector< double > px
x component of momentum direction
Definition: gears.cc:49
void Record()
Record simulated data.
Definition: gears.cc:93
G4OpticalSurface * optic
point to G4OpticalSurface object
Definition: gears.cc:164
vector< double > k
kinetic energy [keV]
Definition: gears.cc:47
vector< double > dt
time elapsed from previous step point [ns]
Definition: gears.cc:39
void BeginOfRunAction(const G4Run *)
enable output if output file name is not empty
Definition: gears.cc:552
Definition: gears.cc:630
~StackingAction()
Definition: gears.cc:612
vector< double > z
z [mm] (origin: center of the world)
Definition: gears.cc:45
Book keeping before and after a run.
Definition: gears.cc:549
void SetSteppingVerbose(int level)
Definition: gears.cc:27
vector< double > py
y component of momentum direction
Definition: gears.cc:50
Construct detector geometry.
Definition: gears.cc:413
vector< double > dl
step length [mm]
Definition: gears.cc:41
Dump simulation results to screen or a file.
Definition: gears.cc:12
vector< double > x
x [mm] (origin: center of the world)
Definition: gears.cc:43
vector< double > et
Total energy deposited in a volume [keV].
Definition: gears.cc:53
G4String v2
name of volume 2
Definition: gears.cc:163
void Reset()
Definition: gears.cc:22
G4VPhysicalVolume * ConstructDetector(const G4tgrVolume *topVol)
Construct detector based on text geometry description.
Definition: gears.cc:360
A link list of G4LogicalBorderSurface.
Definition: gears.cc:159
vector< double > y
y [mm] (origin: center of the world)
Definition: gears.cc:44
Book keeping before and after an event.
Definition: gears.cc:585
vector< double > de
energy deposited [keV]
Definition: gears.cc:40
const G4tgrVolume * ReadDetector()
Read text geometry input.
Definition: gears.cc:348
int GetSteppingVerbose()
Definition: gears.cc:28
vector< double > xx
x [mm] (origin: center of local volume)
Definition: gears.cc:36
G4String name
name of the surface
Definition: gears.cc:161
Output()
use analysis manager to handle output
Definition: gears.cc:58
BorderSurface * next
link to next border surface
Definition: gears.cc:165
Split a radioactive decay chain to different events based on a time window.
Definition: gears.cc:594
vector< double > t
time elapsed from the beginning of an event [ns]
Definition: gears.cc:46
vector< int > vlm
volume copy number
Definition: gears.cc:32
Call Geant4 General Particle Source to generate particles.
Definition: gears.cc:529
void StepInfo()
Information of steps>0.
Definition: gears.cc:20
~Detector()
Definition: gears.cc:417
G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *trk)
send a daughter particle to waiting stack if it appears too late
Definition: gears.cc:613
vector< int > stp
step number
Definition: gears.cc:31
vector< int > trk
track ID
Definition: gears.cc:30
int main(int argc, char **argv)
The main function that calls individual components.
Definition: gears.cc:650
StackingAction()
created macro /grdm/setTimeWindow
Definition: gears.cc:601
vector< int > pro
process ID * 100 + sub-process ID
Definition: gears.cc:33
void NewStage()
save and reset output before processing waiting tracks
Definition: gears.cc:618
vector< double > pz
z component of momentum direction
Definition: gears.cc:51
void SetNewValue(G4UIcommand *cmd, G4String value)
Definition: gears.cc:624
virtual void GeneratePrimaries(G4Event *evt)
add sources to an event
Definition: gears.cc:538
Detector()
Definition: gears.cc:430
void TrackingStarted()
Information of step 0 (initStep)
Definition: gears.cc:18
Extension to default text geometry file line processor.
Definition: gears.cc:173
Construct detector based on text geometry description.
Definition: gears.cc:328
G4VPhysicalVolume * Construct()
called at /run/initialize
Definition: gears.cc:510
void SetNewValue(G4UIcommand *cmd, G4String value)
for G4UI
Definition: gears.cc:468
void SaveAndResetEvent()
save and then reset an event
Definition: gears.cc:568