Geant4 11.4.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 (IN2P3/IJCLab Orsay), 16 October 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) { G4cout << "Scale positive" << G4endl; }
150
152 = new G4PVPlacement(pureTransform3D, LV, name,
153 motherLV, isMany, copyNo, surfCheck);
154
155 G4VPhysicalVolume* pv2 = nullptr;
156 if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
157 {
158 // if mother was reflected
159 // reflect this LV and place it in reflected mother
160
161 pv2 = new G4PVPlacement(fScale * (pureTransform3D * fScale.inverse()),
162 ReflectLV(LV, surfCheck), name, reflMotherLV,
163 isMany, copyNo, surfCheck);
164 }
165
166 return {pv1, pv2};
167 }
168
169 //
170 // reflection IS present in transform3D
171 //
172
173 if (fVerboseLevel>0) { G4cout << "scale negative" << G4endl; }
174
176 = new G4PVPlacement(pureTransform3D, ReflectLV(LV, surfCheck), name,
177 motherLV, isMany, copyNo, surfCheck);
178
179 G4VPhysicalVolume* pv2 = nullptr;
180 if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
181 {
182
183 // if mother was reflected
184 // place the refLV consituent in reflected mother
185
186 pv2 = new G4PVPlacement(fScale * (pureTransform3D * fScale.inverse()),
187 LV, name, reflMotherLV, isMany, copyNo, surfCheck);
188 }
189
190 return {pv1, pv2};
191}
192
193
194//_____________________________________________________________________________
195
198 G4LogicalVolume* LV,
199 G4LogicalVolume* motherLV,
200 EAxis axis,
201 G4int nofReplicas,
202 G4double width,
204{
205 // Creates replica in given mother.
206 // The result is a pair of physical volumes;
207 // the second physical volume is a replica in a reflected mother
208 // - or nullptr if mother LV was not reflected.
209 // ---
210
211 if (fVerboseLevel>0)
212 {
213 G4cout << "Replicate " << name << " lv " << LV << " "
214 << LV->GetName() << G4endl;
215 }
216
218 = new G4PVReplica(name, LV, motherLV, axis, nofReplicas, width, offset);
219
220 G4VPhysicalVolume* pv2 = nullptr;
221 if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
222 {
223 // if mother was reflected
224 // reflect the LV and replicate it in reflected mother
225
226 pv2 = new G4PVReplica(name, ReflectLV(LV), reflMotherLV,
227 axis, nofReplicas, width, offset);
228 }
229
230 return {pv1, pv2};
231}
232
233//_____________________________________________________________________________
234
237 G4LogicalVolume* LV,
238 G4LogicalVolume* motherLV,
239 EAxis axis,
240 G4int nofDivisions,
241 G4double width,
243{
244 // Creates division in the given mother.
245 // The result is a pair of physical volumes;
246 // the second physical volume is a division in a reflected mother
247 // or nullptr if mother LV was not reflected.
248 // ---
249
250 if (fVerboseLevel>0)
251 {
252 G4cout << "Divide " << name << " lv " << LV << " "
253 << LV->GetName() << G4endl;
254 }
255
256 G4VPVDivisionFactory* divisionFactory = GetPVDivisionFactory();
257
258 G4VPhysicalVolume* pv1 = divisionFactory
259 ->CreatePVDivision(name, LV, motherLV, axis, nofDivisions, width, offset);
260
261 G4VPhysicalVolume* pv2 = nullptr;
262 if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
263 {
264 // if mother was reflected
265 // reflect the LV and replicate it in reflected mother
266
267 pv2 = divisionFactory->CreatePVDivision(name, ReflectLV(LV), reflMotherLV,
268 axis, nofDivisions, width, offset);
269 }
270
271 return {pv1, pv2};
272}
273
274
275//_____________________________________________________________________________
276
279 G4LogicalVolume* LV,
280 G4LogicalVolume* motherLV,
281 EAxis axis,
282 G4int nofDivisions,
284{
285 // Creates division in the given mother.
286 // The result is a pair of physical volumes;
287 // the second physical volume is a division in a reflected mother
288 // or nullptr if mother LV was not reflected.
289 // ---
290
291 if (fVerboseLevel>0)
292 {
293 G4cout << "Divide " << name << " lv " << LV << " "
294 << LV->GetName() << G4endl;
295 }
296
297 G4VPVDivisionFactory* divisionFactory = GetPVDivisionFactory();
298
299 G4VPhysicalVolume* pv1 = divisionFactory
300 ->CreatePVDivision(name, LV, motherLV, axis, nofDivisions, offset);
301
302 G4VPhysicalVolume* pv2 = nullptr;
303 if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
304 {
305 // if mother was reflected
306 // reflect the LV and replicate it in reflected mother
307
308 pv2 = divisionFactory->CreatePVDivision(name, ReflectLV(LV), reflMotherLV,
309 axis, nofDivisions, offset);
310 }
311
312 return {pv1, pv2};
313}
314
315
316//_____________________________________________________________________________
317
320 G4LogicalVolume* LV,
321 G4LogicalVolume* motherLV,
322 EAxis axis,
323 G4double width,
325{
326 // Creates division in the given mother.
327 // The result is a pair of physical volumes;
328 // the second physical volume is a division in a reflected mother
329 // or nullptr if mother LV was not reflected.
330 // ---
331
332 if (fVerboseLevel>0)
333 {
334 G4cout << "Divide " << name << " lv " << LV << " "
335 << LV->GetName() << G4endl;
336 }
337
338 G4VPVDivisionFactory* divisionFactory = GetPVDivisionFactory();
339
340 G4VPhysicalVolume* pv1 = divisionFactory
341 -> CreatePVDivision(name, LV, motherLV, axis, width, offset);
342
343 G4VPhysicalVolume* pv2 = nullptr;
344 if (G4LogicalVolume* reflMotherLV = GetReflectedLV(motherLV))
345 {
346 // if mother was reflected
347 // reflect the LV and replicate it in reflected mother
348
349 pv2 = divisionFactory->CreatePVDivision(name, ReflectLV(LV), reflMotherLV,
350 axis, width, offset);
351 }
352
353 return {pv1, pv2};
354}
355
356
357//
358// private methods
359//
360
361//_____________________________________________________________________________
362
363G4LogicalVolume* G4ReflectionFactory::ReflectLV(G4LogicalVolume* LV,
364 G4bool surfCheck)
365{
366 // Gets/creates the reflected solid and logical volume
367 // and copies + transforms LV daughters.
368 // ---
369
370 G4LogicalVolume* refLV = GetReflectedLV(LV);
371
372 if (refLV == nullptr)
373 {
374
375 // create new (reflected) objects
376 //
377 refLV = CreateReflectedLV(LV);
378
379 // process daughters
380 //
381 ReflectDaughters(LV, refLV, surfCheck);
382
383 // check if to be set as root region
384 //
385 if (LV->IsRootRegion())
386 {
387 LV->GetRegion()->AddRootLogicalVolume(refLV);
388 }
389 }
390
391 return refLV;
392}
393
394//_____________________________________________________________________________
395
396G4LogicalVolume* G4ReflectionFactory::CreateReflectedLV(G4LogicalVolume* LV)
397{
398 // Creates the reflected solid and logical volume
399 // and add the logical volumes pair in the maps.
400 // ---
401
402 // consistency check
403 //
404 if (fReflectedLVMap.find(LV) != fReflectedLVMap.cend())
405 {
406 std::ostringstream message;
407 message << "Invalid reflection for volume: "
408 << LV->GetName() << G4endl
409 << "Cannot be applied to a volume already reflected !";
410 G4Exception("G4ReflectionFactory::CreateReflectedLV()",
411 "GeomVol0002", FatalException, message);
412 }
413
414 G4VSolid* refSolid
415 = new G4ReflectedSolid(LV->GetSolid()->GetName() + fNameExtension,
416 LV->GetSolid(), fScale);
417
418 auto refLV
419 = new G4LogicalVolume(refSolid,
420 LV->GetMaterial(),
421 LV->GetName() + fNameExtension,
422 LV->GetFieldManager(),
424 LV->GetUserLimits());
425
426 if (LV->GetVisAttributes() != nullptr) {
427 refLV->SetVisAttributes(*LV->GetVisAttributes()); // vis-attributes
428 }
429 refLV->SetBiasWeight(LV->GetBiasWeight()); // biasing weight
430 if (LV->IsRegion())
431 {
432 refLV->SetRegion(LV->GetRegion()); // set a region in case
433 }
434
435 fConstituentLVMap[LV] = refLV;
436 fReflectedLVMap[refLV] = LV;
437
438 return refLV;
439}
440
441//_____________________________________________________________________________
442
443void G4ReflectionFactory::ReflectDaughters(G4LogicalVolume* LV,
444 G4LogicalVolume* refLV,
445 G4bool surfCheck)
446{
447 // Reflects daughters recursively.
448 // ---
449
450 if (fVerboseLevel>0)
451 {
452 G4cout << "G4ReflectionFactory::ReflectDaughters(): "
453 << LV->GetNoDaughters() << " of " << LV->GetName() << G4endl;
454 }
455
456 for (std::size_t i=0; i<LV->GetNoDaughters(); ++i)
457 {
458 G4VPhysicalVolume* dPV = LV->GetDaughter((G4int)i);
459
460 if (!dPV->IsReplicated())
461 {
462 ReflectPVPlacement(dPV, refLV, surfCheck);
463 }
464 else if (dPV->GetParameterisation() == nullptr)
465 {
466 ReflectPVReplica(dPV, refLV);
467 }
468 else if ((G4VPVDivisionFactory::Instance() != nullptr) &&
469 G4VPVDivisionFactory::Instance()->IsPVDivision(dPV))
470 {
471 ReflectPVDivision(dPV, refLV);
472 }
473 else
474 {
475 ReflectPVParameterised(dPV, refLV, surfCheck);
476 }
477 }
478}
479
480//_____________________________________________________________________________
481
482void G4ReflectionFactory::ReflectPVPlacement(G4VPhysicalVolume* dPV,
483 G4LogicalVolume* refLV,
484 G4bool surfCheck)
485{
486 // Copies and transforms daughter of PVPlacement type of
487 // a constituent volume into a reflected volume.
488 // ---
489
490 G4LogicalVolume* dLV = dPV->GetLogicalVolume();
491
492 // update daughter transformation
493 //
495 dt = fScale * (dt * fScale.inverse());
496
497 G4LogicalVolume* refDLV;
498
499 if (fVerboseLevel>0) { G4cout << "Daughter: " << dPV << " " << dLV->GetName(); }
500
501 if (!IsReflected(dLV))
502 {
503
504 if (fVerboseLevel>0) { G4cout << " will be reflected." << G4endl; }
505
506 // get reflected volume if already created
507 refDLV = GetReflectedLV(dLV);
508
509 if (refDLV == nullptr)
510 {
511 // create new daughter solid and logical volume
512 //
513 refDLV = CreateReflectedLV(dLV);
514
515 // recursive call
516 //
517 ReflectDaughters(dLV, refDLV, surfCheck);
518 }
519
520 // create new daughter physical volume
521 // with updated transformation
522
523 new G4PVPlacement(dt, refDLV, dPV->GetName(), refLV,
524 dPV->IsMany(), dPV->GetCopyNo(), surfCheck);
525
526 }
527 else
528 {
529 if (fVerboseLevel>0) { G4cout << " will be reconstitued." << G4endl; }
530
531 refDLV = GetConstituentLV(dLV);
532
533 new G4PVPlacement(dt, refDLV, dPV->GetName(), refLV,
534 dPV->IsMany(), dPV->GetCopyNo(), surfCheck);
535 }
536}
537
538//_____________________________________________________________________________
539
540void G4ReflectionFactory::ReflectPVReplica(G4VPhysicalVolume* dPV,
541 G4LogicalVolume* refLV)
542{
543 // Copies and transforms daughter of PVReplica type of
544 // a constituent volume into a reflected volume.
545 // ---
546
547 G4LogicalVolume* dLV = dPV->GetLogicalVolume();
548
549 // get replication data
550 //
551 EAxis axis;
552 G4int nofReplicas;
553 G4double width;
555 G4bool consuming;
556
557 dPV->GetReplicationData(axis, nofReplicas, width, offset, consuming);
558
559 G4LogicalVolume* refDLV;
560
561 if (fVerboseLevel>0) { G4cout << "Daughter: " << dPV << " " << dLV->GetName(); }
562
563 if (!IsReflected(dLV))
564 {
565 if (fVerboseLevel>0) { G4cout << " will be reflected." << G4endl; }
566
567 // get reflected volume if already created
568 //
569 refDLV = GetReflectedLV(dLV);
570
571 if (refDLV == nullptr)
572 {
573 // create new daughter solid and logical volume
574 //
575 refDLV = CreateReflectedLV(dLV);
576
577 // recursive call
578 //
579 ReflectDaughters(dLV, refDLV);
580 }
581
582 // create new daughter replica
583 //
584 new G4PVReplica(dPV->GetName(), refDLV, refLV,
585 axis, nofReplicas, width, offset);
586 }
587 else
588 {
589 if (fVerboseLevel>0) { G4cout << " will be reconstitued." << G4endl; }
590
591 refDLV = GetConstituentLV(dLV);
592
593 new G4PVReplica(dPV->GetName(), refDLV, refLV,
594 axis, nofReplicas, width, offset);
595 }
596}
597
598//_____________________________________________________________________________
599
600void G4ReflectionFactory::ReflectPVDivision(G4VPhysicalVolume* dPV,
601 G4LogicalVolume* refLV)
602{
603 // Copies and transforms daughter of PVDivision type of
604 // a constituent volume into a reflected volume.
605 // ---
606
607 G4VPVDivisionFactory* divisionFactory = GetPVDivisionFactory();
608
609 G4LogicalVolume* dLV = dPV->GetLogicalVolume();
610
611 // get parameterisation data
612 //
613 G4VPVParameterisation* param = dPV->GetParameterisation();
614
615 G4LogicalVolume* refDLV;
616
617 if (fVerboseLevel>0) { G4cout << "Daughter: " << dPV << " " << dLV->GetName(); }
618
619 if (!IsReflected(dLV))
620 {
621 if (fVerboseLevel>0) { G4cout << " will be reflected." << G4endl; }
622
623 // get reflected volume if already created
624 //
625 refDLV = GetReflectedLV(dLV);
626
627 if (refDLV == nullptr)
628 {
629 // create new daughter solid and logical volume
630 //
631 refDLV = CreateReflectedLV(dLV);
632
633 // recursive call
634 //
635 ReflectDaughters(dLV, refDLV);
636 }
637
638 // create new daughter replica
639 //
640 divisionFactory->CreatePVDivision(dPV->GetName(), refDLV, refLV, param);
641 }
642 else
643 {
644 if (fVerboseLevel>0) { G4cout << " will be reconstitued." << G4endl; }
645
646 refDLV = GetConstituentLV(dLV);
647
648 divisionFactory->CreatePVDivision(dPV->GetName(), refDLV, refLV, param);
649 }
650}
651
652//_____________________________________________________________________________
653
654void G4ReflectionFactory::ReflectPVParameterised(G4VPhysicalVolume* dPV,
656{
657 // Not implemented.
658 // Should copy and transform daughter of G4PVParameterised type of
659 // a constituent volume into a reflected volume.
660 // ---
661
662 std::ostringstream message;
663 message << "Not yet implemented. Volume: " << dPV->GetName() << G4endl
664 << "Reflection of parameterised volumes is not yet implemented.";
665 G4Exception("G4ReflectionFactory::ReflectPVParameterised()",
666 "GeomVol0001", FatalException, message);
667}
668
669//_____________________________________________________________________________
670
673{
674 // Returns the consituent volume of the given reflected volume,
675 // nullptr if the given reflected volume was not found.
676 // ---
677
678 auto it = fReflectedLVMap.find(reflLV);
679
680 if (it == fReflectedLVMap.cend()) { return nullptr; }
681
682 return (*it).second;
683}
684
685//_____________________________________________________________________________
686
689{
690 // Returns the reflected volume of the given consituent volume,
691 // nullptr if the given volume was not reflected.
692 // ---
693
694 auto it = fConstituentLVMap.find(lv);
695
696 if (it == fConstituentLVMap.cend()) { return nullptr; }
697
698 return (*it).second;
699}
700
701//_____________________________________________________________________________
702
704{
705 // Returns true if the given volume has been already reflected
706 // (is in the map of constituent volumes).
707 // ---
708
709 return (fConstituentLVMap.find(lv) != fConstituentLVMap.cend());
710}
711
712//_____________________________________________________________________________
713
715{
716 // Returns true if the given volume is a reflected volume
717 // (is in the map reflected volumes).
718 // ---
719
720 return (fReflectedLVMap.find(lv) != fReflectedLVMap.cend());
721}
722
723//_____________________________________________________________________________
724
725G4bool G4ReflectionFactory::IsReflection(const G4Scale3D& scale) const
726{
727 // Returns true if the scale is negative, false otherwise.
728 // ---
729
730 return scale(0,0)*scale(1,1)*scale(2,2) < 0.;
731}
732
733//_____________________________________________________________________________
734
737{
738 return fReflectedLVMap;
739}
740
741//_____________________________________________________________________________
742
744{
745 fConstituentLVMap.clear();
746 fReflectedLVMap.clear();
747}
748
749//_____________________________________________________________________________
750
751void G4ReflectionFactory::PrintConstituentLVMap() const
752{
753 // temporary - for debugging purpose
754 // ---
755
756 for (const auto & it : fConstituentLVMap)
757 {
758 G4cout << "lv: " << it.first << " lv_refl: " << it.second << G4endl;
759 }
760 G4cout << G4endl;
761}
762
763//_____________________________________________________________________________
764
765void G4ReflectionFactory::CheckScale(const G4Scale3D& scale) const
766{
767 // Check if scale correspond to fScale,
768 // if not give exception.
769 // ---
770
771 if (!IsReflection(scale)) { return;
772}
773
774 G4double diff = 0.;
775 for (auto i=0; i<4; ++i)
776 {
777 for (auto j=0; j<4; ++j)
778 {
779 diff += std::abs(scale(i,j) - fScale(i,j));
780 }
781 }
782
783 if (diff > fScalePrecision)
784 {
785 std::ostringstream message;
786 message << "Unexpected scale in input !" << G4endl
787 << " Difference: " << diff;
788 G4Exception("G4ReflectionFactory::CheckScale()",
789 "GeomVol0002", FatalException, message);
790 }
791}
792
793//_____________________________________________________________________________
794
795G4VPVDivisionFactory* G4ReflectionFactory::GetPVDivisionFactory() const
796{
797 // Returns the G4PVDivisionFactory instance if it exists,
798 // otherwise gives exception
799 // ---
800
801 G4VPVDivisionFactory* divisionFactory = G4VPVDivisionFactory::Instance();
802 if (divisionFactory == nullptr)
803 {
804 std::ostringstream message;
805 message << "A concrete G4PVDivisionFactory instantiated is required !"
806 << G4endl
807 << " It has been requested to reflect divided volumes."
808 << G4endl
809 << " In this case, it is required to instantiate a concrete"
810 << G4endl
811 << " factory G4PVDivisionFactory in your program -before-"
812 << G4endl
813 << " executing the reflection !";
814 G4Exception("G4ReflectionFactory::GetPVDivisionFactory()",
815 "GeomVol0002", FatalException, message);
816 }
817
818 return divisionFactory;
819}
820
821//_____________________________________________________________________________
822
824{
825 fScalePrecision = scaleValue;
826}
827
828//_____________________________________________________________________________
829
831{
832 return fScalePrecision;
833}
834
835//_____________________________________________________________________________
836
838{
839 fVerboseLevel = verboseLevel;
840}
841
842//_____________________________________________________________________________
843
845{
846 return fVerboseLevel;
847}
848
849//_____________________________________________________________________________
850
852{
853 fNameExtension = nameExtension;
854}
855
856//_____________________________________________________________________________
857
859{
860 return fNameExtension;
861}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ThreadLocal T * G4GeomSplitter< T >::offset
G4PVReplica(const G4String &pName, G4LogicalVolume *pLogical, G4LogicalVolume *pMother, const EAxis pAxis, const G4int nReplicas, const G4double width, const G4double offset=0.)
G4PVReplica represents many touchable detector elements differing only in their positioning....
std::map< G4LogicalVolume *, G4LogicalVolume *, std::less< G4LogicalVolume * > > G4ReflectedVolumesMap
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
HepGeom::Transform3D G4Transform3D
HepGeom::ScaleZ3D G4ScaleZ3D
HepGeom::Translate3D G4Translate3D
HepGeom::Scale3D G4Scale3D
HepGeom::Rotate3D G4Rotate3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
G4VSolid * GetSolid() const
const G4VisAttributes * GetVisAttributes() const
G4VSensitiveDetector * GetSensitiveDetector() const
std::size_t GetNoDaughters() const
void SetRegion(G4Region *reg)
G4bool IsRootRegion() const
G4Region * GetRegion() const
G4double GetBiasWeight() const
G4bool IsRegion() const
G4Material * GetMaterial() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
G4UserLimits * GetUserLimits() const
G4FieldManager * GetFieldManager() const
void SetVisAttributes(const G4VisAttributes *pVA)
void SetBiasWeight(G4double w)
const G4String & GetName() const
G4PVPlacement represents a single volume positioned within and relative to a mother volume.
G4ReflectionFactory provides functions for volumes placements with a general transfomation that can c...
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)
G4ReflectionFactory(const G4ReflectionFactory &)=delete
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:292
G4VPVDivisionFactory is an abstract factory that defines the interfaces for creating volume divisions...
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
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
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
Transform3D inverse() const
EAxis
Definition geomdefs.hh:54
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668
#define G4ThreadLocal
Definition tls.hh:77