Geant4-11
Functions | Variables
G4BuildGeom.cc File Reference
#include <iomanip>
#include <fstream>
#include "globals.hh"
#include "G3toG4.hh"
#include "G3MatTable.hh"
#include "G3MedTable.hh"
#include "G3RotTable.hh"
#include "G3VolTable.hh"
#include "G3PartTable.hh"
#include "G3DetTable.hh"
#include "G3toG4BuildTree.hh"
#include "G4GeometryManager.hh"
#include "G4LogicalVolume.hh"
#include "G4LogicalVolumeStore.hh"
#include "G4PVPlacement.hh"
#include "G4VisAttributes.hh"

Go to the source code of this file.

Functions

void checkVol ()
 
void checkVol (G4LogicalVolume *, G4int)
 
void G3CLRead (G4String &, char *)
 
G4LogicalVolumeG4BuildGeom (G4String &inFile)
 

Variables

std::ofstream ofile
 

Function Documentation

◆ checkVol() [1/2]

void checkVol ( )

Definition at line 109 of file G4BuildGeom.cc.

110{
112 G4LogicalVolume* ll = (*theStore)[0];
113 G4int level=0;
114 checkVol(ll, level);
115}
void checkVol(G4LogicalVolume *, G4int)
Definition: G4BuildGeom.cc:117
int G4int
Definition: G4Types.hh:85
static G4LogicalVolumeStore * GetInstance()

References checkVol(), and G4LogicalVolumeStore::GetInstance().

◆ checkVol() [2/2]

void checkVol ( G4LogicalVolume _lvol,
G4int  level 
)

Definition at line 117 of file G4BuildGeom.cc.

118{
119 G4LogicalVolume* _ldvol;
120 G4VPhysicalVolume* _pdvol;
121 level++;
122
123 G4int ndau = _lvol -> GetNoDaughters();
124
125 G4cout << "G44LogicalVolume " << _lvol->GetName() << " at level " << level
126 << " contains " << ndau << " daughters." << G4endl;
127 for (G4int idau=0; idau<ndau; idau++){
128 _pdvol = _lvol-> GetDaughter(idau);
129 _ldvol = _pdvol -> GetLogicalVolume();
130 G4cout << "G4VPhysical volume " << std::setw(5) << _pdvol -> GetName()
131 << " (G4LogicalVolume " << std::setw(5) << _ldvol->GetName() << ")"
132 << G4endl;
133 checkVol(_ldvol, level);
134 }
135 return;
136}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const

References checkVol(), G4cout, G4endl, and G4LogicalVolume::GetName().

Referenced by checkVol(), and G4BuildGeom().

◆ G3CLRead()

void G3CLRead ( G4String fname,
char *  select 
)

Definition at line 98 of file clparse.cc.

99{
100 //
101 // G3CLRead
102 // Read the call List file, parse the tokens, and pass the token
103 // List to the Geant4 interpreter
104 //
105 // fname: call List filename
106
107 G4String line;
108 G4String tokens[1000];
109
110 const char* ofname = "clparse.out";
111 ofile.open(ofname);
112 ofile << "Output file open\n";
113
114 G4int count = 0;
115 G4int ntokens = 0;
116 std::ifstream istr(fname);
117
118 while (G4StrUtil::readline(istr, line) && ! istr.eof())
119 {
120 count++;
121 ntokens = G3CLTokens(&line,tokens); // tokenize the line
122 for (G4int i=0; i < ntokens; i++)
123 {
124 ofile << tokens[i] << G4endl;
125 }
126
127 // interpret the line as a Geant call
128 //
129 G3CLEval(tokens, select);
130 }
131}
void G3CLEval(G4String *tokens, char *select)
G4int G3CLTokens(G4String *line, G4String *tokens)
std::ofstream ofile
Definition: clparse.cc:44
std::istream & readline(std::istream &is, G4String &str, G4bool skipWhite=true)
Read characters into a G4String from an input stream until end-of-line.
string fname
Definition: test.py:308

References test::fname, G3CLEval(), G3CLTokens(), G4endl, ofile, and G4StrUtil::readline().

Referenced by G4BuildGeom().

◆ G4BuildGeom()

G4LogicalVolume * G4BuildGeom ( G4String inFile)

Definition at line 54 of file G4BuildGeom.cc.

54 {
55
56 G4int irot=0;
57 G4gsrotm(0, 90, 0, 90, 90, 0, 0);
58
59 G4cout << "Instantiated unit rotation matrix irot=" << irot << G4endl;
60
61 // Read the call List and interpret to Generate Geant4 geometry
62
63 G4cout << "Reading the call List file " << inFile << "..." << G4endl;
64
65 G3CLRead(inFile, 0);
66
68
70
72
73 G4cout << "Call List file read completed. Build geometry" << G4endl;
74
75 // Build the geometry
76
78 G4cout << "G3toG4 top level volume is " << topVTE->GetName() << G4endl;
79
80 // modified
81 G3toG4BuildTree(topVTE, 0);
82
83 // Retrieve the top-level G3toG4 logical mother volume pointer
84
85 G4LogicalVolume* topLV = topVTE->GetLV();
86
87 // position the top logical volume
88 // (in Geant3 the top volume is not positioned)
89 //
90 new G4PVPlacement(0, G4ThreeVector(), topLV->GetName(), topLV, 0, false, 0);
91
92 // mark as invisible
93
95
96 G4cout << "Top-level G3toG4 logical volume " << topLV->GetName() << " "
97 << *(topLV->GetVisAttributes()) << G4endl;
98
99 // check the geometry here
100
101 #ifdef G3G4DEBUG
102 G4cout << "scan through G4LogicalVolumeStore:" << G4endl;
103 checkVol();
104 #endif
105
106 return topLV;
107}
G3G4DLL_API G3DetTable G3Det
Definition: clparse.cc:58
void G4gsrotm(G4int irot, G4double theta1, G4double phi1, G4double theta2, G4double phi2, G4double theta3, G4double phi3)
Definition: G4gsrotm.cc:53
G3G4DLL_API G3PartTable G3Part
Definition: clparse.cc:57
G3G4DLL_API G3VolTable G3Vol
Definition: clparse.cc:53
void G3toG4BuildTree(G3VolTableEntry *curVTE, G3VolTableEntry *motherVTE)
void G3CLRead(G4String &, char *)
Definition: clparse.cc:98
CLHEP::Hep3Vector G4ThreeVector
void PrintAll()
Definition: G3DetTable.cc:95
void PrintAll()
Definition: G3PartTable.cc:78
G4LogicalVolume * GetLV()
void PrintAll()
Definition: G3VolTable.cc:60
G3VolTableEntry * GetFirstVTE()
Definition: G3VolTable.cc:106
const G4VisAttributes * GetVisAttributes() const
void SetVisAttributes(const G4VisAttributes *pVA)
static const G4VisAttributes & GetInvisible()

References checkVol(), G3CLRead(), G3Det, G3Part, G3toG4BuildTree(), G3Vol, G4cout, G4endl, G4gsrotm(), G3VolTable::GetFirstVTE(), G4VisAttributes::GetInvisible(), G3VolTableEntry::GetLV(), G3VolTableEntry::GetName(), G4LogicalVolume::GetName(), G4LogicalVolume::GetVisAttributes(), G3DetTable::PrintAll(), G3PartTable::PrintAll(), G3VolTable::PrintAll(), and G4LogicalVolume::SetVisAttributes().

Variable Documentation

◆ ofile

std::ofstream ofile
extern