Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// Class G4PVParameterised implementation
27//
28// Author: Paul Kent (CERN), 29 July 1995 - first non-stub version
29// ----------------------------------------------------------------------
30
31#include "G4PVParameterised.hh"
33#include "G4AffineTransform.hh"
34#include "G4UnitsTable.hh"
35#include "G4VSolid.hh"
36#include "G4LogicalVolume.hh"
37
38// ----------------------------------------------------------------------
39// Constructor
40//
42 G4LogicalVolume* pLogical,
43 G4VPhysicalVolume* pMotherPhysical,
44 const EAxis pAxis,
45 const G4int nReplicas,
47 G4bool pSurfChk )
48: G4PVReplica(pName, nReplicas, pAxis, pLogical,
49 pMotherPhysical != nullptr ? pMotherPhysical->GetLogicalVolume() : nullptr ),
50 fparam(pParam)
51{
52 G4LogicalVolume* motherLogical= pMotherPhysical != nullptr ?
53 pMotherPhysical->GetLogicalVolume() : nullptr;
54
55 SetMotherLogical( motherLogical );
56 if( motherLogical != nullptr )
57 {
58 // Registration moved here to ensure that the volume is recognised as Parameterised
59 motherLogical->AddDaughter(this);
60 }
61
62#ifdef G4VERBOSE
63 if ((pMotherPhysical != nullptr) && (pMotherPhysical->IsParameterised()))
64 {
65 std::ostringstream message;
66 std::ostringstream hint;
67 message << "A parameterised volume is being placed" << G4endl
68 << "inside another parameterised volume !";
69 hint << "To make sure that no overlaps are generated," << G4endl
70 << "you should verify the mother replicated shapes" << G4endl
71 << "are of the same type and dimensions." << G4endl
72 << " Mother physical volume: " << pMotherPhysical->GetName() << G4endl
73 << " Parameterised volume: " << pName << G4endl
74 << " (To switch this warning off, compile with G4_NO_VERBOSE)";
75 G4Exception("G4PVParameterised::G4PVParameterised()", "GeomVol1002",
76 JustWarning, message, G4String(hint.str()));
77 }
78#endif
79 if (pSurfChk) { CheckOverlaps(); }
80}
81
82// ----------------------------------------------------------------------
83// Constructor
84//
86 G4LogicalVolume* pLogical,
87 G4LogicalVolume* pMotherLogical,
88 const EAxis pAxis,
89 const G4int nReplicas,
91 G4bool pSurfChk )
92 : G4PVReplica(pName, nReplicas, pAxis, pLogical, pMotherLogical ),
93 fparam(pParam)
94{
95 SetMotherLogical( pMotherLogical );
96 if( pMotherLogical != nullptr )
97 {
98 // Registration moved here to ensure that the volume is recognised as Parameterised
99 pMotherLogical->AddDaughter(this);
100 }
101 if (pSurfChk) { CheckOverlaps(); }
102}
103
104// ----------------------------------------------------------------------
105// Fake default constructor - sets only member data and allocates memory
106// for usage restricted to object persistency.
107//
109 : G4PVReplica(a)
110{
111}
112
113// ----------------------------------------------------------------------
114// GetParameterisation
115//
120
121// ----------------------------------------------------------------------
122// IsParameterised
123//
125{
126 return true;
127}
128
129// ----------------------------------------------------------------------
130// VolumeType
131//
136
137// ----------------------------------------------------------------------
138// GetReplicationData
139//
141 G4int& nReplicas,
142 G4double& width,
144 G4bool& consuming) const
145{
146 axis = faxis;
147 nReplicas = fnReplicas;
148 width = fwidth;
149 offset = foffset;
150 consuming = false;
151}
152
153// ----------------------------------------------------------------------
154// SetRegularStructureId
155//
157{
158 G4PVReplica::SetRegularStructureId( code );
159 // To undertake additional preparation, a derived volume must
160 // redefine this method, while calling also the above method
161}
162
163
164// ----------------------------------------------------------------------
165// CheckOverlaps
166//
167G4bool
169 G4bool verbose, G4int maxErr)
170{
171 if (res<=0) { return false; }
172
173 G4int trials = 0;
174 G4bool retval = false;
175 G4VSolid* solidA = nullptr;
176 G4VSolid* solidB = nullptr;
177 G4LogicalVolume* motherLog = GetMotherLogical();
178 G4VSolid* motherSolid = motherLog->GetSolid();
179 std::vector<G4ThreeVector> points;
180
181 if (verbose)
182 {
183 G4cout << "Checking overlaps for parameterised volume "
184 << GetName() << " ... ";
185 }
186
187 for (auto i=0; i<GetMultiplicity(); ++i)
188 {
189 solidA = fparam->ComputeSolid(i, this);
190 solidA->ComputeDimensions(fparam, i, this);
191 fparam->ComputeTransformation(i, this);
192
193 // Create the transformation from daughter to mother
194 //
195 G4AffineTransform Tm( GetRotation(), GetTranslation() );
196
197 // Generate random points on surface according to the given resolution,
198 // transform them to the mother's coordinate system and if no overlaps
199 // with the mother volume, cache them in a vector for later use with
200 // the daughters
201 //
202 for (auto n=0; n<res; ++n)
203 {
205
206 // Checking overlaps with the mother volume
207 //
208 if (motherSolid->Inside(mp)==kOutside)
209 {
210 G4double distin = motherSolid->DistanceToIn(mp);
211 if (distin > tol)
212 {
213 ++trials; retval = true;
214 std::ostringstream message;
215 message << "Overlap with mother volume !" << G4endl
216 << " Overlap is detected for volume "
217 << GetName() << ", parameterised instance: " << i << G4endl
218 << " with its mother volume "
219 << motherLog->GetName() << G4endl
220 << " at mother local point " << mp << ", "
221 << "overlapping by at least: "
222 << G4BestUnit(distin, "Length");
223 if (trials>=maxErr)
224 {
225 message << G4endl
226 << "NOTE: Reached maximum fixed number -" << maxErr
227 << "- of overlaps reports for this volume !";
228 }
229 G4Exception("G4PVParameterised::CheckOverlaps()",
230 "GeomVol1002", JustWarning, message);
231 if (trials>=maxErr) { return true; }
232 }
233 }
234 points.push_back(mp);
235 }
236
237 // Checking overlaps with each other parameterised instance
238 //
239 for (auto j=i+1; j<GetMultiplicity(); ++j)
240 {
241 solidB = fparam->ComputeSolid(j,this);
242 solidB->ComputeDimensions(fparam, j, this);
243 fparam->ComputeTransformation(j, this);
244
245 // Create the transformation for daughter volume
246 //
247 G4AffineTransform Td( GetRotation(), GetTranslation() );
248
249 for (const auto & point : points)
250 {
251 // Transform each point according to daughter's frame
252 //
254
255 if (solidB->Inside(md)==kInside)
256 {
257 G4double distout = solidB->DistanceToOut(md);
258 if (distout > tol)
259 {
260 ++trials; retval = true;
261 std::ostringstream message;
262 message << "Overlap within parameterised volumes !" << G4endl
263 << " Overlap is detected for volume "
264 << GetName() << ", parameterised instance: " << i << G4endl
265 << " with parameterised volume instance: " << j
266 << G4endl
267 << " at local point " << md << ", "
268 << "overlapping by at least: "
269 << G4BestUnit(distout, "Length")
270 << ", related to volume instance: " << j << ".";
271 if (trials>=maxErr)
272 {
273 message << G4endl
274 << "NOTE: Reached maximum fixed number -" << maxErr
275 << "- of overlaps reports for this volume !";
276 }
277 G4Exception("G4PVParameterised::CheckOverlaps()",
278 "GeomVol1002", JustWarning, message);
279 if (trials>=maxErr) { return true; }
280 }
281 }
282 }
283 }
284 }
285 if (verbose)
286 {
287 G4cout << "OK! " << G4endl;
288 }
289
290 return retval;
291}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ThreadLocal T * G4GeomSplitter< T >::offset
EAxis faxis
G4int fnReplicas
G4double fwidth
G4double foffset
G4int GetMultiplicity() const override
#define G4BestUnit(a, b)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector InverseTransformPoint(const G4ThreeVector &vec) const
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
G4VSolid * GetSolid() const
void AddDaughter(G4VPhysicalVolume *p)
const G4String & GetName() const
G4bool IsParameterised() const override
void SetRegularStructureId(G4int code) override
G4VPVParameterisation * GetParameterisation() const override
EVolume VolumeType() const final
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1) override
void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const override
G4PVParameterised(const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, G4VPVParameterisation *pParam, G4bool pSurfChk=false)
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
virtual G4bool IsParameterised() const =0
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition G4VSolid.cc:136
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual G4ThreeVector GetPointOnSurface() const
Definition G4VSolid.cc:151
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
EAxis
Definition geomdefs.hh:54
@ kInside
Definition geomdefs.hh:70
@ kOutside
Definition geomdefs.hh:68
EVolume
Definition geomdefs.hh:83
@ kParameterised
Definition geomdefs.hh:86
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668