Geant4  10.02.p01
G4PVPlacement.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: G4PVPlacement.cc 85846 2014-11-05 15:45:28Z gcosmo $
28 //
29 //
30 // class G4PVPlacement Implementation
31 //
32 // ----------------------------------------------------------------------
33 
34 #include "G4PVPlacement.hh"
35 #include "G4AffineTransform.hh"
36 #include "G4UnitsTable.hh"
37 #include "G4LogicalVolume.hh"
38 #include "G4VSolid.hh"
39 
40 // ----------------------------------------------------------------------
41 // Constructor
42 //
44  const G4ThreeVector &tlate,
45  const G4String& pName,
46  G4LogicalVolume *pLogical,
47  G4VPhysicalVolume *pMother,
48  G4bool pMany,
49  G4int pCopyNo,
50  G4bool pSurfChk )
51  : G4VPhysicalVolume(pRot,tlate,pName,pLogical,pMother),
52  fmany(pMany), fallocatedRotM(false), fcopyNo(pCopyNo)
53 {
54  if (pMother)
55  {
56  G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
57  if (pLogical == motherLogical)
58  {
59  G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
60  FatalException, "Cannot place a volume inside itself!");
61  }
62  SetMotherLogical(motherLogical);
63  motherLogical->AddDaughter(this);
64  if (pSurfChk) { CheckOverlaps(); }
65  }
66 }
67 
68 // ----------------------------------------------------------------------
69 // Constructor
70 //
72  const G4String& pName,
73  G4LogicalVolume *pLogical,
74  G4VPhysicalVolume *pMother,
75  G4bool pMany,
76  G4int pCopyNo,
77  G4bool pSurfChk )
78  : G4VPhysicalVolume(NewPtrRotMatrix(Transform3D.getRotation().inverse()),
79  Transform3D.getTranslation(),pName,pLogical,pMother),
80  fmany(pMany), fcopyNo(pCopyNo)
81 {
82  fallocatedRotM = (GetRotation() != 0);
83  if (pMother)
84  {
85  G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
86  if (pLogical == motherLogical)
87  G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
88  FatalException, "Cannot place a volume inside itself!");
89  SetMotherLogical(motherLogical);
90  motherLogical->AddDaughter(this);
91  if (pSurfChk) { CheckOverlaps(); }
92  }
93 }
94 
95 // ----------------------------------------------------------------------
96 // Constructor
97 //
98 // The logical volume of the mother is utilised (not the physical)
99 //
101  const G4ThreeVector &tlate,
102  G4LogicalVolume *pCurrentLogical,
103  const G4String& pName,
104  G4LogicalVolume *pMotherLogical,
105  G4bool pMany,
106  G4int pCopyNo,
107  G4bool pSurfChk )
108  : G4VPhysicalVolume(pRot,tlate,pName,pCurrentLogical,0),
109  fmany(pMany), fallocatedRotM(false), fcopyNo(pCopyNo)
110 {
111  if (pCurrentLogical == pMotherLogical)
112  {
113  G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
114  FatalException, "Cannot place a volume inside itself!");
115  }
116  SetMotherLogical(pMotherLogical);
117  if (pMotherLogical) { pMotherLogical->AddDaughter(this); }
118  if ((pSurfChk) && (pMotherLogical)) { CheckOverlaps(); }
119 }
120 
121 
122 // ----------------------------------------------------------------------
123 // Constructor
124 //
126  G4LogicalVolume *pCurrentLogical,
127  const G4String& pName,
128  G4LogicalVolume *pMotherLogical,
129  G4bool pMany,
130  G4int pCopyNo,
131  G4bool pSurfChk )
132  : G4VPhysicalVolume(0,Transform3D.getTranslation(),pName,pCurrentLogical,0),
133  fmany(pMany), fcopyNo(pCopyNo)
134 {
135  if (pCurrentLogical == pMotherLogical)
136  {
137  G4Exception("G4PVPlacement::G4PVPlacement()", "GeomVol0002",
138  FatalException, "Cannot place a volume inside itself!");
139  }
140  SetRotation( NewPtrRotMatrix(Transform3D.getRotation().inverse()) );
141  fallocatedRotM = (GetRotation() != 0);
142  SetMotherLogical(pMotherLogical);
143  if (pMotherLogical) { pMotherLogical->AddDaughter(this); }
144  if ((pSurfChk) && (pMotherLogical)) { CheckOverlaps(); }
145 }
146 
147 // ----------------------------------------------------------------------
148 // Fake default constructor - sets only member data and allocates memory
149 // for usage restricted to object persistency.
150 //
152  : G4VPhysicalVolume(a), fmany(false), fallocatedRotM(0), fcopyNo(0)
153 {
154 }
155 
156 // ----------------------------------------------------------------------
157 // Destructor
158 //
160 {
161  if( fallocatedRotM ){ delete this->GetRotation() ; }
162 }
163 
164 // ----------------------------------------------------------------------
165 // IsMany
166 //
168 {
169  return fmany;
170 }
171 
172 // ----------------------------------------------------------------------
173 // SetCopyNo
174 //
176 {
177  fcopyNo= newCopyNo;
178 }
179 
180 // ----------------------------------------------------------------------
181 // IsReplicated
182 //
184 {
185  return false;
186 }
187 
188 // ----------------------------------------------------------------------
189 // IsParameterised
190 //
192 {
193  return false;
194 }
195 
196 // ----------------------------------------------------------------------
197 // GetParameterisation
198 //
200 {
201  return 0;
202 }
203 
204 // ----------------------------------------------------------------------
205 // GetReplicationData
206 //
207 void G4PVPlacement::
209 {
210  // No-operations
211 }
212 
213 // ----------------------------------------------------------------------
214 // IsRegularRepeatedStructure
215 //
216 // This is for specialised repeated volumes (replicas, parameterised vol.)
217 //
219 {
220  return false;
221 }
222 
223 // ----------------------------------------------------------------------
224 // IsRegularRepeatedStructure
225 //
226 // This is for specialised repeated volumes (replicas, parameterised vol.)
227 //
229 {
230  return 0;
231 }
232 
233 // ----------------------------------------------------------------------
234 // CheckOverlaps
235 //
237  G4bool verbose, G4int maxErr)
238 {
239  if (res<=0) { return false; }
240 
241  G4VSolid* solid = GetLogicalVolume()->GetSolid();
242  G4LogicalVolume* motherLog = GetMotherLogical();
243  if (!motherLog) { return false; }
244 
245  G4int trials = 0;
246  G4bool retval = false;
247 
248  G4VSolid* motherSolid = motherLog->GetSolid();
249 
250  if (verbose)
251  {
252  G4cout << "Checking overlaps for volume " << GetName() << " ... ";
253  }
254 
255  // Create the transformation from daughter to mother
256  //
258 
259  for (G4int n=0; n<res; n++)
260  {
261  // Generate a random point on the solid's surface
262  //
263  G4ThreeVector point = solid->GetPointOnSurface();
264 
265  // Transform the generated point to the mother's coordinate system
266  //
267  G4ThreeVector mp = Tm.TransformPoint(point);
268 
269  // Checking overlaps with the mother volume
270  //
271  if (motherSolid->Inside(mp)==kOutside)
272  {
273  G4double distin = motherSolid->DistanceToIn(mp);
274  if (distin > tol)
275  {
276  trials++; retval = true;
277  std::ostringstream message;
278  message << "Overlap with mother volume !" << G4endl
279  << " Overlap is detected for volume "
280  << GetName() << G4endl
281  << " with its mother volume "
282  << motherLog->GetName() << G4endl
283  << " at mother local point " << mp << ", "
284  << "overlapping by at least: "
285  << G4BestUnit(distin, "Length");
286  if (trials>=maxErr)
287  {
288  message << G4endl
289  << "NOTE: Reached maximum fixed number -" << maxErr
290  << "- of overlaps reports for this volume !";
291  }
292  G4Exception("G4PVPlacement::CheckOverlaps()",
293  "GeomVol1002", JustWarning, message);
294  if (trials>=maxErr) { return true; }
295  }
296  }
297 
298  // Checking overlaps with each 'sister' volume
299  //
300  for (G4int i=0; i<motherLog->GetNoDaughters(); i++)
301  {
302  G4VPhysicalVolume* daughter = motherLog->GetDaughter(i);
303 
304  if (daughter == this) { continue; }
305 
306  // Create the transformation for daughter volume and transform point
307  //
308  G4AffineTransform Td( daughter->GetRotation(),
309  daughter->GetTranslation() );
310  G4ThreeVector md = Td.Inverse().TransformPoint(mp);
311 
312  G4VSolid* daughterSolid = daughter->GetLogicalVolume()->GetSolid();
313  if (daughterSolid->Inside(md)==kInside)
314  {
315  G4double distout = daughterSolid->DistanceToOut(md);
316  if (distout > tol)
317  {
318  trials++; retval = true;
319  std::ostringstream message;
320  message << "Overlap with volume already placed !" << G4endl
321  << " Overlap is detected for volume "
322  << GetName() << G4endl
323  << " with " << daughter->GetName() << " volume's"
324  << G4endl
325  << " local point " << md << ", "
326  << "overlapping by at least: "
327  << G4BestUnit(distout,"Length");
328  if (trials>=maxErr)
329  {
330  message << G4endl
331  << "NOTE: Reached maximum fixed number -" << maxErr
332  << "- of overlaps reports for this volume !";
333  }
334  G4Exception("G4PVPlacement::CheckOverlaps()",
335  "GeomVol1002", JustWarning, message);
336  if (trials>=maxErr) { return true; }
337  }
338  }
339 
340  // Now checking that 'sister' volume is not totally included and
341  // overlapping. Do it only once, for the first point generated
342  //
343  if (n==0)
344  {
345  // Generate a single point on the surface of the 'sister' volume
346  // and verify that the point is NOT inside the current volume
347 
348  G4ThreeVector dPoint = daughterSolid->GetPointOnSurface();
349 
350  // Transform the generated point to the mother's coordinate system
351  // and finally to current volume's coordinate system
352  //
353  G4ThreeVector mp2 = Td.TransformPoint(dPoint);
354  G4ThreeVector msi = Tm.Inverse().TransformPoint(mp2);
355 
356  if (solid->Inside(msi)==kInside)
357  {
358  trials++; retval = true;
359  std::ostringstream message;
360  message << "Overlap with volume already placed !" << G4endl
361  << " Overlap is detected for volume "
362  << GetName() << G4endl
363  << " apparently fully encapsulating volume "
364  << daughter->GetName() << G4endl
365  << " at the same level !";
366  G4Exception("G4PVPlacement::CheckOverlaps()",
367  "GeomVol1002", JustWarning, message);
368  if (trials>=maxErr) { return true; }
369  }
370  }
371  }
372  }
373 
374  if (verbose)
375  {
376  G4cout << "OK! " << G4endl;
377  }
378 
379  return retval;
380 }
381 
382 // ----------------------------------------------------------------------
383 // NewPtrRotMatrix
384 //
385 // Auxiliary function for 2nd & 4th constructors (those with G4Transform3D)
386 // Creates a new rotation matrix on the heap (using "new") and copies its
387 // argument into it.
388 //
389 // NOTE: Ownership of the returned pointer is left to the caller !
390 // No entity is currently responsible to delete this memory.
391 //
394 {
395  G4RotationMatrix *pRotMatrix;
396  if ( RotMat.isIdentity() )
397  {
398  pRotMatrix = 0;
399  }
400  else
401  {
402  pRotMatrix = new G4RotationMatrix(RotMat);
403  }
404  // fallocatedRotM= ! (RotMat.isIdentity());
405 
406  return pRotMatrix;
407 }
const G4ThreeVector & GetTranslation() const
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
G4AffineTransform Inverse() const
G4bool IsReplicated() const
G4bool fallocatedRotM
G4VPhysicalVolume * GetDaughter(const G4int i) const
G4double a
Definition: TRTMaterials.hh:39
virtual ~G4PVPlacement()
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4bool IsRegularStructure() const
int G4int
Definition: G4Types.hh:78
void SetRotation(G4RotationMatrix *)
G4VPVParameterisation * GetParameterisation() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
bool G4bool
Definition: G4Types.hh:79
G4PVPlacement(G4RotationMatrix *pRot, const G4ThreeVector &tlate, G4LogicalVolume *pCurrentLogical, const G4String &pName, G4LogicalVolume *pMotherLogical, G4bool pMany, G4int pCopyNo, G4bool pSurfChk=false)
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
HepGeom::Transform3D G4Transform3D
const G4int n
G4LogicalVolume * GetMotherLogical() const
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int GetRegularStructureId() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4bool CheckOverlaps(G4int res=1000, G4double tol=0., G4bool verbose=true, G4int maxErr=1)
G4LogicalVolume * GetLogicalVolume() const
G4bool IsMany() const
EAxis
Definition: geomdefs.hh:54
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
const G4RotationMatrix * GetRotation() const
void SetCopyNo(G4int CopyNo)
#define G4endl
Definition: G4ios.hh:61
void SetMotherLogical(G4LogicalVolume *pMother)
G4bool IsParameterised() const
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
void AddDaughter(G4VPhysicalVolume *p)
void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
G4VSolid * GetSolid() const
static G4RotationMatrix * NewPtrRotMatrix(const G4RotationMatrix &RotMat)