LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
MyDetectorConstruction Class Reference

#include "MyDetectorConstruction.hh"

Inheritance diagram for MyDetectorConstruction:

Public Member Functions

 MyDetectorConstruction ()
 
 ~MyDetectorConstruction ()
 
virtual G4VPhysicalVolume * Construct ()
 
void SetSDtoScoreVoxel (G4VSensitiveDetector *asd)
 

Private Attributes

G4LogicalVolume * scoreVoxel
 

Detailed Description

Definition at line 45 of file MyDetectorConstruction.hh.

Constructor & Destructor Documentation

MyDetectorConstruction::MyDetectorConstruction ( )

Definition at line 49 of file MyDetectorConstruction.cc.

50  : scoreVoxel(0)
52 {
53 }
MyDetectorConstruction::~MyDetectorConstruction ( )

Definition at line 56 of file MyDetectorConstruction.cc.

58 {
59 }

Member Function Documentation

G4VPhysicalVolume * MyDetectorConstruction::Construct ( )
virtual

Definition at line 62 of file MyDetectorConstruction.cc.

References scoreVoxel.

64 {
65  G4Material* mate;
66  G4VisAttributes* va;
67 
68  // ==============================================================
69  // world volume
70  // ==============================================================
71  G4Box* areaSolid= new G4Box("area", 25.*cm, 25.*cm, 1.1*m);
72  G4Material* vacuum= G4Material::GetMaterial("Vacuum");
73  G4LogicalVolume* areaLV= new G4LogicalVolume(areaSolid, vacuum, "area");
74  G4PVPlacement* area= new G4PVPlacement(0, G4ThreeVector(), "area",
75  areaLV, 0, false, 0);
76  // vis. attributes
77  va= new G4VisAttributes(G4Color(1.,1.,1.));
78  va-> SetVisibility(false);
79  areaLV-> SetVisAttributes(va);
80 
81  // ==============================================================
82  // phantom
83  // ==============================================================
84  // water phantom
85  const double dxyphantom= 40.*cm;
86  const double dzphantom= 50.*cm;
87 
88  G4Box* sphantom= new G4Box("phantom",
89  dxyphantom/2., dxyphantom/2., dzphantom/2.);
90  G4Material* water= G4Material::GetMaterial("Water");
91  G4LogicalVolume* lphantom= new G4LogicalVolume(sphantom,
92  water, "phantom");
93  G4PVPlacement* phantom= new G4PVPlacement(0,
94  G4ThreeVector(0.,0., dzphantom/2.),
95  lphantom, "phantom",
96  areaLV, false, 0);
97 
98  va= new G4VisAttributes(G4Color(0.,0.1,0.8));
99  lphantom-> SetVisAttributes(va);
100 
101  // score voxels
102  const G4double dvoxel= 2.*mm;
103  const G4double dvoxel_y= 20.*mm;
104 
105  G4Box* svoxel= new G4Box("voxel", dvoxel/2., dvoxel_y/2., dvoxel/2.);
106  scoreVoxel= new G4LogicalVolume(svoxel, water, "voxel");
107  va= new G4VisAttributes(G4Color(0.,0.8,0.8));
108  va-> SetVisibility(false);
109  scoreVoxel-> SetVisAttributes(va);
110 
111  G4int ix, iz;
112  G4int index=0;
113  for (iz=0; iz<200; iz++) {
114  for (ix=-40; ix<=40; ix++) {
115  G4double x0= ix*dvoxel;
116  G4double z0= -dzphantom/2.+(iz+0.5)*dvoxel;
117  G4PVPlacement* pvoxel= new
118  G4PVPlacement(0, G4ThreeVector(x0, 0., z0),
119  scoreVoxel, "voxel", lphantom,
120  false, index);
121  index++;
122  }
123  }
124 
125  return area;
126 }
void MyDetectorConstruction::SetSDtoScoreVoxel ( G4VSensitiveDetector *  asd)

Definition at line 130 of file MyDetectorConstruction.cc.

References scoreVoxel.

Referenced by BOOST_PYTHON_MODULE().

132 {
133  if(scoreVoxel) {
134  scoreVoxel-> SetSensitiveDetector(asd);
135  }
136 }

Member Data Documentation

G4LogicalVolume* MyDetectorConstruction::scoreVoxel
private

Definition at line 47 of file MyDetectorConstruction.hh.

Referenced by Construct(), and SetSDtoScoreVoxel().


The documentation for this class was generated from the following files: