LArSoft  v07_13_02
Liquid Argon Software toolkit - http://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 
13 // LArSoft includes
14 
15 #include "larsim/LArG4/OpParamSD.h"
16 
17 #include "cetlib_except/exception.h"
18 
22 
23 // Geant4 includes
24 
25 #include "Geant4/G4SDManager.hh"
26 #include "Geant4/Randomize.hh"
27 #include "Geant4/G4RandomTools.hh"
28 
29 
30 namespace larg4{
31 
32 
33  OpParamSD::OpParamSD(G4String DetectorUniqueName, std::string ModelName, int Orientation, std::vector<std::vector<double> > ModelParameters)
34  : G4VSensitiveDetector(DetectorUniqueName)
35  {
36  // Register self with sensitive detector manager
37  G4SDManager::GetSDMpointer()->AddNewDetector(this);
38 
39  if(ModelName == "OverlaidWireplanes")
40  fOpa = new OverlaidWireplanesAction(ModelParameters, Orientation);
41 
42  else if(ModelName == "TransparentPlaneAction")
44 
45 // else if(ModelName == "SimpleWireplane")
46 // fOpa = new SimpleWireplaneAction(ModelParameters, Orientation);
47 
48 
49  // else if( your model here )
50 
51  else
52  {
53  throw cet::exception("OpParamSD")<<"Error: Optical parameterization model " << ModelName <<" not found.\n";
54  }
55 
56  }
57 
58 
59  //--------------------------------------------------------
60 
61  G4bool OpParamSD::ProcessHits(G4Step * aStep, G4TouchableHistory *)
62  {
63 
64  const G4Track* aTrack = aStep->GetTrack();
65  const G4DynamicParticle* aParticle = aTrack->GetDynamicParticle();
66 
67 
68  G4ThreeVector mom = aParticle->GetMomentumDirection();
69  G4ThreeVector pos = aStep->GetPostStepPoint()->GetPosition();
70  if(!fPhotonAlreadyCrossed[aTrack->GetTrackID()])
71  {
73  {
74  // photon survives - let it carry on
75  fPhotonAlreadyCrossed[aTrack->GetTrackID()]=true;
76  }
77  else
78  {
79  // photon is absorbed
80  aStep->GetTrack()->SetTrackStatus(fStopAndKill);
81  }
82  }
83  return true;
84  }
85 
86  //--------------------------------------------------------
87 
89  {
90  fPhotonAlreadyCrossed.clear();
91  }
92 
93 }
virtual void Initialize(G4HCofThisEvent *)
Definition: OpParamSD.cxx:88
Geant4 interface.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
contains objects relating to OpDet hits
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: OpParamSD.cxx:61
G4bool G4BooleanRand(const G4double prob) const
Definition: OpParamSD.h:94
std::map< G4int, bool > fPhotonAlreadyCrossed
Definition: OpParamSD.h:87
virtual double GetAttenuationFraction(G4ThreeVector PhotonDirection, G4ThreeVector PhotonPosition)
OpParamSD(G4String name, std::string ModelName, int Orientation, std::vector< std::vector< double > > Parameters)
Definition: OpParamSD.cxx:33
OpParamAction * fOpa
Definition: OpParamSD.h:86
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33