LArSoft  v09_90_00
Liquid Argon Software toolkit - https://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 44 of file MyDetectorConstruction.hh.

Constructor & Destructor Documentation

MyDetectorConstruction::MyDetectorConstruction ( )

Definition at line 48 of file MyDetectorConstruction.cc.

49  : scoreVoxel(0)
51 {
52 }
MyDetectorConstruction::~MyDetectorConstruction ( )

Definition at line 55 of file MyDetectorConstruction.cc.

57 {
58 }

Member Function Documentation

G4VPhysicalVolume * MyDetectorConstruction::Construct ( )
virtual

Definition at line 61 of file MyDetectorConstruction.cc.

References scoreVoxel.

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

Definition at line 129 of file MyDetectorConstruction.cc.

References scoreVoxel.

Referenced by BOOST_PYTHON_MODULE().

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

Member Data Documentation

G4LogicalVolume* MyDetectorConstruction::scoreVoxel
private

Definition at line 46 of file MyDetectorConstruction.hh.

Referenced by Construct(), and SetSDtoScoreVoxel().


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