LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
larg4::OpDetReadoutGeometry Class Reference

#include "OpDetReadoutGeometry.h"

Inheritance diagram for larg4::OpDetReadoutGeometry:

Public Member Functions

 OpDetReadoutGeometry (G4String OpDetSensitiveName, const G4String name="OpDetReadoutGeometry")
 
virtual ~OpDetReadoutGeometry ()
 
virtual void Construct ()
 

Private Member Functions

void FindVolumes (G4VPhysicalVolume *, G4String, std::vector< G4Transform3D >, std::vector< G4LogicalVolume * > &, std::vector< G4Transform3D > &)
 

Private Attributes

std::vector< G4LogicalVolume * > fOpDetVolumes
 
std::vector< G4Transform3D > fOpDetTransformations
 
G4String fOpDetSensitiveName
 

Detailed Description

Definition at line 37 of file OpDetReadoutGeometry.h.

Constructor & Destructor Documentation

larg4::OpDetReadoutGeometry::OpDetReadoutGeometry ( G4String  OpDetSensitiveName,
const G4String  name = "OpDetReadoutGeometry" 
)

Definition at line 24 of file OpDetReadoutGeometry.cxx.

References fOpDetSensitiveName.

24  :
25  G4VUserParallelWorld(name)
26  {
27  fOpDetSensitiveName = OpDetSensitiveName;
28  }
larg4::OpDetReadoutGeometry::~OpDetReadoutGeometry ( )
virtual

Definition at line 30 of file OpDetReadoutGeometry.cxx.

31  {}

Member Function Documentation

void larg4::OpDetReadoutGeometry::Construct ( )
virtual

Definition at line 33 of file OpDetReadoutGeometry.cxx.

References larg4::OpDetLookup::AddPhysicalVolume(), FindVolumes(), fOpDetSensitiveName, fOpDetTransformations, fOpDetVolumes, g4b::DetectorConstruction::GetWorld(), larg4::OpDetLookup::Instance(), sim::LArG4Parameters::OpticalParamModels(), sim::LArG4Parameters::OpticalParamOrientations(), sim::LArG4Parameters::OpticalParamParameters(), sim::LArG4Parameters::OpticalParamVolumes(), and larg4::TheOpDetLookup.

