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 "Par04ParallelFullWorld.hh"
0029
0030 #include "Par04ParallelFullSensitiveDetector.hh"
0031 #include "Par04ParallelMessenger.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 Par04ParallelFullWorld::Par04ParallelFullWorld(G4String aWorldName,
0053 const Par04DetectorConstruction* aMassDetector)
0054 : G4VUserParallelWorld(aWorldName), fMassDetector(aMassDetector)
0055 {
0056 fParallelMessenger = new Par04ParallelMessenger(this);
0057 fNbOfLayers = fMassDetector->GetNbOfLayers();
0058 }
0059
0060
0061
0062 Par04ParallelFullWorld::~Par04ParallelFullWorld() = default;
0063
0064
0065
0066 void Par04ParallelFullWorld::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 G4double sensitiveLayerOffset = 0;
0079 if (fMassDetector->GetAbsorberSensitivity(0))
0080 fLayerThickness = fMassDetector->GetAbsorberThickness(0);
0081 else
0082 sensitiveLayerOffset = fMassDetector->GetAbsorberThickness(0);
0083 if (fMassDetector->GetAbsorberSensitivity(1))
0084 fLayerThickness += fMassDetector->GetAbsorberThickness(1);
0085
0086 fNbOfLayers = fMassDetector->GetNbOfLayers();
0087 G4double detectorRadius = fNbOfLayers * fullLayerThickness;
0088 G4double detectorOuterRadius = detectorInnerRadius + detectorRadius;
0089 G4double rowThickness = detectorLength / fNbOfRows;
0090 G4double full2Pi = 2. * CLHEP::pi * rad;
0091 Print();
0092
0093
0094
0095 auto solidDetector = new G4Tubs("Detector",
0096 detectorInnerRadius,
0097 detectorOuterRadius,
0098 detectorLength / 2.,
0099 0,
0100 full2Pi);
0101 auto logicDetector = new G4LogicalVolume(solidDetector,
0102 dummy,
0103 "Detector");
0104 new G4PVPlacement(0,
0105 G4ThreeVector(0, 0, 0),
0106 logicDetector,
0107 "Detector",
0108 parallelLogicalVolume,
0109 false,
0110 9999,
0111 true);
0112
0113
0114 auto solidRow =
0115 new G4Tubs("Row", detectorInnerRadius, detectorOuterRadius, rowThickness / 2., 0, full2Pi);
0116
0117 auto logicRow = new G4LogicalVolume(solidRow, dummy, "Row");
0118 if (fNbOfRows > 1)
0119 new G4PVReplica("Row", logicRow, logicDetector, kZAxis, fNbOfRows, rowThickness);
0120 else
0121 new G4PVPlacement(0, G4ThreeVector(), logicRow, "Row", logicDetector, false, 0);
0122
0123
0124 G4double cellPhi = full2Pi / fNbOfSlices;
0125 auto solidSlice =
0126 new G4Tubs("Slice", detectorInnerRadius, detectorOuterRadius, rowThickness / 2, 0, cellPhi);
0127 auto logicSlice = new G4LogicalVolume(solidSlice, dummy, "Slice");
0128 if (fNbOfLayers > 1 && fullLayerThickness == fLayerThickness) {
0129 new G4PVReplica("Slice", logicSlice, logicRow, kPhi, fNbOfSlices, cellPhi, -cellPhi);
0130 }
0131 else {
0132
0133 for (int iSlice = 0; iSlice < fNbOfSlices; iSlice++) {
0134 auto rotation = new G4RotationMatrix();
0135 rotation->setPhi((iSlice + 0.5) * cellPhi);
0136 new G4PVPlacement(rotation, G4ThreeVector(), logicSlice, "Slice_" + std::to_string(iSlice),
0137 logicRow, false, iSlice);
0138 }
0139 }
0140
0141
0142 G4VisAttributes attribs;
0143 attribs.SetColour(G4Colour(0, 1, 0, 0.1));
0144 attribs.SetForceSolid(true);
0145 if (fNbOfLayers > 1 && fullLayerThickness == fLayerThickness) {
0146 auto solidCell = new G4Tubs("Cell", detectorInnerRadius + sensitiveLayerOffset,
0147 detectorInnerRadius + sensitiveLayerOffset + fLayerThickness,
0148 rowThickness / 2, 0, cellPhi);
0149 fLogicalCell.push_back(new G4LogicalVolume(solidCell, dummy, "Cell_0"));
0150 new G4PVReplica("Cell", fLogicalCell.back(), logicSlice, kRho, fNbOfLayers, fLayerThickness,
0151 detectorInnerRadius);
0152 }
0153 else {
0154
0155 for (int iLayer = 0; iLayer < fNbOfLayers; iLayer++) {
0156 auto solidCell = new G4Tubs(
0157 "Cell_" + std::to_string(iLayer),
0158 detectorInnerRadius + iLayer * fullLayerThickness + sensitiveLayerOffset,
0159 detectorInnerRadius + iLayer * fullLayerThickness + sensitiveLayerOffset + fLayerThickness,
0160 rowThickness / 2, 0, cellPhi);
0161 fLogicalCell.push_back(
0162 new G4LogicalVolume(solidCell, dummy, "Cell_" + std::to_string(iLayer)));
0163 fLogicalCell.back()->SetVisAttributes(attribs);
0164 new G4PVPlacement(0, G4ThreeVector(), fLogicalCell.back(), "Cell_" + std::to_string(iLayer),
0165 logicSlice, false, iLayer);
0166 }
0167 }
0168 Print();
0169 }
0170
0171
0172
0173 void Par04ParallelFullWorld::ConstructSD()
0174 {
0175
0176 G4SDManager* SDman = G4SDManager::GetSDMpointer();
0177 Par04ParallelFullSensitiveDetector* caloSD =
0178 new Par04ParallelFullSensitiveDetector("parallelFullSD", fNbOfLayers, fNbOfSlices, fNbOfRows);
0179 SDman->AddNewDetector(caloSD);
0180 for (const auto& logicalCell : fLogicalCell)
0181 logicalCell->SetSensitiveDetector(caloSD);
0182 }
0183
0184
0185
0186 void Par04ParallelFullWorld::Print()
0187 {
0188 G4cout << "\n------------------------------------------------------"
0189 << "\n Readout geometry with physics layout is set in parallel geometry:\t"
0190 << "\n Cylindrical detector is divided along radius (layers), phi (slices), and z (rows)."
0191 << "\n Number of layers is determined by number of layers set in detector construction. "
0192 << "\n- Number of layers: " << fNbOfLayers << "\n------- Number of slices: " << fNbOfSlices
0193 << "\n- Number of rows: " << fNbOfRows;
0194 G4cout << "\n Readout will collect energy for full simulation.\n------- Therefore thickness is "
0195 << "only a thickness of sensitive absorbers = " << G4BestUnit(fLayerThickness, "Length")
0196 << "\n-----------------------------------------------------" << G4endl;
0197 }