LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
OpParamSD.cxx
Go to the documentation of this file.
1 //
6 // Implementation of the OpParamSD
7 //
8 // See comments in OpParamSD.h
9 //
10 // Ben Jones, MIT, 2013
11 //
12 
14 #include "cetlib_except/exception.h"
16 
17 #include "Geant4/G4DynamicParticle.hh"
18 #include "Geant4/G4SDManager.hh"
19 #include "Geant4/G4Step.hh"
20 #include "Geant4/G4StepPoint.hh"
21 #include "Geant4/G4ThreeVector.hh"
22 #include "Geant4/G4Track.hh"
23 #include "Geant4/G4TrackStatus.hh"
24 
25 namespace larg4 {
26 
27  OpParamSD::OpParamSD(G4String DetectorUniqueName,
28  std::string ModelName,
29  int Orientation,
30  std::vector<std::vector<double>> ModelParameters)
31  : G4VSensitiveDetector(DetectorUniqueName)
32  {
33  // Register self with sensitive detector manager
34  G4SDManager::GetSDMpointer()->AddNewDetector(this);
35 
36  if (ModelName == "OverlaidWireplanes")
37  fOpa = new OverlaidWireplanesAction(ModelParameters, Orientation);
38 
39  else if (ModelName == "TransparentPlaneAction")
41 
42  // else if(ModelName == "SimpleWireplane")
43  // fOpa = new SimpleWireplaneAction(ModelParameters, Orientation);
44 
45  // else if( your model here )
46 
47  else {
48  throw cet::exception("OpParamSD")
49  << "Error: Optical parameterization model " << ModelName << " not found.\n";
50  }
51  }
52 
53  //--------------------------------------------------------
54 
55  G4bool OpParamSD::ProcessHits(G4Step* aStep, G4TouchableHistory*)
56  {
57 
58  const G4Track* aTrack = aStep->GetTrack();
59  const G4DynamicParticle* aParticle = aTrack->GetDynamicParticle();
60 
61  G4ThreeVector mom = aParticle->GetMomentumDirection();
62  G4ThreeVector pos = aStep->GetPostStepPoint()->GetPosition();
63  if (!fPhotonAlreadyCrossed[aTrack->GetTrackID()]) {
64  if (G4BooleanRand(fOpa->GetAttenuationFraction(mom, pos))) {
65  // photon survives - let it carry on
66  fPhotonAlreadyCrossed[aTrack->GetTrackID()] = true;
67  }
68  else {
69  // photon is absorbed
70  aStep->GetTrack()->SetTrackStatus(fStopAndKill);
71  }
72  }
73  return true;
74  }
75 
76  //--------------------------------------------------------
77 
79  {
80  fPhotonAlreadyCrossed.clear();
81  }
82 
83 }
virtual void Initialize(G4HCofThisEvent *)
Definition: OpParamSD.cxx:78
Geant4 interface.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: OpParamSD.cxx:55
G4bool G4BooleanRand(const G4double prob) const
Definition: OpParamSD.h:78
std::map< G4int, bool > fPhotonAlreadyCrossed
Definition: OpParamSD.h:75
virtual double GetAttenuationFraction(G4ThreeVector PhotonDirection, G4ThreeVector PhotonPosition)
OpParamSD(G4String name, std::string ModelName, int Orientation, std::vector< std::vector< double >> Parameters)
Definition: OpParamSD.cxx:27
OpParamAction * fOpa
Definition: OpParamSD.h:74
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33