source: trunk/source/g3tog4/src/G3toG4BuildTree.cc@ 1059

Last change on this file since 1059 was 965, checked in by garnier, 17 years ago

update g3tog4

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