Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4PVParameterised.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id$
28//
29//
30// class G4PVParameterised
31//
32// Implementation
33//
34// ----------------------------------------------------------------------
35
36#include "G4PVParameterised.hh"
38#include "G4AffineTransform.hh"
39#include "G4UnitsTable.hh"
40#include "G4VSolid.hh"
41#include "G4LogicalVolume.hh"
42
43// ----------------------------------------------------------------------
44// Constructor
45//
47 G4LogicalVolume* pLogical,
48 G4VPhysicalVolume* pMother,
49 const EAxis pAxis,
50 const G4int nReplicas,
52 G4bool pSurfChk )
53 : G4PVReplica(pName, pLogical, pMother, pAxis, nReplicas, 0, 0),
54 fparam(pParam)
55{
56#ifdef G4VERBOSE
57 if ((pMother) && (pMother->IsParameterised()))
58 {
59 std::ostringstream message, hint;
60 message << "A parameterised volume is being placed" << G4endl
61 << "inside another parameterised volume !";
62 hint << "To make sure that no overlaps are generated," << G4endl
63 << "you should verify the mother replicated shapes" << G4endl
64 << "are of the same type and dimensions." << G4endl
65 << " Mother physical volume: " << pMother->GetName() << G4endl
66 << " Parameterised volume: " << pName << G4endl
67 << " (To switch this warning off, compile with G4_NO_VERBOSE)";
68 G4Exception("G4PVParameterised::G4PVParameterised()", "GeomVol1002",
69 JustWarning, message, G4String(hint.str()));
70 }
71#endif
72 if (pSurfChk) { CheckOverlaps(); }
73}
74
75// ----------------------------------------------------------------------
76// Constructor
77//
79 G4LogicalVolume* pLogical,
80 G4LogicalVolume* pMotherLogical,
81 const EAxis pAxis,
82 const G4int nReplicas,
84 G4bool pSurfChk )
85 : G4PVReplica(pName, pLogical, pMotherLogical, pAxis, nReplicas, 0, 0),
86 fparam(pParam)
87{
88 if (pSurfChk) { CheckOverlaps(); }
89}
90
91// ----------------------------------------------------------------------
92// Fake default constructor - sets only member data and allocates memory
93// for usage restricted to object persistency.
94//
96 : G4PVReplica(a), fparam(0)
97{
98}
99
100// ----------------------------------------------------------------------
101// Destructor
102//
104{
105}
106
107// ----------------------------------------------------------------------
108// GetParameterisation
109//
111{
112 return fparam;
113}
114
115// ----------------------------------------------------------------------
116// IsParameterised
117//
119{
120 return true;
121}
122
123// ----------------------------------------------------------------------
124// GetReplicationData
125//
127 G4int& nReplicas,
128 G4double& width,
129 G4double& offset,
130 G4bool& consuming) const
131{
132 axis = faxis;
133 nReplicas = fnReplicas;
134 width = fwidth;
135 offset = foffset;
136 consuming = false;
137}
138
139// ----------------------------------------------------------------------
140// SetRegularStructureId
141//
143{
145 // To undertake additional preparation, a derived volume must
146 // redefine this method, while calling also the above method.
147}
148
149
150// ----------------------------------------------------------------------
151// CheckOverlaps
152//
153G4bool
155{
156 if (res<=0) { return false; }
157
158 G4VSolid *solidA = 0, *solidB = 0;
159 G4LogicalVolume *motherLog = GetMotherLogical();
160 G4VSolid *motherSolid = motherLog->GetSolid();
161 std::vector<G4ThreeVector> points;
162
163 if (verbose)
164 {
165 G4cout << "Checking overlaps for parameterised volume "
166 << GetName() << " ... ";
167 }
168
169 for (G4int i=0; i<GetMultiplicity(); i++)
170 {
171 solidA = fparam->ComputeSolid(i, this);
172 solidA->ComputeDimensions(fparam, i, this);
173 fparam->ComputeTransformation(i, this);
174
175 // Create the transformation from daughter to mother
176 //
178
179 // Generate random points on surface according to the given resolution,
180 // transform them to the mother's coordinate system and if no overlaps
181 // with the mother volume, cache them in a vector for later use with
182 // the daughters
183 //
184 for (G4int n=0; n<res; n++)
185 {
187
188 // Checking overlaps with the mother volume
189 //
190 if (motherSolid->Inside(mp)==kOutside)
191 {
192 G4double distin = motherSolid->DistanceToIn(mp);
193 if (distin > tol)
194 {
195 std::ostringstream message;
196 message << "Overlap with mother volume !" << G4endl
197 << " Overlap is detected for volume "
198 << GetName() << ", parameterised instance: " << i << G4endl
199 << " with its mother volume "
200 << motherLog->GetName() << G4endl
201 << " at mother local point " << mp << ", "
202 << "overlapping by at least: "
203 << G4BestUnit(distin, "Length");
204 G4Exception("G4PVParameterised::CheckOverlaps()",
205 "GeomVol1002", JustWarning, message);
206 return true;
207 }
208 }
209 points.push_back(mp);
210 }
211
212 // Checking overlaps with each other parameterised instance
213 //
214 std::vector<G4ThreeVector>::iterator pos;
215 for (G4int j=i+1; j<GetMultiplicity(); j++)
216 {
217 solidB = fparam->ComputeSolid(j,this);
218 solidB->ComputeDimensions(fparam, j, this);
219 fparam->ComputeTransformation(j, this);
220
221 // Create the transformation for daughter volume
222 //
224
225 for (pos=points.begin(); pos!=points.end(); pos++)
226 {
227 // Transform each point according to daughter's frame
228 //
229 G4ThreeVector md = Td.Inverse().TransformPoint(*pos);
230
231 if (solidB->Inside(md)==kInside)
232 {
233 G4double distout = solidB->DistanceToOut(md);
234 if (distout > tol)
235 {
236 std::ostringstream message;
237 message << "Overlap within parameterised volumes !" << G4endl
238 << " Overlap is detected for volume "
239 << GetName() << ", parameterised instance: " << i << G4endl
240 << " with parameterised volume instance: " << j
241 << G4endl
242 << " at local point " << md << ", "
243 << "overlapping by at least: "
244 << G4BestUnit(distout, "Length")
245 << ", related to volume instance: " << j << ".";
246 G4Exception("G4PVParameterised::CheckOverlaps()",
247 "GeomVol1002", JustWarning, message);
248 return true;
249 }
250 }
251 }
252 }
253 }
254 if (verbose)
255 {
256 G4cout << "OK! " << G4endl;
257 }
258
259 return false;
260}
@ JustWarning
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4AffineTransform Inverse() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4VSolid * GetSolid() const
G4String GetName() const
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true)
G4bool IsParameterised() const
void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
G4VPVParameterisation * GetParameterisation() const
virtual void SetRegularStructureId(G4int Code)
G4PVParameterised(const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, G4VPVParameterisation *pParam, G4bool pSurfChk=false)
G4double fwidth
Definition: G4PVReplica.hh:136
virtual G4int GetMultiplicity() const
Definition: G4PVReplica.cc:212
G4int fnReplicas
Definition: G4PVReplica.hh:135
G4double foffset
Definition: G4PVReplica.hh:136
virtual void SetRegularStructureId(G4int Code)
Definition: G4PVReplica.cc:242
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0
G4LogicalVolume * GetMotherLogical() const
const G4RotationMatrix * GetRotation() const
const G4String & GetName() const
const G4ThreeVector & GetTranslation() const
virtual G4bool IsParameterised() const =0
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
#define Code
Definition: deflate.h:70
EAxis
Definition: geomdefs.hh:54
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41