Geant4  10.00.p03
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: G4PVParameterised.cc 73250 2013-08-22 13:22:23Z gcosmo $
28 //
29 //
30 // class G4PVParameterised
31 //
32 // Implementation
33 //
34 // ----------------------------------------------------------------------
35 
36 #include "G4PVParameterised.hh"
37 #include "G4VPVParameterisation.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,
51  G4VPVParameterisation *pParam,
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,
83  G4VPVParameterisation *pParam,
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 //
153 G4bool
155  G4bool verbose, G4int maxErr)
156 {
157  if (res<=0) { return false; }
158 
159  G4int trials = 0;
160  G4bool retval = false;
161  G4VSolid *solidA = 0, *solidB = 0;
162  G4LogicalVolume *motherLog = GetMotherLogical();
163  G4VSolid *motherSolid = motherLog->GetSolid();
164  std::vector<G4ThreeVector> points;
165 
166  if (verbose)
167  {
168  G4cout << "Checking overlaps for parameterised volume "
169  << GetName() << " ... ";
170  }
171 
172  for (G4int i=0; i<GetMultiplicity(); i++)
173  {
174  solidA = fparam->ComputeSolid(i, this);
175  solidA->ComputeDimensions(fparam, i, this);
176  fparam->ComputeTransformation(i, this);
177 
178  // Create the transformation from daughter to mother
179  //
181 
182  // Generate random points on surface according to the given resolution,
183  // transform them to the mother's coordinate system and if no overlaps
184  // with the mother volume, cache them in a vector for later use with
185  // the daughters
186  //
187  for (G4int n=0; n<res; n++)
188  {
190 
191  // Checking overlaps with the mother volume
192  //
193  if (motherSolid->Inside(mp)==kOutside)
194  {
195  G4double distin = motherSolid->DistanceToIn(mp);
196  if (distin > tol)
197  {
198  trials++; retval = true;
199  std::ostringstream message;
200  message << "Overlap with mother volume !" << G4endl
201  << " Overlap is detected for volume "
202  << GetName() << ", parameterised instance: " << i << G4endl
203  << " with its mother volume "
204  << motherLog->GetName() << G4endl
205  << " at mother local point " << mp << ", "
206  << "overlapping by at least: "
207  << G4BestUnit(distin, "Length");
208  if (trials>=maxErr)
209  {
210  message << G4endl
211  << "NOTE: Reached maximum fixed number -" << maxErr
212  << "- of overlaps reports for this volume !";
213  }
214  G4Exception("G4PVParameterised::CheckOverlaps()",
215  "GeomVol1002", JustWarning, message);
216  if (trials>=maxErr) { return true; }
217  }
218  }
219  points.push_back(mp);
220  }
221 
222  // Checking overlaps with each other parameterised instance
223  //
224  std::vector<G4ThreeVector>::iterator pos;
225  for (G4int j=i+1; j<GetMultiplicity(); j++)
226  {
227  solidB = fparam->ComputeSolid(j,this);
228  solidB->ComputeDimensions(fparam, j, this);
229  fparam->ComputeTransformation(j, this);
230 
231  // Create the transformation for daughter volume
232  //
234 
235  for (pos=points.begin(); pos!=points.end(); pos++)
236  {
237  // Transform each point according to daughter's frame
238  //
239  G4ThreeVector md = Td.Inverse().TransformPoint(*pos);
240 
241  if (solidB->Inside(md)==kInside)
242  {
243  G4double distout = solidB->DistanceToOut(md);
244  if (distout > tol)
245  {
246  trials++; retval = true;
247  std::ostringstream message;
248  message << "Overlap within parameterised volumes !" << G4endl
249  << " Overlap is detected for volume "
250  << GetName() << ", parameterised instance: " << i << G4endl
251  << " with parameterised volume instance: " << j
252  << G4endl
253  << " at local point " << md << ", "
254  << "overlapping by at least: "
255  << G4BestUnit(distout, "Length")
256  << ", related to volume instance: " << j << ".";
257  if (trials>=maxErr)
258  {
259  message << G4endl
260  << "NOTE: Reached maximum fixed number -" << maxErr
261  << "- of overlaps reports for this volume !";
262  }
263  G4Exception("G4PVParameterised::CheckOverlaps()",
264  "GeomVol1002", JustWarning, message);
265  if (trials>=maxErr) { return true; }
266  }
267  }
268  }
269  }
270  }
271  if (verbose)
272  {
273  G4cout << "OK! " << G4endl;
274  }
275 
276  return retval;
277 }
const G4ThreeVector & GetTranslation() const
G4String GetName() const
CLHEP::Hep3Vector G4ThreeVector
virtual void SetRegularStructureId(G4int Code)
Definition: G4PVReplica.cc:259
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4VSolid.cc:137
G4AffineTransform Inverse() const
virtual G4VSolid * ComputeSolid(const G4int, G4VPhysicalVolume *)
#define width
G4double fwidth
Definition: G4PVReplica.hh:194
G4bool IsParameterised() const
G4double a
Definition: TRTMaterials.hh:39
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
G4double foffset
Definition: G4PVReplica.hh:194
G4VPVParameterisation * GetParameterisation() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
G4VPVParameterisation * fparam
G4int fnReplicas
Definition: G4PVReplica.hh:193
bool G4bool
Definition: G4Types.hh:79
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
const G4int n
G4LogicalVolume * GetMotherLogical() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4bool IsParameterised() const =0
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
virtual void SetRegularStructureId(G4int Code)
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1)
void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
EAxis
Definition: geomdefs.hh:54
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual void ComputeTransformation(const G4int, G4VPhysicalVolume *) const =0
virtual G4int GetMultiplicity() const
Definition: G4PVReplica.cc:229
const G4RotationMatrix * GetRotation() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4PVParameterised(const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMotherLogical, const EAxis pAxis, const G4int nReplicas, G4VPVParameterisation *pParam, G4bool pSurfChk=false)
static const G4double pos
G4VSolid * GetSolid() const