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