Geant4-11
Functions
G3toG4BuildTree.hh File Reference
#include "G3VolTableEntry.hh"

Go to the source code of this file.

Functions

void G3toG4BuildLVTree (G3VolTableEntry *curVTE, G3VolTableEntry *motherVTE)
 
void G3toG4BuildPVTree (G3VolTableEntry *curVTE)
 
void G3toG4BuildTree (G3VolTableEntry *curVTE, G3VolTableEntry *motherVTE)
 

Function Documentation

◆ G3toG4BuildLVTree()

void G3toG4BuildLVTree ( G3VolTableEntry curVTE,
G3VolTableEntry motherVTE 
)

Definition at line 48 of file G3toG4BuildTree.cc.

49{
50 // check existence of the solid
51 if (curVTE->GetSolid()) {
52 G4LogicalVolume* curLog = curVTE->GetLV();
53 if (!curLog) {
54 // skip creating logical volume
55 // in case it already exists
56
57 // material
59 G3MedTableEntry* mte = G3Med.get(curVTE->GetNmed());
60 if (mte) material = mte->GetMaterial();
61 if (!material) {
62 G4String err_message = "VTE " + curVTE->GetName()
63 + " has not defined material!!";
64 G4Exception("G3toG4BuildLVTree()", "G3toG40001",
65 FatalException, err_message);
66 return;
67 }
68
69 // logical volume
70 curLog =
71 new G4LogicalVolume(curVTE->GetSolid(), material, curVTE->GetName());
72 curVTE->SetLV(curLog);
73
74 // insert logical volume to G3SensVol vector
75 // in case it is sensitive
76 if (mte->GetISVOL()) G3SensVol.push_back(curLog);
77 }
78 }
79 else {
80 if ( !(curVTE->GetDivision() && motherVTE->GetMasterClone() == motherVTE &&
81 motherVTE->GetNoClones()>1)) {
82 // ignore dummy vte's
83 // (this should be the only case when the vte is dummy but
84 // is present in mother <-> daughters tree
85 G4String err_message = "VTE " + curVTE->GetName()
86 + " has not defined solid!!";
87 G4Exception("G3toG4BuildLVTree()", "G3toG40002",
88 FatalException, err_message);
89 return;
90 }
91 }
92
93 // process daughters
94 G4int Ndau = curVTE->GetNoDaughters();
95 for (int Idau=0; Idau<Ndau; Idau++){
96 G3toG4BuildLVTree(curVTE->GetDaughter(Idau), curVTE);
97 }
98}
G3G4DLL_API G3MedTable G3Med
Definition: clparse.cc:55
G3G4DLL_API G3SensVolVector G3SensVol
Definition: clparse.cc:60
void G3toG4BuildLVTree(G3VolTableEntry *curVTE, G3VolTableEntry *motherVTE)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
int G4int
Definition: G4Types.hh:85
G4Material * GetMaterial() const
G4int GetISVOL() const
G3MedTableEntry * get(G4int id) const
Definition: G3MedTable.cc:43
G3VolTableEntry * GetMasterClone()
G4VSolid * GetSolid()
G3Division * GetDivision()
G4LogicalVolume * GetLV()
void SetLV(G4LogicalVolume *lv)
G3VolTableEntry * GetDaughter(G4int i)
string material
Definition: eplot.py:19

References FatalException, G3Med, G3SensVol, G3toG4BuildLVTree(), G4Exception(), G3MedTable::get(), G3VolTableEntry::GetDaughter(), G3VolTableEntry::GetDivision(), G3MedTableEntry::GetISVOL(), G3VolTableEntry::GetLV(), G3VolTableEntry::GetMasterClone(), G3MedTableEntry::GetMaterial(), G3VolTableEntry::GetName(), G3VolTableEntry::GetNmed(), G3VolTableEntry::GetNoClones(), G3VolTableEntry::GetNoDaughters(), G3VolTableEntry::GetSolid(), eplot::material, and G3VolTableEntry::SetLV().

Referenced by G3toG4BuildLVTree(), and G3toG4BuildTree().

◆ G3toG4BuildPVTree()

void G3toG4BuildPVTree ( G3VolTableEntry curVTE)

Definition at line 100 of file G3toG4BuildTree.cc.

101{
102 // check existence of the solid
103 if (curVTE->GetSolid()) {
104 G4LogicalVolume* curLog = curVTE->GetLV();
105
106 // positions in motherVTE
107 for (G4int i=0; i<curVTE->NPCopies(); i++){
108
109 G3Pos* theG3Pos = curVTE->GetG3PosCopy(i);
110 if (theG3Pos) {
111
112 // loop over all mothers
113 for (G4int im=0; im<curVTE->GetNoMothers(); im++) {
114
115 G3VolTableEntry* motherVTE = curVTE->GetMother(im);
116 if (theG3Pos->GetMotherName() == motherVTE->GetMasterClone()->GetName()) {
117
118 // get mother logical volume
119 G4LogicalVolume* mothLV=0;
120 G4String motherName = motherVTE->GetName();
121 if (!curVTE->FindMother(motherName)) continue;
122 if (curVTE->FindMother(motherName)->GetName() != motherName) {
123 // check consistency - tbr
124 G4String err_message =
125 "G3toG4BuildTree: Inconsistent mother <-> daughter !!";
126 G4Exception("G3toG4BuildPVTree()", "G3toG40003",
127 FatalException, err_message);
128 return;
129 }
130 mothLV = motherVTE->GetLV();
131
132 // copy number
133 // (in G3 numbering starts from 1 but in G4 from 0)
134 G4int copyNo = theG3Pos->GetCopy() - 1;
135
136 // position it if not top-level volume
137
138 if (mothLV != 0) {
139
140 // transformation
141 G4int irot = theG3Pos->GetIrot();
142 G4RotationMatrix* theMatrix = 0;
143 if (irot>0) theMatrix = G3Rot.Get(irot);
144 G4Rotate3D rotation;
145 if (theMatrix) {
146 rotation = G4Rotate3D(*theMatrix);
147 }
148
149 #ifndef G3G4_NO_REFLECTION
150 G4Translate3D translation(*(theG3Pos->GetPos()));
151 G4Transform3D transform3D = translation * (rotation.inverse());
152
154 ->Place(transform3D, // transformation
155 curVTE->GetName(), // PV name
156 curLog, // its logical volume
157 mothLV, // mother logical volume
158 false, // only
159 copyNo); // copy
160 #else
161 new G4PVPlacement(theMatrix, // rotation matrix
162 *(theG3Pos->GetPos()), // its position
163 curLog, // its LogicalVolume
164 curVTE->GetName(), // PV name
165 mothLV, // Mother LV
166 0, // only
167 copyNo); // copy
168 #endif
169
170 // verbose
171 #ifdef G3G4DEBUG
172 G4cout << "PV: " << i << "th copy of " << curVTE->GetName()
173 << " in " << motherVTE->GetName() << " copyNo: "
174 << copyNo << " irot: " << irot << " pos: "
175 << *(theG3Pos->GetPos()) << G4endl;
176 #endif
177 }
178 }
179 }
180 // clear this position
181 curVTE->ClearG3PosCopy(i);
182 i--;
183 }
184 }
185
186 // divisions
187 if (curVTE->GetDivision()) {
188 curVTE->GetDivision()->CreatePVReplica();
189 // verbose
190 #ifdef G3G4DEBUG
191 G4cout << "CreatePVReplica: " << curVTE->GetName()
192 << " in " << curVTE->GetMother()->GetName() << G4endl;
193 #endif
194
195 // clear this divison
196 curVTE->ClearDivision();
197 }
198 }
199
200 // process daughters
201 G4int Ndau = curVTE->GetNoDaughters();
202 for (int Idau=0; Idau<Ndau; Idau++){
203 G3toG4BuildPVTree(curVTE->GetDaughter(Idau));
204 }
205}
G3G4DLL_API G3RotTable G3Rot
Definition: clparse.cc:56
void G3toG4BuildPVTree(G3VolTableEntry *curVTE)
HepGeom::Rotate3D G4Rotate3D
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void CreatePVReplica()
Definition: G3Division.cc:127
Definition: G3Pos.hh:43
G4ThreeVector * GetPos()
Definition: G3Pos.cc:72
G4int GetCopy()
Definition: G3Pos.cc:67
G4String & GetMotherName()
Definition: G3Pos.cc:57
G4int GetIrot()
Definition: G3Pos.cc:62
G4RotationMatrix * Get(G4int id) const
Definition: G3RotTable.cc:43
G3VolTableEntry * FindMother(const G4String &vname)
G3VolTableEntry * GetMother(G4int i)
void ClearG3PosCopy(G4int copy)
G3Pos * GetG3PosCopy(G4int copy=0)
static G4ReflectionFactory * Instance()
G4PhysicalVolumesPair Place(const G4Transform3D &transform3D, const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, G4bool isMany, G4int copyNo, G4bool surfCheck=false)
Transform3D inverse() const
Definition: Transform3D.cc:141

References G3VolTableEntry::ClearDivision(), G3VolTableEntry::ClearG3PosCopy(), G3Division::CreatePVReplica(), FatalException, G3VolTableEntry::FindMother(), G3Rot, G3toG4BuildPVTree(), G4cout, G4endl, G4Exception(), G3RotTable::Get(), G3Pos::GetCopy(), G3VolTableEntry::GetDaughter(), G3VolTableEntry::GetDivision(), G3VolTableEntry::GetG3PosCopy(), G3Pos::GetIrot(), G3VolTableEntry::GetLV(), G3VolTableEntry::GetMasterClone(), G3VolTableEntry::GetMother(), G3Pos::GetMotherName(), G3VolTableEntry::GetName(), G3VolTableEntry::GetNoDaughters(), G3VolTableEntry::GetNoMothers(), G3Pos::GetPos(), G3VolTableEntry::GetSolid(), G4ReflectionFactory::Instance(), HepGeom::Transform3D::inverse(), G3VolTableEntry::NPCopies(), and G4ReflectionFactory::Place().

Referenced by G3toG4BuildPVTree(), and G3toG4BuildTree().

◆ G3toG4BuildTree()

void G3toG4BuildTree ( G3VolTableEntry curVTE,
G3VolTableEntry motherVTE 
)

Definition at line 42 of file G3toG4BuildTree.cc.

43{
44 G3toG4BuildLVTree(curVTE, motherVTE);
45 G3toG4BuildPVTree(curVTE);
46}

References G3toG4BuildLVTree(), and G3toG4BuildPVTree().

Referenced by G4BuildGeom().