source: trunk/source/digits_hits/utils/test/test2/src/Test2PSConstruction.cc@ 1355

Last change on this file since 1355 was 1350, checked in by garnier, 15 years ago

update to last version 4.9.4

File size: 5.6 KB
Line 
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#include "Test2PSConstruction.hh"
29
30#include "G4LogicalVolume.hh"
31
32Test2PSConstruction::Test2PSConstruction(const G4String& name, G4int Segment[3])
33 :fname(name)
34{
35 nSegment[0] = Segment[0];
36 nSegment[1] = Segment[1];
37 nSegment[2] = Segment[2];
38}
39
40Test2PSConstruction::~Test2PSConstruction()
41{;}
42
43
44#include "G4SDManager.hh"
45
46
47#include "G4SDManager.hh"
48#include "G4MultiFunctionalDetector.hh"
49#include "G4PSEnergyDeposit3D.hh"
50#include "G4PSTrackLength3D.hh"
51#include "G4PSPassageTrackLength3D.hh"
52#include "G4PSDoseDeposit3D.hh"
53#include "G4PSFlatSurfaceCurrent3D.hh"
54#include "G4PSSphereSurfaceCurrent3D.hh"
55#include "G4PSPassageCellCurrent3D.hh"
56#include "G4PSFlatSurfaceFlux3D.hh"
57#include "G4PSCellFlux3D.hh"
58#include "G4PSPassageCellFlux3D.hh"
59#include "G4PSMinKinEAtGeneration3D.hh"
60#include "G4PSNofSecondary3D.hh"
61#include "G4PSCellCharge3D.hh"
62#include "G4PSNofStep3D.hh"
63#include "G4SDParticleFilter.hh"
64
65void Test2PSConstruction::SetupSensitivity(G4LogicalVolume* logVol) {
66 //
67 // sensitive detectors
68 //
69 G4SDManager * sdManager = G4SDManager::GetSDMpointer();
70 sdManager->SetVerboseLevel(1);
71
72 G4MultiFunctionalDetector* det = new G4MultiFunctionalDetector(fname);
73
74 // filters
75 G4String filterName, particleName;
76 G4SDParticleFilter* gammaFilter
77 = new G4SDParticleFilter(filterName="gammaFilter",particleName="gamma");
78 G4SDParticleFilter* electronFilter
79 = new G4SDParticleFilter(filterName="electronFilter",particleName="e-");
80 G4SDParticleFilter* positronFilter
81 = new G4SDParticleFilter(filterName="positronFilter",particleName="e+");
82
83 // primitive scorers
84 G4VPrimitiveScorer* primitive;
85 primitive =
86 new G4PSEnergyDeposit3D("eDep",nSegment[0],nSegment[1],nSegment[2]);
87 det->RegisterPrimitive(primitive);
88 primitive =
89 new G4PSTrackLength3D("trackLengthGamma",nSegment[0],nSegment[1],nSegment[2]);
90 primitive->SetFilter(gammaFilter);
91 det->RegisterPrimitive(primitive);
92 primitive =
93 new G4PSTrackLength3D("trackLengthElec",nSegment[0],nSegment[1],nSegment[2]);
94 primitive->SetFilter(electronFilter);
95 det->RegisterPrimitive(primitive);
96 primitive =
97 new G4PSTrackLength3D("trackLengthPosi",nSegment[0],nSegment[1],nSegment[2]);
98 primitive->SetFilter(positronFilter);
99 det->RegisterPrimitive(primitive);
100 primitive =
101 new G4PSNofStep3D("nStepGamma",nSegment[0],nSegment[1],nSegment[2]);
102 primitive->SetFilter(gammaFilter);
103 det->RegisterPrimitive(primitive);
104 primitive =
105 new G4PSNofStep3D("nStepElec",nSegment[0],nSegment[1],nSegment[2]);
106 primitive->SetFilter(electronFilter);
107 det->RegisterPrimitive(primitive);
108 primitive =
109 new G4PSNofStep3D("nStepPosi",nSegment[0],nSegment[1],nSegment[2]);
110 primitive->SetFilter(positronFilter);
111 det->RegisterPrimitive(primitive);
112
113 //-----
114 primitive = new G4PSPassageTrackLength3D("PassageTrackLength",nSegment[0],nSegment[1],nSegment[2]);
115 det->RegisterPrimitive(primitive);
116
117 primitive = new G4PSDoseDeposit3D("DoseDeposit",nSegment[0],nSegment[1],nSegment[2]);
118 det->RegisterPrimitive(primitive);
119
120 primitive = new G4PSFlatSurfaceCurrent3D("FlatSurfaceCurrent",fCurrent_In);
121 det->RegisterPrimitive(primitive);
122
123 primitive = new G4PSPassageCellCurrent3D("PassageCellCurrent",nSegment[0],nSegment[1],nSegment[2]);
124 det->RegisterPrimitive(primitive);
125
126 primitive = new G4PSFlatSurfaceFlux3D("FlatSurfaceFlux",fFlux_In);
127 det->RegisterPrimitive(primitive);
128
129 primitive = new G4PSCellFlux3D("CellFlux",nSegment[0],nSegment[1],nSegment[2]);
130 det->RegisterPrimitive(primitive);
131
132 primitive = new G4PSPassageCellFlux3D("PassageCellFlux",nSegment[0],nSegment[1],nSegment[2]);
133 det->RegisterPrimitive(primitive);
134
135 primitive = new G4PSNofSecondary3D("NofSecondary",nSegment[0],nSegment[1],nSegment[2]);
136 det->RegisterPrimitive(primitive);
137
138 primitive = new G4PSCellCharge3D("CellCharge",nSegment[0],nSegment[1],nSegment[2]);
139 det->RegisterPrimitive(primitive);
140
141 sdManager->AddNewDetector(det);
142 logVol->SetSensitiveDetector(det);
143 sdManager->SetVerboseLevel(0);
144
145}
Note: See TracBrowser for help on using the repository browser.