File indexing completed on 2025-02-23 09:22:35
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 #include "Par04ParallelFastWorld.hh"
0029
0030 #include "Par04ParallelFastSensitiveDetector.hh"
0031 #include "Par04ParallelFullWorld.hh"
0032
0033
0034 #include "G4AutoDelete.hh"
0035 #include "G4Colour.hh"
0036 #include "G4LogicalVolume.hh"
0037 #include "G4Material.hh"
0038 #include "G4NistManager.hh"
0039 #include "G4PVPlacement.hh"
0040 #include "G4PVReplica.hh"
0041 #include "G4SDManager.hh"
0042 #include "G4SystemOfUnits.hh"
0043 #include "G4ThreeVector.hh"
0044 #include "G4Tubs.hh"
0045 #include "G4UnitsTable.hh"
0046 #include "G4VPhysicalVolume.hh"
0047 #include "G4VisAttributes.hh"
0048 #include "globals.hh"
0049
0050
0051
0052 Par04ParallelFastWorld::Par04ParallelFastWorld(G4String aWorldName,
0053 const Par04DetectorConstruction* aMassDetector,
0054 const Par04ParallelFullWorld* aParallelFull)
0055 : G4VUserParallelWorld(aWorldName), fMassDetector(aMassDetector), fParallelFull(aParallelFull)
0056 {
0057 fNbOfLayers = fMassDetector->GetNbOfLayers();
0058 }
0059
0060
0061
0062 Par04ParallelFastWorld::~Par04ParallelFastWorld() = default;
0063
0064
0065
0066 void Par04ParallelFastWorld::Construct()
0067 {
0068
0069 G4Material* dummy = nullptr;
0070
0071
0072 auto parallelLogicalVolume = GetWorld()->GetLogicalVolume();
0073
0074 G4double detectorInnerRadius = fMassDetector->GetInnerRadius();
0075 G4double detectorLength = fMassDetector->GetLength();
0076 G4double fullLayerThickness =
0077 fMassDetector->GetAbsorberThickness(0) + fMassDetector->GetAbsorberThickness(1);
0078 fLayerThickness = fullLayerThickness;
0079
0080 fNbOfLayers = fMassDetector->GetNbOfLayers();
0081 fNbOfSlices = fParallelFull->GetNbOfSlices();
0082 fNbOfRows = fParallelFull->GetNbOfRows();
0083 G4double detectorRadius = fNbOfLayers * fullLayerThickness;
0084 G4double detectorOuterRadius = detectorInnerRadius + detectorRadius;
0085 G4double rowThickness = detectorLength / fNbOfRows;
0086 G4double full2Pi = 2. * CLHEP::pi * rad;
0087 Print();
0088
0089
0090
0091 auto solidDetector = new G4Tubs("Detector",
0092 detectorInnerRadius,
0093 detectorOuterRadius,
0094 detectorLength / 2.,
0095 0,
0096 full2Pi);
0097 auto logicDetector = new G4LogicalVolume(solidDetector,
0098 dummy,
0099 "Detector");
0100 new G4PVPlacement(0,
0101 G4ThreeVector(0, 0, 0),
0102 logicDetector,
0103 "Detector",
0104 parallelLogicalVolume,
0105 false,
0106 9999,
0107 true);
0108
0109
0110 auto solidRow =
0111 new G4Tubs("Row", detectorInnerRadius, detectorOuterRadius, rowThickness / 2., 0, full2Pi);
0112
0113 auto logicRow = new G4LogicalVolume(solidRow, dummy, "Row");
0114 if (fNbOfRows > 1)
0115 new G4PVReplica("Row", logicRow, logicDetector, kZAxis, fNbOfRows, rowThickness);
0116 else
0117 new G4PVPlacement(0, G4ThreeVector(), logicRow, "Row", logicDetector, false, 0);
0118
0119
0120 G4double cellPhi = full2Pi / fNbOfSlices;
0121 auto solidSlice =
0122 new G4Tubs("Slice", detectorInnerRadius, detectorOuterRadius, rowThickness / 2, 0, cellPhi);
0123 auto logicSlice = new G4LogicalVolume(solidSlice, dummy, "Slice");
0124 if (fNbOfSlices > 1) {
0125 new G4PVReplica("Slice", logicSlice, logicRow, kPhi, fNbOfSlices, cellPhi, -cellPhi);
0126 }
0127 else {
0128 new G4PVPlacement(0, G4ThreeVector(), logicSlice, "Slice", logicRow, false, 0);
0129 }
0130
0131
0132 G4VisAttributes attribs;
0133 attribs.SetColour(G4Colour(0, 1, 0, 0.1));
0134 attribs.SetForceSolid(true);
0135 if (fNbOfLayers > 1) {
0136 auto solidCell = new G4Tubs("Cell", detectorInnerRadius, detectorInnerRadius + fLayerThickness,
0137 rowThickness / 2, 0, cellPhi);
0138 fLogicalCell.push_back(new G4LogicalVolume(solidCell, dummy, "Cell_0"));
0139 new G4PVReplica("Cell", fLogicalCell.back(), logicSlice, kRho, fNbOfLayers, fLayerThickness,
0140 detectorInnerRadius);
0141 }
0142 else {
0143 auto solidCell = new G4Tubs("Cell", detectorInnerRadius, detectorInnerRadius + fLayerThickness,
0144 rowThickness / 2, 0, cellPhi);
0145 fLogicalCell.push_back(new G4LogicalVolume(solidCell, dummy, "Cell"));
0146 fLogicalCell.back()->SetVisAttributes(attribs);
0147 new G4PVPlacement(0, G4ThreeVector(), fLogicalCell.back(), "Cell", logicSlice, false, 0);
0148 }
0149 Print();
0150 }
0151
0152
0153
0154 void Par04ParallelFastWorld::ConstructSD()
0155 {
0156
0157 G4SDManager* SDman = G4SDManager::GetSDMpointer();
0158 Par04ParallelFastSensitiveDetector* caloSD =
0159 new Par04ParallelFastSensitiveDetector("parallelFastSD", fNbOfLayers, fNbOfSlices);
0160 SDman->AddNewDetector(caloSD);
0161 for (const auto& logicalCell : fLogicalCell)
0162 logicalCell->SetSensitiveDetector(caloSD);
0163 }
0164
0165
0166
0167 void Par04ParallelFastWorld::Print()
0168 {
0169 G4cout << "\n------------------------------------------------------"
0170 << "\n Readout geometry with physics layout is set in parallel geometry:\t"
0171 << "\n Cylindrical detector is divided along radius (layers), phi (slices), and z (rows)."
0172 << "\n Number of layers is determined by number of layers set in detector construction. "
0173 << "\n- Number of layers: " << fNbOfLayers << "\n------- Number of slices: " << fNbOfSlices
0174 << "\n- Number of rows: " << fNbOfRows;
0175 G4cout << "\n Readout will collect energy from fast simulation.\n------- Thickness is "
0176 << "a sum of all absorbers"
0177 << " = " << G4BestUnit(fLayerThickness, "Length")
0178 << "\n-----------------------------------------------------" << G4endl;
0179 }