LArSoft  v09_90_00
Liquid Argon Software toolkit - https://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 48 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 27 of file OpParamSD.cxx.

References fOpa.

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

Definition at line 55 of file OpParamSD.h.

References Initialize().

55 {}

Member Function Documentation

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

Definition at line 62 of file OpParamSD.h.

References ProcessHits().

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

Definition at line 68 of file OpParamSD.h.

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

Definition at line 59 of file OpParamSD.h.

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

Definition at line 78 of file OpParamSD.h.

Referenced by PrintAll(), and ProcessHits().

79  {
80  /* Returns a random boolean variable with the specified probability */
81  return (G4UniformRand() < prob);
82  }
void larg4::OpParamSD::Initialize ( G4HCofThisEvent )
virtual

Definition at line 78 of file OpParamSD.cxx.

References fPhotonAlreadyCrossed.

Referenced by ~OpParamSD().

79  {
80  fPhotonAlreadyCrossed.clear();
81  }
std::map< G4int, bool > fPhotonAlreadyCrossed
Definition: OpParamSD.h:75
virtual void larg4::OpParamSD::PrintAll ( )
inlinevirtual

Definition at line 69 of file OpParamSD.h.

References G4BooleanRand().

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

Definition at line 55 of file OpParamSD.cxx.

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

Referenced by clear().

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

Member Data Documentation

OpParamAction* larg4::OpParamSD::fOpa
private

Definition at line 74 of file OpParamSD.h.

Referenced by OpParamSD(), and ProcessHits().

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

Definition at line 75 of file OpParamSD.h.

Referenced by Initialize(), and ProcessHits().


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