File indexing completed on 2025-02-23 09:20:44
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032 #include "B03ImportanceDetectorConstruction.hh"
0033
0034 #include "G4LogicalVolume.hh"
0035 #include "G4Material.hh"
0036 #include "G4PVPlacement.hh"
0037 #include "G4PhysicalConstants.hh"
0038 #include "G4SystemOfUnits.hh"
0039 #include "G4ThreeVector.hh"
0040 #include "G4Tubs.hh"
0041 #include "globals.hh"
0042
0043 #include <sstream>
0044
0045
0046 #include "G4MultiFunctionalDetector.hh"
0047 #include "G4PSNofCollision.hh"
0048 #include "G4PSPopulation.hh"
0049 #include "G4PSTrackCounter.hh"
0050 #include "G4PSTrackLength.hh"
0051 #include "G4SDManager.hh"
0052 #include "G4SDParticleFilter.hh"
0053
0054
0055
0056 B03ImportanceDetectorConstruction::B03ImportanceDetectorConstruction(G4String worldName)
0057 : G4VUserParallelWorld(worldName), fLogicalVolumeVector()
0058 {
0059
0060 }
0061
0062
0063
0064 B03ImportanceDetectorConstruction::~B03ImportanceDetectorConstruction()
0065 {
0066 fLogicalVolumeVector.clear();
0067 }
0068
0069
0070
0071 void B03ImportanceDetectorConstruction::Construct()
0072 {
0073 G4cout << " constructing parallel world " << G4endl;
0074
0075 G4Material* dummyMat = 0;
0076
0077
0078
0079 fGhostWorld = GetWorld();
0080 G4cout << " B03ImportanceDetectorConstruction:: ghostWorldName = " << fGhostWorld->GetName()
0081 << G4endl;
0082 G4LogicalVolume* worldLogical = fGhostWorld->GetLogicalVolume();
0083 fLogicalVolumeVector.push_back(worldLogical);
0084
0085 G4String name("none");
0086
0087 fPVolumeStore.AddPVolume(G4GeometryCell(*fGhostWorld, 0));
0088
0089
0090
0091 G4double innerRadiusShield = 0 * cm;
0092 G4double outerRadiusShield = 100 * cm;
0093 G4double heightShield = 5 * cm;
0094 G4double startAngleShield = 0 * deg;
0095 G4double spanningAngleShield = 360 * deg;
0096
0097 G4Tubs* aShield = new G4Tubs("aShield", innerRadiusShield, outerRadiusShield, heightShield,
0098 startAngleShield, spanningAngleShield);
0099
0100
0101
0102 G4LogicalVolume* aShield_log_imp = new G4LogicalVolume(aShield, dummyMat, "aShield_log_imp");
0103 fLogicalVolumeVector.push_back(aShield_log_imp);
0104
0105
0106
0107 G4int i = 1;
0108 G4double startz = -85 * cm;
0109
0110 for (i = 1; i <= 18; i++) {
0111 name = GetCellName(i);
0112
0113 G4double pos_x = 0 * cm;
0114 G4double pos_y = 0 * cm;
0115 G4double pos_z = startz + (i - 1) * (2 * heightShield);
0116 G4VPhysicalVolume* pvol = new G4PVPlacement(0, G4ThreeVector(pos_x, pos_y, pos_z),
0117 aShield_log_imp, name, worldLogical, false, i);
0118
0119 G4GeometryCell cell(*pvol, i);
0120
0121 fPVolumeStore.AddPVolume(cell);
0122 }
0123
0124
0125
0126
0127 innerRadiusShield = 0 * cm;
0128
0129 outerRadiusShield = 100 * cm;
0130
0131 heightShield = 5 * cm;
0132 startAngleShield = 0 * deg;
0133 spanningAngleShield = 360 * deg;
0134
0135 G4Tubs* aRest = new G4Tubs("Rest", innerRadiusShield, outerRadiusShield, heightShield,
0136 startAngleShield, spanningAngleShield);
0137
0138 G4LogicalVolume* aRest_log = new G4LogicalVolume(aRest, dummyMat, "aRest_log");
0139
0140 fLogicalVolumeVector.push_back(aRest_log);
0141
0142 name = GetCellName(19);
0143
0144 G4double pos_x = 0 * cm;
0145 G4double pos_y = 0 * cm;
0146
0147 G4double pos_z = 95 * cm;
0148 G4VPhysicalVolume* pvol = new G4PVPlacement(0, G4ThreeVector(pos_x, pos_y, pos_z), aRest_log,
0149 name, worldLogical, false, 19);
0150
0151 G4GeometryCell cell(*pvol, 19);
0152
0153 fPVolumeStore.AddPVolume(cell);
0154
0155 SetSensitive();
0156 }
0157
0158
0159
0160 const G4VPhysicalVolume&
0161 B03ImportanceDetectorConstruction::GetPhysicalVolumeByName(const G4String& name) const
0162 {
0163 return *fPVolumeStore.GetPVolume(name);
0164 }
0165
0166
0167
0168 G4String B03ImportanceDetectorConstruction::ListPhysNamesAsG4String()
0169 {
0170 G4String names(fPVolumeStore.GetPNames());
0171 return names;
0172 }
0173
0174
0175
0176 G4String B03ImportanceDetectorConstruction::GetCellName(G4int i)
0177 {
0178 std::ostringstream os;
0179 os << "cell_";
0180 if (i < 10) {
0181 os << "0";
0182 }
0183 os << i;
0184 G4String name = os.str();
0185 return name;
0186 }
0187
0188
0189
0190 G4GeometryCell B03ImportanceDetectorConstruction::GetGeometryCell(G4int i)
0191 {
0192 G4String name(GetCellName(i));
0193 const G4VPhysicalVolume* p = 0;
0194 p = fPVolumeStore.GetPVolume(name);
0195 if (p) {
0196 return G4GeometryCell(*p, 0);
0197 }
0198 else {
0199 G4cout << "B03ImportanceDetectorConstruction::GetGeometryCell: " << G4endl
0200 << " couldn't get G4GeometryCell" << G4endl;
0201 return G4GeometryCell(*fGhostWorld, -2);
0202 }
0203 }
0204
0205
0206
0207 G4VPhysicalVolume& B03ImportanceDetectorConstruction::GetWorldVolumeAddress() const
0208 {
0209 return *fGhostWorld;
0210 }
0211
0212
0213
0214 G4VPhysicalVolume* B03ImportanceDetectorConstruction::GetWorldVolume()
0215 {
0216 return fGhostWorld;
0217 }
0218
0219
0220
0221 void B03ImportanceDetectorConstruction::SetSensitive()
0222 {
0223
0224
0225
0226
0227
0228
0229
0230
0231
0232
0233
0234
0235
0236
0237 }
0238
0239
0240 void B03ImportanceDetectorConstruction::ConstructSD()
0241 {
0242 G4SDManager* SDman = G4SDManager::GetSDMpointer();
0243
0244
0245 G4String concreteSDname = "ConcreteSD";
0246
0247
0248
0249
0250
0251
0252 G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
0253 SDman->AddNewDetector(MFDet);
0254
0255 G4String fltName, particleName;
0256 G4SDParticleFilter* neutronFilter =
0257 new G4SDParticleFilter(fltName = "neutronFilter", particleName = "neutron");
0258
0259 MFDet->SetFilter(neutronFilter);
0260
0261 for (std::vector<G4LogicalVolume*>::iterator it = fLogicalVolumeVector.begin();
0262 it != fLogicalVolumeVector.end(); it++)
0263 {
0264
0265 SetSensitiveDetector((*it)->GetName(), MFDet);
0266 }
0267
0268 G4String psName;
0269 G4PSNofCollision* scorer0 = new G4PSNofCollision(psName = "Collisions");
0270 MFDet->RegisterPrimitive(scorer0);
0271
0272 G4PSNofCollision* scorer1 = new G4PSNofCollision(psName = "CollWeight");
0273 scorer1->Weighted(true);
0274 MFDet->RegisterPrimitive(scorer1);
0275
0276 G4PSPopulation* scorer2 = new G4PSPopulation(psName = "Population");
0277 MFDet->RegisterPrimitive(scorer2);
0278
0279 G4PSTrackCounter* scorer3 = new G4PSTrackCounter(psName = "TrackEnter", fCurrent_In);
0280 MFDet->RegisterPrimitive(scorer3);
0281
0282 G4PSTrackLength* scorer4 = new G4PSTrackLength(psName = "SL");
0283 MFDet->RegisterPrimitive(scorer4);
0284
0285 G4PSTrackLength* scorer5 = new G4PSTrackLength(psName = "SLW");
0286 scorer5->Weighted(true);
0287 MFDet->RegisterPrimitive(scorer5);
0288
0289 G4PSTrackLength* scorer6 = new G4PSTrackLength(psName = "SLWE");
0290 scorer6->Weighted(true);
0291 scorer6->MultiplyKineticEnergy(true);
0292 MFDet->RegisterPrimitive(scorer6);
0293
0294 G4PSTrackLength* scorer7 = new G4PSTrackLength(psName = "SLW_V");
0295 scorer7->Weighted(true);
0296 scorer7->DivideByVelocity(true);
0297 MFDet->RegisterPrimitive(scorer7);
0298
0299 G4PSTrackLength* scorer8 = new G4PSTrackLength(psName = "SLWE_V");
0300 scorer8->Weighted(true);
0301 scorer8->MultiplyKineticEnergy(true);
0302 scorer8->DivideByVelocity(true);
0303 MFDet->RegisterPrimitive(scorer8);
0304 }
0305
0306