34  {
35  mf::LogInfo("OpDetReadoutGeometry") << "constructing parallel world, looking for "
37 
38  // Get an empty parallel world volume
39  G4VPhysicalVolume * ParallelWorld = GetWorld();
40 
41  // Start with empty vectors
42  std::vector<G4LogicalVolume*> OpDetVolumes;
43  std::vector<G4Transform3D> OpDetTransformations;
44 
45  // Get the primary world volume
46  G4VPhysicalVolume * WorldPhysical = g4b::DetectorConstruction::GetWorld();
47 
48  // Find the OpDet volumes
49  std::vector<G4Transform3D> EmptyVector;
50  EmptyVector.clear();
51  FindVolumes(WorldPhysical, fOpDetSensitiveName, EmptyVector, OpDetVolumes, OpDetTransformations );
52 
53  fOpDetVolumes = OpDetVolumes;
54  fOpDetTransformations = OpDetTransformations;
55 
56  // Get the OpDet Lookup Table
57  OpDetLookup * TheOpDetLookup = OpDetLookup::Instance();
58 
59  // Create sensitive detector
60  OpDetSensitiveDetector * TheSD = new OpDetSensitiveDetector("OpDetSensitiveDetector");
61 
62 
63  if(OpDetVolumes.size()>0)
64  {
65 
66 
67 
68  // Make placements
69  for(unsigned int i=0; i!=OpDetVolumes.size(); i++)
70  {
71  std::stringstream VolumeName;
72  VolumeName.flush();
73  VolumeName.str("OpDetVolume_");
74  VolumeName<<i;
75 
76  G4Transform3D TheTransform = OpDetTransformations.at(i);
77 
78  G4VSolid * TheSolid = OpDetVolumes.at(i)->GetSolid();
79  G4Material * TheMaterial = OpDetVolumes.at(i)->GetMaterial();
80  G4LogicalVolume * TheLogVolume = new G4LogicalVolume(TheSolid,TheMaterial,VolumeName.str().c_str());
81 
82  TheLogVolume->SetSensitiveDetector(TheSD );
83 
84 
85  G4PVPlacement * ThePlacement
86  = new G4PVPlacement( TheTransform,
87  VolumeName.str().c_str(),
88  TheLogVolume,
89  ParallelWorld,
90  false,
91  0);
92 
93  // CLHEP::Hep3Vector trans = ThePlacement->GetTranslation();
94 
95  TheOpDetLookup->AddPhysicalVolume(ThePlacement);
96 
97  }
98  }
99 
100 
101  // Now add any optically parameterized volumes
102 
103  std::vector<G4LogicalVolume*> OpParamVolumesFound;
104  std::vector<G4Transform3D> OpParamTransformationsFound;
105 
107 
108  std::vector<std::string> OpParamModels = lgp->OpticalParamModels();
109  std::vector<std::string> OpParamVolumes = lgp->OpticalParamVolumes();
110  std::vector<int> OpParamOrientations = lgp->OpticalParamOrientations();;
111  std::vector<std::vector<std::vector<double > > > OpParamParameters = lgp->OpticalParamParameters();
112 
113  if((OpParamModels.size()!=OpParamVolumes.size())||
114  (OpParamModels.size()!=OpParamOrientations.size())||
115  (OpParamModels.size()!=OpParamParameters.size()))
116  {
117  throw cet::exception("OpDetReadoutGeometry")<<"sizes of OpParam specification vectors do not match\n";
118  }
119 
120  for(size_t imodel=0; imodel!=OpParamVolumes.size(); ++imodel)
121  {
122  EmptyVector.clear();
123  FindVolumes(WorldPhysical, OpParamVolumes.at(imodel), EmptyVector, OpParamVolumesFound, OpParamTransformationsFound );
124  mf::LogInfo("OpDetReadoutGeometry")<< "Found " << OpParamVolumesFound.size()<< " volumes of name " << OpParamVolumes.at(imodel)<< " to attach optical parameterization " << OpParamModels.at(imodel)<<std::endl;
125 
126 
127  // Since the same named model may be instantiated more than once with
128  // different parameters, create a unique sensitive detector name for this
129  // instance
130 
131  std::stringstream SDName("");
132  SDName<<OpParamModels.at(imodel)<<"_"<<imodel;
133 
134  OpParamSD * ParamSD = new OpParamSD(SDName.str().c_str(), OpParamModels.at(imodel), OpParamOrientations.at(imodel), OpParamParameters.at(imodel));
135 
136 
137  if(OpParamVolumesFound.size()>0)
138  {
139  for(unsigned int ivol=0; ivol!=OpParamVolumes.size(); ivol++)
140  {
141  std::stringstream VolumeName;
142  VolumeName.flush();
143  VolumeName.str("OpParamVolume_");
144  VolumeName<<ivol;
145 
146  G4Transform3D TheTransform = OpParamTransformationsFound.at(ivol);
147 
148  G4VSolid * TheSolid = OpParamVolumesFound.at(ivol)->GetSolid();
149  G4Material * TheMaterial = OpParamVolumesFound.at(ivol)->GetMaterial();
150  G4LogicalVolume * TheLogVolume = new G4LogicalVolume(TheSolid,TheMaterial,VolumeName.str().c_str());
151 
152 
153  G4PVPlacement * ThePlacement
154  = new G4PVPlacement( TheTransform,
155  VolumeName.str().c_str(),
156  TheLogVolume,
157  ParallelWorld,
158  false,
159  0);
160 
161  TheLogVolume->SetSensitiveDetector(ParamSD );
162 
163 
164  // This line just suppressed a compiler warning
165  ThePlacement->GetTranslation();
166 
167 
168  }
169 
170  }
171  }
172  }
void FindVolumes(G4VPhysicalVolume *, G4String, std::vector< G4Transform3D >, std::vector< G4LogicalVolume * > &, std::vector< G4Transform3D > &)
const std::vector< int > OpticalParamOrientations() const
const std::vector< std::string > OpticalParamModels() const
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
static G4VPhysicalVolume * GetWorld()
std::vector< G4Transform3D > fOpDetTransformations
std::vector< G4LogicalVolume * > fOpDetVolumes
const std::vector< std::vector< std::vector< double > > > OpticalParamParameters() const
static OpDetLookup * Instance()
Definition: OpDetLookup.cxx:31
OpDetLookup * TheOpDetLookup
Definition: OpDetLookup.cxx:22
const std::vector< std::string > OpticalParamVolumes() const
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
void larg4::OpDetReadoutGeometry::FindVolumes ( G4VPhysicalVolume *  PhysicalVolume,
G4String  OpDetName,
std::vector< G4Transform3D >  TransformSoFar,
std::vector< G4LogicalVolume * > &  OpDetVolumes,
std::vector< G4Transform3D > &  OpDetTransformations 
)
private

