VGM Version 5.3
Loading...
Searching...
No Matches
Placement.cxx
Go to the documentation of this file.
1// $Id$
2
3// -----------------------------------------------------------------------
4// The Geant4GM package of the Virtual Geometry Model
5// Copyright (C) 2007, Ivana Hrivnacova
6// All rights reserved.
7//
8// For the licensing terms see vgm/LICENSE.
9// Contact: ivana@ipno.in2p3.fr
10// -----------------------------------------------------------------------
11
12//
13// Class Placement
14// -------------------
15// VGM implementation for Geant4 positions of volumes.
16//
17// Author: Ivana Hrivnacova; IPN Orsay
18
19#include "VGM/common/Axis.h"
20#include "VGM/solids/ICons.h"
23#include "VGM/solids/ITubs.h"
24#include "VGM/volumes/IVolume.h"
25
26#include "ClhepVGM/Units.h"
27#include "ClhepVGM/transform.h"
28
32
33#include "G4LogicalVolume.hh"
34#include "G4Material.hh"
35#include "G4PVDivision.hh"
36#include "G4PVParameterised.hh"
37#include "G4PVPlacement.hh"
38#include "G4PVReplica.hh"
39#include "G4ReflectionFactory.hh"
40#include "G4ReplicatedSlice.hh"
41#include "G4VDivisionParameterisation.hh"
42#include "G4VPhysicalVolume.hh"
43#include "G4VSolid.hh"
44
45//_____________________________________________________________________________
47 VGM::IVolume* volume, VGM::IVolume* motherVolume, G4VPhysicalVolume* pv)
48 : VGM::IPlacement(),
49 BaseVGM::VPlacement(volume, motherVolume),
50 fPhysicalVolume(pv)
51{
53
54 // Register physical volume in the map
55 Geant4GM::PlacementMap::Instance()->AddPlacement(this, fPhysicalVolume);
56}
57
58//_____________________________________________________________________________
59Geant4GM::Placement::Placement() : VGM::IPlacement(), BaseVGM::VPlacement()
60{
62}
63
64//_____________________________________________________________________________
66 : VGM::IPlacement(rhs), BaseVGM::VPlacement(rhs)
67{
69}
70
71//_____________________________________________________________________________
76
77//
78// static public functions
79//
80
81//_____________________________________________________________________________
83{
84 switch (axis) {
85 case VGM::kXAxis:
86 return kXAxis;
87 break;
88 case VGM::kYAxis:
89 return kYAxis;
90 break;
91 case VGM::kZAxis:
92 return kZAxis;
93 break;
94 case VGM::kRho:
95 return kRho;
96 break;
97 case VGM::kRadial3D:
98 return kRadial3D;
99 break;
100 case VGM::kPhi:
101 return kPhi;
102 break;
104 return kUndefined;
105 break;
106 default:
107 return kUndefined;
108 break;
109 }
110}
111
112//_____________________________________________________________________________
114{
115 switch (axis) {
116 case kXAxis:
117 return VGM::kXAxis;
118 break;
119 case kYAxis:
120 return VGM::kYAxis;
121 break;
122 case kZAxis:
123 return VGM::kZAxis;
124 break;
125 case kRho:
126 return VGM::kRho;
127 break;
128 case kRadial3D:
129 return VGM::kRadial3D;
130 break;
131 case kPhi:
132 return VGM::kPhi;
133 break;
134 case kUndefined:
135 return VGM::kUnknownAxis;
136 break;
137 default:
138 return VGM::kUnknownAxis;
139 break;
140 }
141}
142
143//
144// public functions
145//
146
147//_____________________________________________________________________________
149{
150 // paremeterised volumes with user defined parameterisation
151 if (dynamic_cast<G4PVParameterised*>(fPhysicalVolume) &&
152 (! dynamic_cast<G4PVDivision*>(fPhysicalVolume)))
153 return VGM::kParameterised;
154
155 // divisions and replicas
156 if (fPhysicalVolume->IsParameterised() || fPhysicalVolume->IsReplicated())
158
160}
161
162//_____________________________________________________________________________
163std::string Geant4GM::Placement::Name() const
164{
165 //
166 return fPhysicalVolume->GetName();
167}
168
169//_____________________________________________________________________________
171{
172 //
173 return fPhysicalVolume->GetCopyNo();
174}
175
176//_____________________________________________________________________________
178{
179 //
180 return ClhepVGM::Transform(*fPhysicalVolume->GetObjectRotation(),
181 fPhysicalVolume->GetObjectTranslation());
182}
183
184//_____________________________________________________________________________
186 double& width, double& offset, double& halfGap) const
187{
188 // Fill data only if multiple placement
189 if (Type() != VGM::kMultiplePlacement) return false;
190
191 // Different get functions for PVReplica and PVDivisions
192 EAxis g4Axis = kUndefined;
193 if (dynamic_cast<G4PVDivision*>(fPhysicalVolume) ||
194 dynamic_cast<G4ReplicatedSlice*>(fPhysicalVolume)) {
195 // G4PVDivision derives from G4PVReplica (since Geant4 11.0)
196 // - it must be tested first
197 G4VDivisionParameterisation* param =
198 dynamic_cast<G4VDivisionParameterisation*>(
199 fPhysicalVolume->GetParameterisation());
200
201 if (!param) {
202 std::cerr << " Geant4GM::Placement::MultiplePlacementData: "
203 << std::endl;
204 std::cerr << " Incorrect parameterisation type for "
205 "G4PVDivision/G4ReplicatedSlice"
206 << std::endl;
207 std::cerr << " (G4VDivisionParameterisation type was expected.)"
208 << std::endl;
209 std::cerr << "*** Error: Aborting execution ***" << std::endl;
210 exit(1);
211 }
212
213 g4Axis = param->GetAxis();
214 nofItems = param->GetNoDiv();
215 width = param->GetWidth();
216 offset = param->GetOffset();
217 halfGap = param->GetHalfGap();
218 }
219 else if (dynamic_cast<G4PVReplica*>(fPhysicalVolume)) {
220 bool consuming;
221
222 G4double offset0;
223 fPhysicalVolume->GetReplicationData(
224 g4Axis, nofItems, width, offset0, consuming);
225
226 // Different meaning of offset in R/Phi division
227 // for replica than for division:
228 // the offset includes start value (rmin, sphi)
229
230 G4double xlo = 0.;
231 VGM::ISolid* solid = Volume()->Solid();
232
233 if (g4Axis == kRho) {
234 if (solid->Type() == VGM::kTubs)
235 xlo = dynamic_cast<VGM::ITubs*>(solid)->InnerRadius();
236 }
237
238 if (g4Axis == kPhi) {
239 // Different meaning of offset in R/Phi division
240 // for replica than for division:
241 // the offset includes start value (rmin, sphi)
242
243 if (solid->Type() == VGM::kCons)
244 xlo = dynamic_cast<VGM::ICons*>(solid)->StartPhi();
245
246 if (solid->Type() == VGM::kTubs)
247 xlo = dynamic_cast<VGM::ITubs*>(solid)->StartPhi();
248
249 if (solid->Type() == VGM::kPolycone)
250 xlo = dynamic_cast<VGM::IPolycone*>(solid)->StartPhi();
251
252 if (solid->Type() == VGM::kPolyhedra)
253 xlo = dynamic_cast<VGM::IPolyhedra*>(solid)->StartPhi();
254 }
255
256 offset = offset0 - xlo;
257 }
258 axis = GetAxis(g4Axis);
259
260 // Convert units
261 offset *= ClhepVGM::Units::AxisUnit(axis);
262 width *= ClhepVGM::Units::AxisUnit(axis);
263 halfGap *= ClhepVGM::Units::AxisUnit(axis);
264
265 return true;
266}
267
268//_____________________________________________________________________________
270 std::vector<VGM::Transform>& transforms,
271 std::vector<VGM::IVolume*>& volumes) const
272{
273 // Fill data only if parameterised placement
274 if (Type() != VGM::kParameterised) return false;
275
276 G4PVParameterised* paraPhysicalVolume =
277 dynamic_cast<G4PVParameterised*>(fPhysicalVolume);
278 if (!paraPhysicalVolume) {
279 std::cerr << " Geant4GM::Placement::ParameterisedPlacementData: "
280 << std::endl;
281 std::cerr << " Physical volume is not parameterised..." << std::endl;
282 std::cerr << "*** Error: Aborting execution ***" << std::endl;
283 exit(1);
284 }
285
286 // Get parameterised volume elements saved in the map with import
287 const auto& g4ParamVolumes = Geant4GM::VolumeMap::Instance()->
288 GetParamVolumes(paraPhysicalVolume->GetLogicalVolume());
289
290 G4VPVParameterisation* pPara = paraPhysicalVolume->GetParameterisation();
291 for (int i = 0; i < paraPhysicalVolume->GetMultiplicity(); ++i) {
292 // transformations
293 pPara->ComputeTransformation(i, paraPhysicalVolume);
294 transforms.emplace_back(
295 ClhepVGM::Transform(*paraPhysicalVolume->GetObjectRotation(),
296 paraPhysicalVolume->GetObjectTranslation()));
297
298 // volumes
299 auto paramVolume = Geant4GM::VolumeMap::Instance()->GetVolume(g4ParamVolumes[i]);
300 volumes.emplace_back(paramVolume);
301 }
302
303 return true;
304}
virtual VGM::ISolid * Solid() const
Return the associated solid.
Definition VVolume.cxx:73
static double AxisUnit(VGM::Axis axis)
Convert CLHEP default unit for given axis type in VGM units.
Definition Units.cxx:63
static PlacementMap * Instance()
void AddPlacement(VGM::IPlacement *, G4VPhysicalVolume *)
VGM implementation for Geant4 positions of volumes.
Definition Placement.h:40
virtual VGM::Transform Transformation() const
Return the 3D transformation (if simple placement)
virtual VGM::PlacementType Type() const
Return the type of this placement.
static EAxis GetAxis(VGM::Axis axis)
Definition Placement.cxx:82
virtual bool MultiplePlacementData(VGM::Axis &axis, int &nofItems, double &width, double &offset, double &halfGap) const
Fill the multiple placement data if relevant and return true; return false if not multiple placement.
virtual bool ParameterisedPlacementData(std::vector< VGM::Transform > &transforms, std::vector< VGM::IVolume * > &volumes) const
Fill the parameterised placement data if relevant and return true; return false if not parameterised ...
virtual int CopyNo() const
Return the copy number of this placement.
virtual std::string Name() const
Return the name of this placement.
static VolumeMap * Instance()
Definition VolumeMap.cxx:28
G4LogicalVolume * GetVolume(VGM::IVolume *iVolume) const
Definition VolumeMap.cxx:95
VGM implementation for Geant4 volume.
Definition Volume.h:36
VGM Axis enumeration.
The VGM interface to cons solids.
Definition ICons.h:30
virtual double StartPhi() const =0
Return the starting phi angle of the segment in deg.
The VGM interface to polycone solids.
Definition IPolycone.h:30
virtual double StartPhi() const =0
Return starting phi angle of the segment in deg.
The VGM interface to polyhedra solids.
Definition IPolyhedra.h:30
virtual double StartPhi() const =0
Return starting phi angle of the segment in deg.
The VGM interface to solids.
Definition ISolid.h:58
virtual SolidType Type() const =0
Return the type of this solid.
The VGM interface to tubs solids.
Definition ITubs.h:30
virtual double StartPhi() const =0
Return the starting angle of the segment in deg.
The VGM interface to volumes.
Definition IVolume.h:32
BaseVGM utilities.
Definition utilities.h:23
VGM::Transform Transform(const CLHEP::HepRotation &rotation, const CLHEP::Hep3Vector &translation)
Definition transform.cxx:26
VGM interfaces.
Definition VMedium.h:28
std::vector< double > Transform
Definition Transform.h:40
@ kRho
Definition Axis.h:38
@ kPhi
Definition Axis.h:40
@ kYAxis
Definition Axis.h:36
@ kRadial3D
Definition Axis.h:39
@ kXAxis
Definition Axis.h:35
@ kUnknownAxis
Definition Axis.h:42
@ kZAxis
Definition Axis.h:37
PlacementType
Definition IPlacement.h:32
@ kParameterised
Definition IPlacement.h:35
@ kSimplePlacement
Definition IPlacement.h:33
@ kMultiplePlacement
Definition IPlacement.h:34
@ kCons
Definition ISolid.h:32
@ kTubs
Definition ISolid.h:47
@ kPolycone
Definition ISolid.h:40
@ kPolyhedra
Definition ISolid.h:41