Geant4  10.02.p01
G4ReflectionFactory.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: G4ReflectionFactory.cc 66872 2013-01-15 01:25:57Z japost $
28 //
29 //
30 // Class G4ReflectionFactory Implementation
31 //
32 // Decomposition of a general transformation
33 // that can include reflection in a "reflection-free" transformation:
34 //
35 // x(inM') = TG*x(inM) TG - general transformation
36 // = T*(R*x(inM)) T - "reflection-free" transformation
37 // = T* x(inReflM)
38 //
39 // Daughters transformation:
40 // When a volume V containing daughter D with transformation TD
41 // is placed in mother M with a general tranformation TGV,
42 // the TGV is decomposed,
43 // new reflected volume ReflV containing a new daughter ReflD
44 // with reflected transformation ReflTD is created:
45 //
46 // x(inV) = TD * x(inD);
47 // x(inM) = TGV * x(inV)
48 // = TV * R * x(inV)
49 // = TV * R * TD * x(inD)
50 // = TV * R*TD*R-1 * R*x(inD)
51 // = TV * ReflTD * x(inReflD)
52 
53 // Author: Ivana Hrivnacova, 16.10.2001 (Ivana.Hrivnacova@cern.ch)
54 // --------------------------------------------------------------------
55 
56 #include "G4ReflectionFactory.hh"
57 #include "G4ReflectedSolid.hh"
58 #include "G4Region.hh"
59 #include "G4LogicalVolume.hh"
60 #include "G4PVPlacement.hh"
61 #include "G4PVReplica.hh"
62 #include "G4VPVDivisionFactory.hh"
63 #include "G4GeometryTolerance.hh"
64 
68 
69 //_____________________________________________________________________________
70 
72 {
73  // Static singleton access method.
74  // ---
75 
76  if (!fInstance) { fInstance = new G4ReflectionFactory(); }
77 
78  return fInstance;
79 }
80 
81 //_____________________________________________________________________________
82 
84  : fVerboseLevel(0),
85  fNameExtension(fDefaultNameExtension)
86 {
87  // Protected singleton constructor.
88  // ---
89 
90  fScalePrecision = 10.
92  fInstance = this;
93 }
94 
95 //_____________________________________________________________________________
96 
98 {
99  delete fInstance;
100 }
101 
102 //
103 // public methods
104 //
105 
106 //_____________________________________________________________________________
107 
110  const G4String& name,
111  G4LogicalVolume* LV,
112  G4LogicalVolume* motherLV,
113  G4bool isMany,
114  G4int copyNo,
115  G4bool surfCheck)
116 {
117  // Evaluates the passed transformation; if it contains reflection
118  // it performs its decomposition, creates new reflected solid and
119  // logical volume (or retrieves them from a map if the reflected
120  // objects were already created), transforms the daughters (if present)
121  // and place it in the given mother.
122  // The result is a pair of physical volumes;
123  // the second physical volume is a placement in a reflected mother
124  // - or 0 if mother LV was not reflected.
125  // ---
126 
127  if (fVerboseLevel>0)
128  {
129  G4cout << "Place " << name << " lv " << LV << " "
130  << LV->GetName() << G4endl;
131  }
132 
133  // decompose transformation
134  G4Scale3D scale;
135  G4Rotate3D rotation;
136  G4Translate3D translation;
137 
138  transform3D.getDecomposition(scale, rotation, translation);
139  G4Transform3D pureTransform3D = translation * rotation;
140 
141  //PrintTransform(transform3D);
142  //PrintTransform(pureTransform3D);
143 
144  // check that scale correspond to fScale
145  //
146  CheckScale(scale);
147 
148  //
149  // reflection IS NOT present in transform3D
150  //
151 
152  if (! IsReflection(scale))
153  {
154  if (fVerboseLevel>0)
155  G4cout << "Scale positive" << G4endl;
156 
157  G4VPhysicalVolume* pv1
158  = new G4PVPlacement(pureTransform3D, LV, name,
159  motherLV, isMany, copyNo, surfCheck);
160 
161  G4VPhysicalVolume* pv2 = 0;
162  if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
163  {
164  // if mother was reflected
165  // reflect this LV and place it in reflected mother
166 
167  pv2 = new G4PVPlacement(fScale * (pureTransform3D * fScale.inverse()),
168  ReflectLV(LV, surfCheck), name, reflMotherLV,
169  isMany, copyNo, surfCheck);
170  }
171 
172  return G4PhysicalVolumesPair(pv1, pv2);
173  }
174 
175  //
176  // reflection IS present in transform3D
177  //
178 
179  if (fVerboseLevel>0)
180  G4cout << "scale negative" << G4endl;
181 
182  G4VPhysicalVolume* pv1
183  = new G4PVPlacement(pureTransform3D, ReflectLV(LV, surfCheck), name,
184  motherLV, isMany, copyNo, surfCheck);
185 
186  G4VPhysicalVolume* pv2 = 0;
187  if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
188  {
189 
190  // if mother was reflected
191  // place the refLV consituent in reflected mother
192 
193  pv2 = new G4PVPlacement(fScale * (pureTransform3D * fScale.inverse()),
194  LV, name, reflMotherLV, isMany, copyNo, surfCheck);
195  }
196 
197  return G4PhysicalVolumesPair(pv1, pv2);
198 }
199 
200 
201 //_____________________________________________________________________________
202 
205  G4LogicalVolume* LV,
206  G4LogicalVolume* motherLV,
207  EAxis axis,
208  G4int nofReplicas,
209  G4double width,
210  G4double offset)
211 {
212  // Creates replica in given mother.
213  // The result is a pair of physical volumes;
214  // the second physical volume is a replica in a reflected mother
215  // - or 0 if mother LV was not reflected.
216  // ---
217 
218  if (fVerboseLevel>0)
219  {
220  G4cout << "Replicate " << name << " lv " << LV << " "
221  << LV->GetName() << G4endl;
222  }
223 
224  G4VPhysicalVolume* pv1
225  = new G4PVReplica(name, LV, motherLV, axis, nofReplicas, width, offset);
226 
227  G4VPhysicalVolume* pv2 = 0;
228  if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
229  {
230  // if mother was reflected
231  // reflect the LV and replicate it in reflected mother
232 
233  pv2 = new G4PVReplica(name, ReflectLV(LV), reflMotherLV,
234  axis, nofReplicas, width, offset);
235  }
236 
237  return G4PhysicalVolumesPair(pv1, pv2);
238 }
239 
240 //_____________________________________________________________________________
241 
244  G4LogicalVolume* LV,
245  G4LogicalVolume* motherLV,
246  EAxis axis,
247  G4int nofDivisions,
248  G4double width,
249  G4double offset)
250 {
251  // Creates division in the given mother.
252  // The result is a pair of physical volumes;
253  // the second physical volume is a division in a reflected mother
254  // or 0 if mother LV was not reflected.
255  // ---
256 
257  if (fVerboseLevel>0)
258  {
259  G4cout << "Divide " << name << " lv " << LV << " "
260  << LV->GetName() << G4endl;
261  }
262 
263  G4VPVDivisionFactory* divisionFactory = GetPVDivisionFactory();
264 
265  G4VPhysicalVolume* pv1 = divisionFactory
266  ->CreatePVDivision(name, LV, motherLV, axis, nofDivisions, width, offset);
267 
268  G4VPhysicalVolume* pv2 = 0;
269  if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
270  {
271  // if mother was reflected
272  // reflect the LV and replicate it in reflected mother
273 
274  pv2 = divisionFactory->CreatePVDivision(name, ReflectLV(LV), reflMotherLV,
275  axis, nofDivisions, width, offset);
276  }
277 
278  return G4PhysicalVolumesPair(pv1, pv2);
279 }
280 
281 
282 //_____________________________________________________________________________
283 
286  G4LogicalVolume* LV,
287  G4LogicalVolume* motherLV,
288  EAxis axis,
289  G4int nofDivisions,
290  G4double offset)
291 {
292  // Creates division in the given mother.
293  // The result is a pair of physical volumes;
294  // the second physical volume is a division in a reflected mother
295  // or 0 if mother LV was not reflected.
296  // ---
297 
298  if (fVerboseLevel>0)
299  {
300  G4cout << "Divide " << name << " lv " << LV << " "
301  << LV->GetName() << G4endl;
302  }
303 
304  G4VPVDivisionFactory* divisionFactory = GetPVDivisionFactory();
305 
306  G4VPhysicalVolume* pv1 = divisionFactory
307  ->CreatePVDivision(name, LV, motherLV, axis, nofDivisions, offset);
308 
309  G4VPhysicalVolume* pv2 = 0;
310  if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
311  {
312  // if mother was reflected
313  // reflect the LV and replicate it in reflected mother
314 
315  pv2 = divisionFactory->CreatePVDivision(name, ReflectLV(LV), reflMotherLV,
316  axis, nofDivisions, offset);
317  }
318 
319  return G4PhysicalVolumesPair(pv1, pv2);
320 }
321 
322 
323 //_____________________________________________________________________________
324 
327  G4LogicalVolume* LV,
328  G4LogicalVolume* motherLV,
329  EAxis axis,
330  G4double width,
331  G4double offset)
332 {
333  // Creates division in the given mother.
334  // The result is a pair of physical volumes;
335  // the second physical volume is a division in a reflected mother
336  // or 0 if mother LV was not reflected.
337  // ---
338 
339  if (fVerboseLevel>0)
340  {
341  G4cout << "Divide " << name << " lv " << LV << " "
342  << LV->GetName() << G4endl;
343  }
344 
345  G4VPVDivisionFactory* divisionFactory = GetPVDivisionFactory();
346 
347  G4VPhysicalVolume* pv1 = divisionFactory
348  -> CreatePVDivision(name, LV, motherLV, axis, width, offset);
349 
350  G4VPhysicalVolume* pv2 = 0;
351  if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
352  {
353  // if mother was reflected
354  // reflect the LV and replicate it in reflected mother
355 
356  pv2 = divisionFactory->CreatePVDivision(name, ReflectLV(LV), reflMotherLV,
357  axis, width, offset);
358  }
359 
360  return G4PhysicalVolumesPair(pv1, pv2);
361 }
362 
363 
364 //
365 // private methods
366 //
367 
368 //_____________________________________________________________________________
369 
371  G4bool surfCheck)
372 {
373  // Gets/creates the reflected solid and logical volume
374  // and copies + transforms LV daughters.
375  // ---
376 
377  G4LogicalVolume* refLV = GetReflectedLV(LV);
378 
379  if (!refLV)
380  {
381 
382  // create new (reflected) objects
383  //
384  refLV = CreateReflectedLV(LV);
385 
386  // process daughters
387  //
388  ReflectDaughters(LV, refLV, surfCheck);
389 
390  // check if to be set as root region
391  //
392  if (LV->IsRootRegion())
393  {
394  LV->GetRegion()->AddRootLogicalVolume(refLV);
395  }
396  }
397 
398  return refLV;
399 }
400 
401 //_____________________________________________________________________________
402 
404 {
405  // Creates the reflected solid and logical volume
406  // and add the logical volumes pair in the maps.
407  // ---
408 
409  // consistency check
410  //
411  if (fReflectedLVMap.find(LV) != fReflectedLVMap.end())
412  {
413  std::ostringstream message;
414  message << "Invalid reflection for volume: "
415  << LV->GetName() << G4endl
416  << "Cannot be applied to a volume already reflected !";
417  G4Exception("G4ReflectionFactory::CreateReflectedLV()",
418  "GeomVol0002", FatalException, message);
419  }
420 
421  G4VSolid* refSolid
423  LV->GetSolid(), fScale);
424 
425  G4LogicalVolume* refLV
426  = new G4LogicalVolume(refSolid,
427  LV->GetMaterial(),
428  LV->GetName() + fNameExtension,
429  LV->GetFieldManager(),
430  LV->GetSensitiveDetector(),
431  LV->GetUserLimits());
432  refLV->SetVisAttributes(LV->GetVisAttributes()); // vis-attributes
433  refLV->SetBiasWeight(LV->GetBiasWeight()); // biasing weight
434  if (LV->IsRegion())
435  {
436  refLV->SetRegion(LV->GetRegion()); // set a region in case
437  }
438 
439  fConstituentLVMap[LV] = refLV;
440  fReflectedLVMap[refLV] = LV;
441 
442  return refLV;
443 }
444 
445 //_____________________________________________________________________________
446 
448  G4LogicalVolume* refLV,
449  G4bool surfCheck)
450 {
451  // Reflects daughters recursively.
452  // ---
453 
454  if (fVerboseLevel>0)
455  {
456  G4cout << "G4ReflectionFactory::ReflectDaughters(): "
457  << LV->GetNoDaughters() << " of " << LV->GetName() << G4endl;
458  }
459 
460  for (G4int i=0; i<LV->GetNoDaughters(); i++)
461  {
462  G4VPhysicalVolume* dPV = LV->GetDaughter(i);
463 
464  if (! dPV->IsReplicated())
465  {
466  ReflectPVPlacement(dPV, refLV, surfCheck);
467  }
468  else if (! dPV->GetParameterisation())
469  {
470  ReflectPVReplica(dPV, refLV);
471  }
472  else if (G4VPVDivisionFactory::Instance() &&
473  G4VPVDivisionFactory::Instance()->IsPVDivision(dPV))
474  {
475  ReflectPVDivision(dPV, refLV);
476  }
477  else
478  {
479  ReflectPVParameterised(dPV, refLV, surfCheck);
480  }
481  }
482 }
483 
484 //_____________________________________________________________________________
485 
487  G4LogicalVolume* refLV,
488  G4bool surfCheck)
489 {
490  // Copies and transforms daughter of PVPlacement type of
491  // a constituent volume into a reflected volume.
492  // ---
493 
494  G4LogicalVolume* dLV = dPV->GetLogicalVolume();
495 
496  // update daughter transformation
497  //
499  dt = fScale * (dt * fScale.inverse());
500 
501  G4LogicalVolume* refDLV;
502 
503  if (fVerboseLevel>0)
504  G4cout << "Daughter: " << dPV << " " << dLV->GetName();
505 
506  if (!IsReflected(dLV))
507  {
508 
509  if (fVerboseLevel>0)
510  G4cout << " will be reflected." << G4endl;
511 
512  // get reflected volume if already created //
513  refDLV = GetReflectedLV(dLV);
514 
515  if (!refDLV)
516  {
517  // create new daughter solid and logical volume
518  //
519  refDLV = CreateReflectedLV(dLV);
520 
521  // recursive call
522  //
523  ReflectDaughters(dLV, refDLV, surfCheck);
524  }
525 
526  // create new daughter physical volume
527  // with updated transformation
528 
529  new G4PVPlacement(dt, refDLV, dPV->GetName(), refLV,
530  dPV->IsMany(), dPV->GetCopyNo(), surfCheck);
531 
532  }
533  else
534  {
535  if (fVerboseLevel>0)
536  G4cout << " will be reconstitued." << G4endl;
537 
538  refDLV = GetConstituentLV(dLV);
539 
540  new G4PVPlacement(dt, refDLV, dPV->GetName(), refLV,
541  dPV->IsMany(), dPV->GetCopyNo(), surfCheck);
542  }
543 }
544 
545 //_____________________________________________________________________________
546 
548  G4LogicalVolume* refLV)
549 {
550  // Copies and transforms daughter of PVReplica type of
551  // a constituent volume into a reflected volume.
552  // ---
553 
554  G4LogicalVolume* dLV = dPV->GetLogicalVolume();
555 
556  // get replication data
557  //
558  EAxis axis;
559  G4int nofReplicas;
560  G4double width;
561  G4double offset;
562  G4bool consuming;
563 
564  dPV->GetReplicationData(axis, nofReplicas, width, offset, consuming);
565 
566  G4LogicalVolume* refDLV;
567 
568  if (fVerboseLevel>0)
569  G4cout << "Daughter: " << dPV << " " << dLV->GetName();
570 
571  if (!IsReflected(dLV))
572  {
573  if (fVerboseLevel>0)
574  G4cout << " will be reflected." << G4endl;
575 
576  // get reflected volume if already created
577  //
578  refDLV = GetReflectedLV(dLV);
579 
580  if (!refDLV)
581  {
582  // create new daughter solid and logical volume
583  //
584  refDLV = CreateReflectedLV(dLV);
585 
586  // recursive call
587  //
588  ReflectDaughters(dLV, refDLV);
589  }
590 
591  // create new daughter replica
592  //
593  new G4PVReplica(dPV->GetName(), refDLV, refLV,
594  axis, nofReplicas, width, offset);
595  }
596  else
597  {
598  if (fVerboseLevel>0)
599  G4cout << " will be reconstitued." << G4endl;
600 
601  refDLV = GetConstituentLV(dLV);
602 
603  new G4PVReplica(dPV->GetName(), refDLV, refLV,
604  axis, nofReplicas, width, offset);
605  }
606 }
607 
608 //_____________________________________________________________________________
609 
611  G4LogicalVolume* refLV)
612 {
613  // Copies and transforms daughter of PVDivision type of
614  // a constituent volume into a reflected volume.
615  // ---
616 
617  G4VPVDivisionFactory* divisionFactory = GetPVDivisionFactory();
618 
619  G4LogicalVolume* dLV = dPV->GetLogicalVolume();
620 
621  // get parameterisation data
622  //
624 
625  G4LogicalVolume* refDLV;
626 
627  if (fVerboseLevel>0)
628  G4cout << "Daughter: " << dPV << " " << dLV->GetName();
629 
630  if (!IsReflected(dLV))
631  {
632  if (fVerboseLevel>0)
633  G4cout << " will be reflected." << G4endl;
634 
635  // get reflected volume if already created
636  //
637  refDLV = GetReflectedLV(dLV);
638 
639  if (!refDLV)
640  {
641  // create new daughter solid and logical volume
642  //
643  refDLV = CreateReflectedLV(dLV);
644 
645  // recursive call
646  //
647  ReflectDaughters(dLV, refDLV);
648  }
649 
650  // create new daughter replica
651  //
652  divisionFactory->CreatePVDivision(dPV->GetName(), refDLV, refLV, param);
653  }
654  else
655  {
656  if (fVerboseLevel>0)
657  G4cout << " will be reconstitued." << G4endl;
658 
659  refDLV = GetConstituentLV(dLV);
660 
661  divisionFactory->CreatePVDivision(dPV->GetName(), refDLV, refLV, param);
662  }
663 }
664 
665 //_____________________________________________________________________________
666 
669 {
670  // Not implemented.
671  // Should copy and transform daughter of PVReplica type of
672  // a constituent volume into a reflected volume.
673  // ---
674 
675  std::ostringstream message;
676  message << "Not yet implemented. Volume: " << dPV->GetName() << G4endl
677  << "Reflection of parameterised volumes is not yet implemented.";
678  G4Exception("G4ReflectionFactory::ReflectPVParameterised()",
679  "GeomVol0001", FatalException, message);
680 }
681 
682 //_____________________________________________________________________________
683 
686 {
687  // Returns the consituent volume of the given reflected volume,
688  // 0 if the given reflected volume was not found.
689  // ---
690 
691  LogicalVolumesMapIterator it = fReflectedLVMap.find(reflLV);
692 
693  if (it == fReflectedLVMap.end()) return 0;
694 
695  return (*it).second;
696 }
697 
698 //_____________________________________________________________________________
699 
702 {
703  // Returns the reflected volume of the given consituent volume,
704  // 0 if the given volume was not reflected.
705  // ---
706 
708 
709  if (it == fConstituentLVMap.end()) return 0;
710 
711  return (*it).second;
712 }
713 
714 //_____________________________________________________________________________
715 
717 {
718  // Returns true if the given volume has been already reflected
719  // (is in the map of constituent volumes).
720  // ---
721 
722  return (fConstituentLVMap.find(lv) != fConstituentLVMap.end());
723 }
724 
725 //_____________________________________________________________________________
726 
728 {
729  // Returns true if the given volume is a reflected volume
730  // (is in the map reflected volumes).
731  // ---
732 
733  return (fReflectedLVMap.find(lv) != fReflectedLVMap.end());
734 }
735 
736 //_____________________________________________________________________________
737 
739 {
740  // Returns true if the scale is negative, false otherwise.
741  // ---
742 
743  if (scale(0,0)*scale(1,1)*scale(2,2) < 0.)
744  return true;
745  else
746  return false;
747 }
748 
749 //_____________________________________________________________________________
750 
753 {
754  return fReflectedLVMap;
755 }
756 
757 //_____________________________________________________________________________
758 
759 void
761 {
762  fConstituentLVMap.~map();
763  fReflectedLVMap.~map();
764 }
765 
766 //_____________________________________________________________________________
767 
769 {
770  // temporary - for debugging purpose
771  // ---
772 
774  for (it = fConstituentLVMap.begin(); it != fConstituentLVMap.end(); it++)
775  {
776  G4cout << "lv: " << (*it).first << " lv_refl: " << (*it).second << G4endl;
777  }
778  G4cout << G4endl;
779 }
780 
781 //_____________________________________________________________________________
782 
784 {
785  // Check if scale correspond to fScale,
786  // if not give exception.
787  // ---
788 
789  if (!IsReflection(scale)) return;
790 
791  G4double diff = 0.;
792  for (G4int i=0; i<4; i++)
793  for (G4int j=0; j<4; j++)
794  diff += std::abs(scale(i,j) - fScale(i,j));
795 
796  if (diff > fScalePrecision)
797  {
798  std::ostringstream message;
799  message << "Unexpected scale in input !" << G4endl
800  << " Difference: " << diff;
801  G4Exception("G4ReflectionFactory::CheckScale()",
802  "GeomVol0002", FatalException, message);
803  }
804 }
805 
806 //_____________________________________________________________________________
807 
809 {
810  // Returns the G4PVDivisionFactory instance if it exists,
811  // otherwise gives exception
812  // ---
813 
815  if (!divisionFactory)
816  {
817  std::ostringstream message;
818  message << "A concrete G4PVDivisionFactory instantiated is required !"
819  << G4endl
820  << " It has been requested to reflect divided volumes."
821  << G4endl
822  << " In this case, it is required to instantiate a concrete"
823  << G4endl
824  << " factory G4PVDivisionFactory in your program -before-"
825  << G4endl
826  << " executing the reflection !";
827  G4Exception("G4ReflectionFactory::GetPVDivisionFactory()",
828  "GeomVol0002", FatalException, message);
829  }
830 
831  return divisionFactory;
832 }
833 
834 //_____________________________________________________________________________
835 
837 {
838  fScalePrecision = scaleValue;
839 }
840 
841 //_____________________________________________________________________________
842 
844 {
845  return fScalePrecision;
846 }
847 
848 //_____________________________________________________________________________
849 
851 {
852  fVerboseLevel = verboseLevel;
853 }
854 
855 //_____________________________________________________________________________
856 
858 {
859  return fVerboseLevel;
860 }
861 
862 //_____________________________________________________________________________
863 
865 {
866  fNameExtension = nameExtension;
867 }
868 
869 //_____________________________________________________________________________
870 
872 {
873  return fNameExtension;
874 }
875 
876 /*
877  // placement with decomposed transformation
878 
879  G4VPhysicalVolume* pv1
880  = new G4PVPlacement(new G4RotationMatrix(rotation.getRotation().inverse()),
881  translation.getTranslation(),
882  refLV, name, motherLV, isMany, copyNo);
883 */
G4PhysicalVolumesPair Divide(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofDivisions, G4double width, G4double offset)
G4String GetName() const
G4ReflectedVolumesMap fReflectedLVMap
void ReflectPVParameterised(G4VPhysicalVolume *PV, G4LogicalVolume *refLV, G4bool surfCheck=false)
const G4String & GetVolumesNameExtension() const
void AddRootLogicalVolume(G4LogicalVolume *lv)
Definition: G4Region.cc:254
virtual G4bool IsReplicated() const =0
G4Material * GetMaterial() const
G4String name
Definition: TRTMaterials.hh:40
void SetScalePrecision(G4double scaleValue)
G4double GetSurfaceTolerance() const
G4VPhysicalVolume * GetDaughter(const G4int i) const
G4LogicalVolume * GetConstituentLV(G4LogicalVolume *reflLV) const
#define width
G4PhysicalVolumesPair Replicate(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofReplicas, G4double width, G4double offset=0)
G4Region * GetRegion() const
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
void CheckScale(const G4Scale3D &scale) const
G4bool IsReflected(G4LogicalVolume *lv) const
static G4ReflectionFactory * Instance()
static G4VPVDivisionFactory * Instance()
G4double GetBiasWeight() const
G4FieldManager * GetFieldManager() const
G4GLOB_DLL std::ostream G4cout
G4PhysicalVolumesPair Place(const G4Transform3D &transform3D, const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, G4bool isMany, G4int copyNo, G4bool surfCheck=false)
G4ReflectedVolumesMap fConstituentLVMap
const G4String & GetName() const
virtual G4bool IsMany() const =0
G4ReflectedVolumesMap::const_iterator LogicalVolumesMapIterator
const G4ReflectedVolumesMap & GetReflectedVolumesMap() const
bool G4bool
Definition: G4Types.hh:79
G4bool IsRootRegion() const
G4VPVDivisionFactory * GetPVDivisionFactory() const
virtual G4VPVParameterisation * GetParameterisation() const =0
G4RotationMatrix GetObjectRotationValue() const
std::map< G4LogicalVolume *, G4LogicalVolume *, std::less< G4LogicalVolume * > > G4ReflectedVolumesMap
G4bool IsRegion() const
const G4VisAttributes * GetVisAttributes() const
static const G4Scale3D fScale
void SetVolumesNameExtension(const G4String &nameExtension)
HepGeom::Transform3D G4Transform3D
HepGeom::Scale3D G4Scale3D
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4UserLimits * GetUserLimits() const
static const G4String fDefaultNameExtension
void SetVerboseLevel(G4int verboseLevel)
virtual G4VPhysicalVolume * CreatePVDivision(const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMother, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double offset)=0
HepGeom::Rotate3D G4Rotate3D
G4double GetScalePrecision() const
G4LogicalVolume * GetLogicalVolume() const
G4int GetVerboseLevel() const
EAxis
Definition: geomdefs.hh:54
G4bool IsReflection(const G4Scale3D &scale) const
G4LogicalVolume * GetReflectedLV(G4LogicalVolume *lv) const
virtual G4int GetCopyNo() const =0
G4LogicalVolume * ReflectLV(G4LogicalVolume *LV, G4bool surfCheck=false)
HepGeom::ScaleZ3D G4ScaleZ3D
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
HepGeom::Translate3D G4Translate3D
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
#define G4endl
Definition: G4ios.hh:61
static G4ThreadLocal G4ReflectionFactory * fInstance
void ReflectPVPlacement(G4VPhysicalVolume *PV, G4LogicalVolume *refLV, G4bool surfCheck=false)
G4ThreeVector GetObjectTranslation() const
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
void ReflectDaughters(G4LogicalVolume *LV, G4LogicalVolume *refLV, G4bool surfCheck=false)
void ReflectPVDivision(G4VPhysicalVolume *PV, G4LogicalVolume *refLV)
G4VSensitiveDetector * GetSensitiveDetector() const
G4LogicalVolume * CreateReflectedLV(G4LogicalVolume *LV)
G4bool IsConstituent(G4LogicalVolume *lv) const
void SetVisAttributes(const G4VisAttributes *pVA)
static G4GeometryTolerance * GetInstance()
G4VSolid * GetSolid() const
void ReflectPVReplica(G4VPhysicalVolume *PV, G4LogicalVolume *refLV)