Geant4-11
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// 29.07.95, P.Kent - 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 ? pMotherPhysical->GetLogicalVolume() : nullptr ),
50 fparam(pParam)
51{
52 G4LogicalVolume* motherLogical= pMotherPhysical ?
53 pMotherPhysical->GetLogicalVolume() : nullptr;
54
55 SetMotherLogical( motherLogical );
56 if( motherLogical )
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, hint;
66 message << "A parameterised volume is being placed" << G4endl
67 << "inside another parameterised volume !";
68 hint << "To make sure that no overlaps are generated," << G4endl
69 << "you should verify the mother replicated shapes" << G4endl
70 << "are of the same type and dimensions." << G4endl
71 << " Mother physical volume: " << pMotherPhysical->GetName() << G4endl
72 << " Parameterised volume: " << pName << G4endl
73 << " (To switch this warning off, compile with G4_NO_VERBOSE)";
74 G4Exception("G4PVParameterised::G4PVParameterised()", "GeomVol1002",
75 JustWarning, message, G4String(hint.str()));
76 }
77#endif
78 if (pSurfChk) { CheckOverlaps(); }
79}
80
81// ----------------------------------------------------------------------
82// Constructor
83//
85 G4LogicalVolume* pLogical,
86 G4LogicalVolume* pMotherLogical,
87 const EAxis pAxis,
88 const G4int nReplicas,
90 G4bool pSurfChk )
91 : G4PVReplica(pName, nReplicas, pAxis, pLogical, pMotherLogical ),
92 fparam(pParam)
93{
94 SetMotherLogical( pMotherLogical );
95 if( pMotherLogical )
96 {
97 // Registration moved here to ensure that the volume is recognised as Parameterised
98 pMotherLogical->AddDaughter(this);
99 }
100 if (pSurfChk) { CheckOverlaps(); }
101}
102
103// ----------------------------------------------------------------------
104// Fake default constructor - sets only member data and allocates memory
105// for usage restricted to object persistency.
106//
108 : G4PVReplica(a)
109{
110}
111
112// ----------------------------------------------------------------------
113// Destructor
114//
116{
117}
118
119// ----------------------------------------------------------------------
120// GetParameterisation
121//
123{
124 return fparam;
125}
126
127// ----------------------------------------------------------------------
128// IsParameterised
129//
131{
132 return true;
133}
134
135// ----------------------------------------------------------------------
136// VolumeType
137//
139{
140 return kParameterised;
141}
142
143// ----------------------------------------------------------------------
144// GetReplicationData
145//
147 G4int& nReplicas,
148 G4double& width,
149 G4double& offset,
150 G4bool& consuming) const
151{
152 axis = faxis;
153 nReplicas = fnReplicas;
154 width = fwidth;
155 offset = foffset;
156 consuming = false;
157}
158
159// ----------------------------------------------------------------------
160// SetRegularStructureId
161//
163{
165 // To undertake additional preparation, a derived volume must
166 // redefine this method, while calling also the above method
167}
168
169
170// ----------------------------------------------------------------------
171// CheckOverlaps
172//
173G4bool
175 G4bool verbose, G4int maxErr)
176{
177 if (res<=0) { return false; }
178
179 G4int trials = 0;
180 G4bool retval = false;
181 G4VSolid *solidA = nullptr, *solidB = nullptr;
182 G4LogicalVolume* motherLog = GetMotherLogical();
183 G4VSolid *motherSolid = motherLog->GetSolid();
184 std::vector<G4ThreeVector> points;
185
186 if (verbose)
187 {
188 G4cout << "Checking overlaps for parameterised volume "
189 << GetName() << " ... ";
190 }
191
192 for (auto i=0; i<GetMultiplicity(); ++i)
193 {
194 solidA = fparam->ComputeSolid(i, this);
195 solidA->ComputeDimensions(fparam, i, this);
197
198 // Create the transformation from daughter to mother
199 //
201
202 // Generate random points on surface according to the given resolution,
203 // transform them to the mother's coordinate system and if no overlaps
204 // with the mother volume, cache them in a vector for later use with
205 // the daughters
206 //
207 for (auto n=0; n<res; ++n)
208 {
210
211 // Checking overlaps with the mother volume
212 //
213 if (motherSolid->Inside(mp)==kOutside)
214 {
215 G4double distin = motherSolid->DistanceToIn(mp);
216 if (distin > tol)
217 {
218 ++trials; retval = true;
219 std::ostringstream message;
220 message << "Overlap with mother volume !" << G4endl
221 << " Overlap is detected for volume "
222 << GetName() << ", parameterised instance: " << i << G4endl
223 << " with its mother volume "
224 << motherLog->GetName() << G4endl
225 << " at mother local point " << mp << ", "
226 << "overlapping by at least: "
227 << G4BestUnit(distin, "Length");
228 if (trials>=maxErr)
229 {
230 message << G4endl
231 << "NOTE: Reached maximum fixed number -" << maxErr
232 << "- of overlaps reports for this volume !";
233 }
234 G4Exception("G4PVParameterised::CheckOverlaps()",
235 "GeomVol1002", JustWarning, message);
236 if (trials>=maxErr) { return true; }
237 }
238 }
239 points.push_back(mp);
240 }
241
242 // Checking overlaps with each other parameterised instance
243 //
244 for (auto j=i+1; j<GetMultiplicity(); ++j)
245 {
246 solidB = fparam->ComputeSolid(j,this);
247 solidB->ComputeDimensions(fparam, j, this);
249
250 // Create the transformation for daughter volume
251 //
253
254 for (auto pos=points.cbegin(); pos!=points.cend(); ++pos)
255 {
256 // Transform each point according to daughter's frame
257 //
259
260 if (solidB->Inside(md)==kInside)
261 {
262 G4double distout = solidB->DistanceToOut(md);
263 if (distout > tol)
264 {
265 ++trials; retval = true;
266 std::ostringstream message;
267 message << "Overlap within parameterised volumes !" << G4endl
268 << " Overlap is detected for volume "
269 << GetName() << ", parameterised instance: " << i << G4endl
270 << " with parameterised volume instance: " << j
271 << G4endl
272 << " at local point " << md << ", "
273 << "overlapping by at least: "
274 << G4BestUnit(distout, "Length")
275 << ", related to volume instance: " << j << ".";
276 if (trials>=maxErr)
277 {
278 message << G4endl
279 << "NOTE: Reached maximum fixed number -" << maxErr
280 << "- of overlaps reports for this volume !";
281 }
282 G4Exception("G4PVParameterised::CheckOverlaps()",
283 "GeomVol1002", JustWarning, message);
284 if (trials>=maxErr) { return true; }
285 }
286 }
287 }
288 }
289 }
290 if (verbose)
291 {
292 G4cout << "OK! " << G4endl;
293 }
294
295 return retval;
296}
static const G4double pos
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector InverseTransformPoint(const G4ThreeVector &vec) const
G4VSolid * GetSolid() const
void AddDaughter(G4VPhysicalVolume *p)
const G4String & GetName() const
G4VPVParameterisation * fparam
G4bool IsParameterised() const
void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
G4VPVParameterisation * GetParameterisation() const
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1)
EVolume VolumeType() const final
virtual void SetRegularStructureId(G4int code)
G4PVParameterised(const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, G4VPVParameterisation *pParam, G4bool pSurfChk=false)
virtual void SetRegularStructureId(G4int code)
Definition: G4PVReplica.cc:335
G4double fwidth
Definition: G4PVReplica.hh:187
virtual G4int GetMultiplicity() const
Definition: G4PVReplica.cc:297
G4int fnReplicas
Definition: G4PVReplica.hh:186
G4double foffset
Definition: G4PVReplica.hh:187
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0
G4LogicalVolume * GetMotherLogical() const
const G4RotationMatrix * GetRotation() const
const G4ThreeVector GetTranslation() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
virtual G4bool IsParameterised() const =0
void SetMotherLogical(G4LogicalVolume *pMother)
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
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
Definition: inftrees.h:24