LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
larg4::OpParamSD Class Reference

#include "OpParamSD.h"

Inheritance diagram for larg4::OpParamSD:

Public Member Functions

 OpParamSD (G4String name, std::string ModelName, int Orientation, std::vector< std::vector< double > > Parameters)
 
virtual ~OpParamSD ()
 
virtual void Initialize (G4HCofThisEvent *)
 
virtual void EndOfEvent (G4HCofThisEvent *)
 
virtual void clear ()
 
virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)
 
virtual void DrawAll ()
 
virtual void PrintAll ()
 

Private Member Functions

G4bool G4BooleanRand (const G4double prob) const
 

Private Attributes

OpParamActionfOpa
 
std::map< G4int, bool > fPhotonAlreadyCrossed
 

Detailed Description

Definition at line 58 of file OpParamSD.h.

Constructor & Destructor Documentation

larg4::OpParamSD::OpParamSD ( G4String  name,
std::string  ModelName,
int  Orientation,
std::vector< std::vector< double > >  Parameters 
)

Definition at line 33 of file OpParamSD.cxx.

References fOpa.

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")
43  fOpa = new 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  }
OpParamAction * fOpa
Definition: OpParamSD.h:86
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
virtual larg4::OpParamSD::~OpParamSD ( )
inlinevirtual

Definition at line 64 of file OpParamSD.h.

References Initialize().

64 {}

Member Function Documentation

virtual void larg4::OpParamSD::clear ( void  )
inlinevirtual

Definition at line 72 of file OpParamSD.h.

72 {}
virtual void larg4::OpParamSD::DrawAll ( )
inlinevirtual

Definition at line 79 of file OpParamSD.h.

79 {}
virtual void larg4::OpParamSD::EndOfEvent ( G4HCofThisEvent )
inlinevirtual

Definition at line 69 of file OpParamSD.h.

69 {}
G4bool larg4::OpParamSD::G4BooleanRand ( const G4double  prob) const
inlineprivate

Definition at line 94 of file OpParamSD.h.

Referenced by ProcessHits().

95  {
96  /* Returns a random boolean variable with the specified probability */
97  return (G4UniformRand() < prob);
98  }
void larg4::OpParamSD::Initialize ( G4HCofThisEvent )
virtual

Definition at line 88 of file OpParamSD.cxx.

References fPhotonAlreadyCrossed.

89  {
90  fPhotonAlreadyCrossed.clear();
91  }
std::map< G4int, bool > fPhotonAlreadyCrossed
Definition: OpParamSD.h:87
virtual void larg4::OpParamSD::PrintAll ( )
inlinevirtual

Definition at line 80 of file OpParamSD.h.

80 {}
G4bool larg4::OpParamSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
virtual

Definition at line 61 of file OpParamSD.cxx.

References fOpa, fPhotonAlreadyCrossed, G4BooleanRand(), and larg4::OpParamAction::GetAttenuationFraction().

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  }
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)
OpParamAction * fOpa
Definition: OpParamSD.h:86

Member Data Documentation

OpParamAction* larg4::OpParamSD::fOpa
private

Definition at line 86 of file OpParamSD.h.

Referenced by OpParamSD(), and ProcessHits().

std::map<G4int, bool> larg4::OpParamSD::fPhotonAlreadyCrossed
private

Definition at line 87 of file OpParamSD.h.

Referenced by Initialize(), and ProcessHits().


The documentation for this class was generated from the following files: