Geant4
9.6.p02
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
geant4_9_6_p02
source
g3tog4
src
G3toG4BuildTree.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$
28
//
29
// modified by I. Hrivnacova, 2.8.99
30
31
#include "
globals.hh
"
32
#include "
G3toG4BuildTree.hh
"
33
#include "
G3RotTable.hh
"
34
#include "
G3MedTable.hh
"
35
#include "
G3VolTable.hh
"
36
#include "
G3SensVolVector.hh
"
37
#include "
G3Pos.hh
"
38
#include "
G4LogicalVolume.hh
"
39
#include "
G4PVPlacement.hh
"
40
#include "
G4ReflectionFactory.hh
"
41
#include "
G4Transform3D.hh
"
42
43
void
G3toG4BuildTree
(
G3VolTableEntry
* curVTE,
G3VolTableEntry
* motherVTE)
44
{
45
G3toG4BuildLVTree
(curVTE, motherVTE);
46
G3toG4BuildPVTree
(curVTE);
47
}
48
49
void
G3toG4BuildLVTree
(
G3VolTableEntry
* curVTE,
G3VolTableEntry
* motherVTE)
50
{
51
// check existence of the solid
52
if
(curVTE->
GetSolid
()) {
53
G4LogicalVolume
* curLog = curVTE->
GetLV
();
54
if
(!curLog) {
55
// skip creating logical volume
56
// in case it already exists
57
58
// material
59
G4Material
*
material
= 0;
60
G3MedTableEntry
* mte =
G3Med
.
get
(curVTE->
GetNmed
());
61
if
(mte) material = mte->
GetMaterial
();
62
if
(!material) {
63
G4String
err_message =
"VTE "
+ curVTE->
GetName
()
64
+
" has not defined material!!"
;
65
G4Exception
(
"G3toG4BuildLVTree()"
,
"G3toG40001"
,
66
FatalException
, err_message);
67
return
;
68
}
69
70
// logical volume
71
curLog =
72
new
G4LogicalVolume
(curVTE->
GetSolid
(),
material
, curVTE->
GetName
());
73
curVTE->
SetLV
(curLog);
74
75
// insert logical volume to G3SensVol vector
76
// in case it is sensitive
77
if
(mte->
GetISVOL
())
G3SensVol
.push_back(curLog);
78
}
79
}
80
else
{
81
if
( !(curVTE->
GetDivision
() && motherVTE->
GetMasterClone
() == motherVTE &&
82
motherVTE->
GetNoClones
()>1)) {
83
// ignore dummy vte's
84
// (this should be the only case when the vte is dummy but
85
// is present in mother <-> daughters tree
86
G4String
err_message =
"VTE "
+ curVTE->
GetName
()
87
+
" has not defined solid!!"
;
88
G4Exception
(
"G3toG4BuildLVTree()"
,
"G3toG40002"
,
89
FatalException
, err_message);
90
return
;
91
}
92
}
93
94
// process daughters
95
G4int
Ndau = curVTE->
GetNoDaughters
();
96
for
(
int
Idau=0; Idau<Ndau; Idau++){
97
G3toG4BuildLVTree
(curVTE->
GetDaughter
(Idau), curVTE);
98
}
99
}
100
101
void
G3toG4BuildPVTree
(
G3VolTableEntry
* curVTE)
102
{
103
// check existence of the solid
104
if
(curVTE->
GetSolid
()) {
105
G4LogicalVolume
* curLog = curVTE->
GetLV
();
106
107
// positions in motherVTE
108
for
(
G4int
i=0; i<curVTE->
NPCopies
(); i++){
109
110
G3Pos
* theG3Pos = curVTE->
GetG3PosCopy
(i);
111
if
(theG3Pos) {
112
113
// loop over all mothers
114
for
(
G4int
im=0; im<curVTE->
GetNoMothers
(); im++) {
115
116
G3VolTableEntry
* motherVTE = curVTE->
GetMother
(im);
117
if
(theG3Pos->
GetMotherName
() == motherVTE->
GetMasterClone
()->
GetName
()) {
118
119
// get mother logical volume
120
G4LogicalVolume
* mothLV=0;
121
G4String
motherName = motherVTE->
GetName
();
122
if
(!curVTE->
FindMother
(motherName))
continue
;
123
if
(curVTE->
FindMother
(motherName)->
GetName
() != motherName) {
124
// check consistency - tbr
125
G4String
err_message =
126
"G3toG4BuildTree: Inconsistent mother <-> daughter !!"
;
127
G4Exception
(
"G3toG4BuildPVTree()"
,
"G3toG40003"
,
128
FatalException
, err_message);
129
return
;
130
}
131
mothLV = motherVTE->
GetLV
();
132
133
// copy number
134
// (in G3 numbering starts from 1 but in G4 from 0)
135
G4int
copyNo = theG3Pos->
GetCopy
() - 1;
136
137
// position it if not top-level volume
138
139
if
(mothLV != 0) {
140
141
// transformation
142
G4int
irot = theG3Pos->
GetIrot
();
143
G4RotationMatrix
* theMatrix = 0;
144
if
(irot>0) theMatrix =
G3Rot
.
Get
(irot);
145
G4Rotate3D
rotation;
146
if
(theMatrix) {
147
rotation =
G4Rotate3D
(*theMatrix);
148
}
149
150
#ifndef G3G4_NO_REFLECTION
151
G4Translate3D
translation(*(theG3Pos->
GetPos
()));
152
G4Transform3D
transform3D = translation * (rotation.
inverse
());
153
154
G4ReflectionFactory::Instance
()
155
->
Place
(transform3D,
// transformation
156
curVTE->
GetName
(),
// PV name
157
curLog,
// its logical volume
158
mothLV,
// mother logical volume
159
false
,
// only
160
copyNo);
// copy
161
#else
162
new
G4PVPlacement
(theMatrix,
// rotation matrix
163
*(theG3Pos->
GetPos
()),
// its position
164
curLog,
// its LogicalVolume
165
curVTE->
GetName
(),
// PV name
166
mothLV,
// Mother LV
167
0,
// only
168
copyNo);
// copy
169
#endif
170
171
// verbose
172
#ifdef G3G4DEBUG
173
G4cout
<<
"PV: "
<< i <<
"th copy of "
<< curVTE->
GetName
()
174
<<
" in "
<< motherVTE->
GetName
() <<
" copyNo: "
175
<< copyNo <<
" irot: "
<< irot <<
" pos: "
176
<< *(theG3Pos->
GetPos
()) <<
G4endl
;
177
#endif
178
}
179
}
180
}
181
// clear this position
182
curVTE->
ClearG3PosCopy
(i);
183
i--;
184
}
185
}
186
187
// divisions
188
if
(curVTE->
GetDivision
()) {
189
curVTE->
GetDivision
()->
CreatePVReplica
();
190
// verbose
191
#ifdef G3G4DEBUG
192
G4cout
<<
"CreatePVReplica: "
<< curVTE->
GetName
()
193
<<
" in "
<< curVTE->
GetMother
()->
GetName
() <<
G4endl
;
194
#endif
195
196
// clear this divison
197
curVTE->
ClearDivision
();
198
}
199
}
200
201
// process daughters
202
G4int
Ndau = curVTE->
GetNoDaughters
();
203
for
(
int
Idau=0; Idau<Ndau; Idau++){
204
G3toG4BuildPVTree
(curVTE->
GetDaughter
(Idau));
205
}
206
}
207
Generated on Sat May 25 2013 14:33:09 for Geant4 by
1.8.4