Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G3toG4MANY.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//
28// By I. Hrivnacova, 22.10.01
29
30//#define G3G4DEBUG 1
31
32#include "globals.hh"
33#include "G3toG4MANY.hh"
34#include "G3Pos.hh"
35#include "G3RotTable.hh"
36#include "G4SubtractionSolid.hh"
37
39{
40 if (curVTE->GetNoOverlaps() > 0) {
41
42 // check consistency
43 if (!curVTE->HasMANY()) {
44 G4String text = "G3toG4MANY: volume ";
45 text = text + curVTE->GetName() + " has specified overlaps \n";
46 text = text + " but is not defined as MANY.";
47 G4Exception("G3toG4MANY()", "G3toG40009",
48 FatalException, text);
49 return;
50 }
51
52 // only MANY volumes with one position are supported
53 if (curVTE->NPCopies() != 1) {
54 G4String text = "G3toG4MANY: volume ";
55 text = text + curVTE->GetName() + " which has MANY has not just one position.";
56 G4Exception("G3toG4MANY()", "G3toG40010",
57 FatalException, text);
58 return;
59 }
60
61 #ifdef G3G4DEBUG
62 G4cout << "G3toG4MANY " << curVTE->GetName() << " boolean" << G4endl;
63 #endif
64
65 G4Transform3D transform = GetTransform3D(curVTE->GetG3PosCopy(0));
66
67 MakeBooleanSolids(curVTE, curVTE->GetOverlaps(), transform.inverse());
68 }
69
70 // process daughters
71 for (G4int i=0; i<curVTE->GetNoDaughters(); i++)
72 G3toG4MANY(curVTE->GetDaughter(i));
73}
74
76 const G4Transform3D& transform)
77{
78 // loop over overlap VTEs
79 for (size_t i=0; i<overlaps->size(); i++){
80
81 G3VolTableEntry* overlapVTE = (*overlaps)[i];
82
83 // loop over clone VTEs
84 for (G4int ij=0; ij<overlapVTE->GetMasterClone()->GetNoClones(); ij++){
85
86 G3VolTableEntry* cloneVTE = overlapVTE->GetMasterClone()->GetClone(ij);
87
88 // loop over clone positions
89 for (G4int j=0; j<cloneVTE->NPCopies(); j++){
90
91 #ifdef G3G4DEBUG
92 G4cout << "From '" << curVTE->GetName() << "' "
93 << "cut '" << cloneVTE->GetName() << "' :"
94 << i << "th overlap (from " << overlaps->size() << ") "
95 << ij << "th clone (from " << overlapVTE->GetMasterClone()->GetNoClones() << ") "
96 << j << "th copy (from " << cloneVTE->NPCopies() << ") "
97 << G4endl;
98 #endif
99
100 SubstractSolids(curVTE, cloneVTE, j, transform);
101 }
102 }
103 }
104}
105
107 G4int copy, const G4Transform3D& transform)
108{
109 // vte2 transformation
110 G4Transform3D transform2 = GetTransform3D(vte2->GetG3PosCopy(copy));
111
112 // compose new name
113 G4String newName = vte1->GetSolid()->GetName();
114 newName = newName + "-" + vte2->GetSolid()->GetName();
115
116 #ifdef G3G4DEBUG
117 G4cout << " " << newName << G4endl;
118 #endif
119
120 G4VSolid* newSolid
121 = new G4SubtractionSolid(newName, vte1->GetSolid(), vte2->GetSolid(),
122 transform*transform2);
123
124 // update vte1
125 vte1->SetSolid(newSolid);
126
127 // process daughters
128 for (G4int k=0; k<vte1->GetNoDaughters(); k++){
129
130 G3VolTableEntry* dVTE = vte1->GetDaughter(k);
131
132 if (dVTE->NPCopies() != 1) {
133 G4String text = "G3toG4MANY: volume ";
134 text = text + dVTE->GetName() + " which has MANY has not just one position.";
135 G4Exception("G3toG4MANY()", "G3toG40011",
136 FatalException, text);
137 return;
138 }
139
141 SubstractSolids(dVTE, vte2, copy, dt.inverse()*transform);
142 }
143}
144
146{
147 G4int irot = g3pos->GetIrot();
148 G4RotationMatrix* theMatrix = 0;
149 if (irot>0) theMatrix = G3Rot.Get(irot);
150
151 G4Rotate3D rotation;
152 if (theMatrix) {
153 rotation = G4Rotate3D(*theMatrix);
154 }
155
156 G4Translate3D translation(*(g3pos->GetPos()));
157 G4Transform3D transform3D = translation * (rotation.inverse());
158
159 return transform3D;
160}
161
162
G3G4DLL_API G3RotTable G3Rot
Definition: clparse.cc:56
void G3toG4MANY(G3VolTableEntry *curVTE)
Definition: G3toG4MANY.cc:38
G4Transform3D GetTransform3D(G3Pos *g3pos)
Definition: G3toG4MANY.cc:145
void MakeBooleanSolids(G3VolTableEntry *curVTE, G3VolTableEntryVector *overlaps, const G4Transform3D &transform)
Definition: G3toG4MANY.cc:75
void SubstractSolids(G3VolTableEntry *vte1, G3VolTableEntry *vte2, G4int copy, const G4Transform3D &transform)
Definition: G3toG4MANY.cc:106
std::vector< G3VolTableEntry * > G3VolTableEntryVector
Definition: G3toG4MANY.hh:68
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
HepGeom::Rotate3D G4Rotate3D
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Definition: G3Pos.hh:43
G4ThreeVector * GetPos()
Definition: G3Pos.cc:72
G4int GetIrot()
Definition: G3Pos.cc:62
G4RotationMatrix * Get(G4int id) const
Definition: G3RotTable.cc:43
std::vector< G3VolTableEntry * > * GetOverlaps()
G3VolTableEntry * GetMasterClone()
G4VSolid * GetSolid()
G3VolTableEntry * GetClone(G4int i)
void SetSolid(G4VSolid *solid)
G3Pos * GetG3PosCopy(G4int copy=0)
G3VolTableEntry * GetDaughter(G4int i)
G4String GetName() const
Transform3D inverse() const
Definition: Transform3D.cc:141