Definition at line 176 of file OpDetReadoutGeometry.cxx.

Referenced by Construct().

177  {
178 
179  // Add the next layer of transformation to the vector
180  G4ThreeVector Translation = PhysicalVolume->GetObjectTranslation();
181  G4RotationMatrix Rotation = PhysicalVolume->GetObjectRotationValue();
182  G4Transform3D NextTransform( Rotation, Translation );
183 
184  TransformSoFar.push_back(NextTransform);
185 
186 
187  // Check if this volume is a OpDet
188  G4String OpDetNameUnderscore = OpDetName+"_";
189  G4String VolumeName = PhysicalVolume->GetName();
190  if( ( VolumeName == OpDetName ) ||
191  ( VolumeName.find( OpDetNameUnderscore,0,OpDetNameUnderscore.length() )==0 )
192  )
193  {
194 
195  // We found a OpDet! Store its volume and global transformation
196  G4ThreeVector Trans(0,0,0);
197  G4RotationMatrix Rot(0,0,0);
198  G4Transform3D TotalTransform(Rot,Trans);
199  //for ( std::vector<G4Transform3D>::reverse_iterator it = TransformSoFar.rbegin();
200  // it != TransformSoFar.rend(); ++it )
201 
202 
203 
204  for ( std::vector<G4Transform3D>::iterator it = TransformSoFar.begin();
205  it != TransformSoFar.end(); ++it )
206  {
207  CLHEP::Hep3Vector trans = (*it).getTranslation();
208  CLHEP::HepRotation rot = (*it).getRotation();
209  TotalTransform = TotalTransform * (*it);
210  }
211 
212  OpDetVolumes.push_back(PhysicalVolume->GetLogicalVolume());
213  OpDetTransformations.push_back(TotalTransform);
214 
215  }
216  else
217  {
218  // We did not find a OpDet. Keep looking through daughters.
219  G4LogicalVolume * LogicalVolume = PhysicalVolume->GetLogicalVolume();
220 
221  // Loop through the daughters of the volume
222  G4int NumberDaughters = LogicalVolume->GetNoDaughters();
223  for(G4int i=0; i!=NumberDaughters; ++i)
224  {
225  // Get the ith daughter volume
226  G4VPhysicalVolume * Daughter = LogicalVolume->GetDaughter(i);
227 
228  // Recursively step into this volume
229  FindVolumes(Daughter, OpDetName, TransformSoFar, OpDetVolumes, OpDetTransformations);
230  }
231  }
232  }
void FindVolumes(G4VPhysicalVolume *, G4String, std::vector< G4Transform3D >, std::vector< G4LogicalVolume * > &, std::vector< G4Transform3D > &)
intermediate_table::iterator iterator

Member Data Documentation

G4String larg4::OpDetReadoutGeometry::fOpDetSensitiveName
private

Definition at line 48 of file OpDetReadoutGeometry.h.

Referenced by Construct(), and OpDetReadoutGeometry().

std::vector<G4Transform3D> larg4::OpDetReadoutGeometry::fOpDetTransformations
private

Definition at line 47 of file OpDetReadoutGeometry.h.

Referenced by Construct().

std::vector<G4LogicalVolume*> larg4::OpDetReadoutGeometry::fOpDetVolumes
private

Definition at line 46 of file OpDetReadoutGeometry.h.

Referenced by Construct().


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