Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G3toG4MANY.cc File Reference
#include "globals.hh"
#include "G3toG4MANY.hh"
#include "G3Pos.hh"
#include "G3RotTable.hh"
#include "G4SubtractionSolid.hh"

Go to the source code of this file.

Functions

void G3toG4MANY (G3VolTableEntry *curVTE)
 
void MakeBooleanSolids (G3VolTableEntry *curVTE, G3VolTableEntryVector *overlaps, const G4Transform3D &transform)
 
void SubstractSolids (G3VolTableEntry *vte1, G3VolTableEntry *vte2, G4int copy, const G4Transform3D &transform)
 
G4Transform3D GetTransform3D (G3Pos *g3pos)
 

Function Documentation

◆ G3toG4MANY()

void G3toG4MANY ( G3VolTableEntry curVTE)

Definition at line 39 of file G3toG4MANY.cc.

40{
41 if (curVTE->GetNoOverlaps() > 0) {
42
43 // check consistency
44 if (!curVTE->HasMANY()) {
45 G4String text = "G3toG4MANY: volume ";
46 text = text + curVTE->GetName() + " has specified overlaps \n";
47 text = text + " but is not defined as MANY.";
48 G4Exception("G3toG4MANY()", "G3toG40009",
49 FatalException, text);
50 return;
51 }
52
53 // only MANY volumes with one position are supported
54 if (curVTE->NPCopies() != 1) {
55 G4String text = "G3toG4MANY: volume ";
56 text = text + curVTE->GetName() + " which has MANY has not just one position.";
57 G4Exception("G3toG4MANY()", "G3toG40010",
58 FatalException, text);
59 return;
60 }
61
62 #ifdef G3G4DEBUG
63 G4cout << "G3toG4MANY " << curVTE->GetName() << " boolean" << G4endl;
64 #endif
65
66 G4Transform3D transform = GetTransform3D(curVTE->GetG3PosCopy(0));
67
68 MakeBooleanSolids(curVTE, curVTE->GetOverlaps(), transform.inverse());
69 }
70
71 // process daughters
72 for (G4int i=0; i<curVTE->GetNoDaughters(); i++)
73 G3toG4MANY(curVTE->GetDaughter(i));
74}
void G3toG4MANY(G3VolTableEntry *curVTE)
Definition: G3toG4MANY.cc:39
G4Transform3D GetTransform3D(G3Pos *g3pos)
Definition: G3toG4MANY.cc:146
void MakeBooleanSolids(G3VolTableEntry *curVTE, G3VolTableEntryVector *overlaps, const G4Transform3D &transform)
Definition: G3toG4MANY.cc:76
@ FatalException
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
std::vector< G3VolTableEntry * > * GetOverlaps()
G3Pos * GetG3PosCopy(G4int copy=0)
G3VolTableEntry * GetDaughter(G4int i)
Transform3D inverse() const
Definition: Transform3D.cc:142
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Referenced by G3toG4MANY().

◆ GetTransform3D()

G4Transform3D GetTransform3D ( G3Pos g3pos)

Definition at line 146 of file G3toG4MANY.cc.

147{
148 G4int irot = g3pos->GetIrot();
149 G4RotationMatrix* theMatrix = 0;
150 if (irot>0) theMatrix = G3Rot.Get(irot);
151
152 G4Rotate3D rotation;
153 if (theMatrix) {
154 rotation = G4Rotate3D(*theMatrix);
155 }
156
157 G4Translate3D translation(*(g3pos->GetPos()));
158 G4Transform3D transform3D = translation * (rotation.inverse());
159
160 return transform3D;
161}
G3G4DLL_API G3RotTable G3Rot
Definition: clparse.cc:57
HepGeom::Rotate3D G4Rotate3D
G4ThreeVector * GetPos()
Definition: G3Pos.cc:73
G4int GetIrot()
Definition: G3Pos.cc:63
G4RotationMatrix * Get(G4int id) const
Definition: G3RotTable.cc:44

Referenced by G3toG4MANY(), and SubstractSolids().

◆ MakeBooleanSolids()

void MakeBooleanSolids ( G3VolTableEntry curVTE,
G3VolTableEntryVector overlaps,
const G4Transform3D transform 
)

Definition at line 76 of file G3toG4MANY.cc.

78{
79 // loop over overlap VTEs
80 for (size_t i=0; i<overlaps->size(); i++){
81
82 G3VolTableEntry* overlapVTE = (*overlaps)[i];
83
84 // loop over clone VTEs
85 for (G4int ij=0; ij<overlapVTE->GetMasterClone()->GetNoClones(); ij++){
86
87 G3VolTableEntry* cloneVTE = overlapVTE->GetMasterClone()->GetClone(ij);
88
89 // loop over clone positions
90 for (G4int j=0; j<cloneVTE->NPCopies(); j++){
91
92 #ifdef G3G4DEBUG
93 G4cout << "From '" << curVTE->GetName() << "' "
94 << "cut '" << cloneVTE->GetName() << "' :"
95 << i << "th overlap (from " << overlaps->size() << ") "
96 << ij << "th clone (from " << overlapVTE->GetMasterClone()->GetNoClones() << ") "
97 << j << "th copy (from " << cloneVTE->NPCopies() << ") "
98 << G4endl;
99 #endif
100
101 SubstractSolids(curVTE, cloneVTE, j, transform);
102 }
103 }
104 }
105}
void SubstractSolids(G3VolTableEntry *vte1, G3VolTableEntry *vte2, G4int copy, const G4Transform3D &transform)
Definition: G3toG4MANY.cc:107
G3VolTableEntry * GetMasterClone()
G3VolTableEntry * GetClone(G4int i)

Referenced by G3toG4MANY().

◆ SubstractSolids()

void SubstractSolids ( G3VolTableEntry vte1,
G3VolTableEntry vte2,
G4int  copy,
const G4Transform3D transform 
)

Definition at line 107 of file G3toG4MANY.cc.

109{
110 // vte2 transformation
111 G4Transform3D transform2 = GetTransform3D(vte2->GetG3PosCopy(copy));
112
113 // compose new name
114 G4String newName = vte1->GetSolid()->GetName();
115 newName = newName + "-" + vte2->GetSolid()->GetName();
116
117 #ifdef G3G4DEBUG
118 G4cout << " " << newName << G4endl;
119 #endif
120
121 G4VSolid* newSolid
122 = new G4SubtractionSolid(newName, vte1->GetSolid(), vte2->GetSolid(),
123 transform*transform2);
124
125 // update vte1
126 vte1->SetSolid(newSolid);
127
128 // process daughters
129 for (G4int k=0; k<vte1->GetNoDaughters(); k++){
130
131 G3VolTableEntry* dVTE = vte1->GetDaughter(k);
132
133 if (dVTE->NPCopies() != 1) {
134 G4String text = "G3toG4MANY: volume ";
135 text = text + dVTE->GetName() + " which has MANY has not just one position.";
136 G4Exception("G3toG4MANY()", "G3toG40011",
137 FatalException, text);
138 return;
139 }
140
142 SubstractSolids(dVTE, vte2, copy, dt.inverse()*transform);
143 }
144}
G4VSolid * GetSolid()
void SetSolid(G4VSolid *solid)
G4String GetName() const

Referenced by MakeBooleanSolids(), and SubstractSolids().