Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-03-30 07:49:58

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 /// \file B03ImportanceDetectorConstruction.cc
0027 /// \brief Implementation of the B03ImportanceDetectorConstruction class
0028 
0029 #include "B03ImportanceDetectorConstruction.hh"
0030 
0031 #include "G4LogicalVolume.hh"
0032 #include "G4Material.hh"
0033 #include "G4PVPlacement.hh"
0034 #include "G4PhysicalConstants.hh"
0035 #include "G4SystemOfUnits.hh"
0036 #include "G4ThreeVector.hh"
0037 #include "G4Tubs.hh"
0038 #include "globals.hh"
0039 
0040 #include <sstream>
0041 
0042 // For Primitive Scorers
0043 #include "G4MultiFunctionalDetector.hh"
0044 #include "G4PSNofCollision.hh"
0045 #include "G4PSPopulation.hh"
0046 #include "G4PSTrackCounter.hh"
0047 #include "G4PSTrackLength.hh"
0048 #include "G4SDManager.hh"
0049 #include "G4SDParticleFilter.hh"
0050 
0051 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0052 
0053 B03ImportanceDetectorConstruction::B03ImportanceDetectorConstruction(G4String worldName)
0054   : G4VUserParallelWorld(worldName), fLogicalVolumeVector()
0055 {
0056   //  Construct();
0057 }
0058 
0059 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0060 
0061 B03ImportanceDetectorConstruction::~B03ImportanceDetectorConstruction()
0062 {
0063   fLogicalVolumeVector.clear();
0064 }
0065 
0066 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0067 
0068 void B03ImportanceDetectorConstruction::Construct()
0069 {
0070   G4cout << " constructing parallel world " << G4endl;
0071 
0072   G4Material* dummyMat = 0;
0073 
0074   // GetWorld methods create a clone of the mass world to the parallel world (!)
0075   //  via the transportation manager
0076   fGhostWorld = GetWorld();
0077   G4cout << " B03ImportanceDetectorConstruction:: ghostWorldName = " << fGhostWorld->GetName()
0078          << G4endl;
0079   G4LogicalVolume* worldLogical = fGhostWorld->GetLogicalVolume();
0080   fLogicalVolumeVector.push_back(worldLogical);
0081 
0082   G4String name("none");
0083   //  fPVolumeStore.AddPVolume(G4GeometryCell(*pWorldVolume, 0));
0084   fPVolumeStore.AddPVolume(G4GeometryCell(*fGhostWorld, 0));
0085 
0086   // creating 18 slobs of 10 cm thicknes
0087 
0088   G4double innerRadiusShield = 0 * cm;
0089   G4double outerRadiusShield = 100 * cm;
0090   G4double heightShield = 5 * cm;
0091   G4double startAngleShield = 0 * deg;
0092   G4double spanningAngleShield = 360 * deg;
0093 
0094   G4Tubs* aShield = new G4Tubs("aShield", innerRadiusShield, outerRadiusShield, heightShield,
0095                                startAngleShield, spanningAngleShield);
0096 
0097   // logical parallel cells
0098 
0099   G4LogicalVolume* aShield_log_imp = new G4LogicalVolume(aShield, dummyMat, "aShield_log_imp");
0100   fLogicalVolumeVector.push_back(aShield_log_imp);
0101 
0102   // physical parallel cells
0103 
0104   G4int i = 1;
0105   G4double startz = -85 * cm;
0106   //  for (i=1; i<=18; ++i) {
0107   for (i = 1; i <= 18; i++) {
0108     name = GetCellName(i);
0109 
0110     G4double pos_x = 0 * cm;
0111     G4double pos_y = 0 * cm;
0112     G4double pos_z = startz + (i - 1) * (2 * heightShield);
0113     G4VPhysicalVolume* pvol = new G4PVPlacement(0, G4ThreeVector(pos_x, pos_y, pos_z),
0114                                                 aShield_log_imp, name, worldLogical, false, i);
0115     //                        0);
0116     G4GeometryCell cell(*pvol, i);
0117     //    G4GeometryCell cell(*pvol, 0);
0118     fPVolumeStore.AddPVolume(cell);
0119   }
0120 
0121   // filling the rest of the world volumr behind the concrete with
0122   // another slob which should get the same importance value as the
0123   // last slob
0124   innerRadiusShield = 0 * cm;
0125   //  outerRadiusShield = 110*cm; exceeds world volume!!!!
0126   outerRadiusShield = 100 * cm;
0127   //  heightShield       = 10*cm;
0128   heightShield = 5 * cm;
0129   startAngleShield = 0 * deg;
0130   spanningAngleShield = 360 * deg;
0131 
0132   G4Tubs* aRest = new G4Tubs("Rest", innerRadiusShield, outerRadiusShield, heightShield,
0133                              startAngleShield, spanningAngleShield);
0134 
0135   G4LogicalVolume* aRest_log = new G4LogicalVolume(aRest, dummyMat, "aRest_log");
0136 
0137   fLogicalVolumeVector.push_back(aRest_log);
0138 
0139   name = GetCellName(19);
0140 
0141   G4double pos_x = 0 * cm;
0142   G4double pos_y = 0 * cm;
0143   //  G4double pos_z = 100*cm;
0144   G4double pos_z = 95 * cm;
0145   G4VPhysicalVolume* pvol = new G4PVPlacement(0, G4ThreeVector(pos_x, pos_y, pos_z), aRest_log,
0146                                               name, worldLogical, false, 19);
0147   //                      0);
0148   G4GeometryCell cell(*pvol, 19);
0149   //  G4GeometryCell cell(*pvol, 0);
0150   fPVolumeStore.AddPVolume(cell);
0151 
0152   SetSensitive();
0153 }
0154 
0155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0156 
0157 const G4VPhysicalVolume&
0158 B03ImportanceDetectorConstruction::GetPhysicalVolumeByName(const G4String& name) const
0159 {
0160   return *fPVolumeStore.GetPVolume(name);
0161 }
0162 
0163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0164 
0165 G4String B03ImportanceDetectorConstruction::ListPhysNamesAsG4String()
0166 {
0167   G4String names(fPVolumeStore.GetPNames());
0168   return names;
0169 }
0170 
0171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0172 
0173 G4String B03ImportanceDetectorConstruction::GetCellName(G4int i)
0174 {
0175   std::ostringstream os;
0176   os << "cell_";
0177   if (i < 10) {
0178     os << "0";
0179   }
0180   os << i;
0181   G4String name = os.str();
0182   return name;
0183 }
0184 
0185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0186 
0187 G4GeometryCell B03ImportanceDetectorConstruction::GetGeometryCell(G4int i)
0188 {
0189   G4String name(GetCellName(i));
0190   const G4VPhysicalVolume* p = 0;
0191   p = fPVolumeStore.GetPVolume(name);
0192   if (p) {
0193     return G4GeometryCell(*p, 0);
0194   }
0195   else {
0196     G4cout << "B03ImportanceDetectorConstruction::GetGeometryCell: " << G4endl
0197            << " couldn't get G4GeometryCell" << G4endl;
0198     return G4GeometryCell(*fGhostWorld, -2);
0199   }
0200 }
0201 
0202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0203 
0204 G4VPhysicalVolume& B03ImportanceDetectorConstruction::GetWorldVolumeAddress() const
0205 {
0206   return *fGhostWorld;
0207 }
0208 
0209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0210 
0211 G4VPhysicalVolume* B03ImportanceDetectorConstruction::GetWorldVolume()
0212 {
0213   return fGhostWorld;
0214 }
0215 
0216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0217 
0218 void B03ImportanceDetectorConstruction::SetSensitive()
0219 {
0220   //  -------------------------------------------------
0221   //   The collection names of defined Primitives are
0222   //   0       ConcreteSD/Collisions
0223   //   1       ConcreteSD/CollWeight
0224   //   2       ConcreteSD/Population
0225   //   3       ConcreteSD/TrackEnter
0226   //   4       ConcreteSD/SL
0227   //   5       ConcreteSD/SLW
0228   //   6       ConcreteSD/SLWE
0229   //   7       ConcreteSD/SLW_V
0230   //   8       ConcreteSD/SLWE_V
0231   //  -------------------------------------------------
0232 
0233   // now moved to ConstructSD()
0234 }
0235 
0236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0237 void B03ImportanceDetectorConstruction::ConstructSD()
0238 {
0239   G4SDManager* SDman = G4SDManager::GetSDMpointer();
0240   //
0241   // Sensitive Detector Name
0242   G4String concreteSDname = "ConcreteSD";
0243 
0244   //------------------------
0245   // MultiFunctionalDetector
0246   //------------------------
0247   //
0248   // Define MultiFunctionalDetector with name.
0249   G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
0250   SDman->AddNewDetector(MFDet);  // Register SD to SDManager
0251 
0252   G4String fltName, particleName;
0253   G4SDParticleFilter* neutronFilter =
0254     new G4SDParticleFilter(fltName = "neutronFilter", particleName = "neutron");
0255 
0256   MFDet->SetFilter(neutronFilter);
0257 
0258   for (std::vector<G4LogicalVolume*>::iterator it = fLogicalVolumeVector.begin();
0259        it != fLogicalVolumeVector.end(); it++)
0260   {
0261     //      (*it)->SetSensitiveDetector(MFDet);
0262     SetSensitiveDetector((*it)->GetName(), MFDet);
0263   }
0264 
0265   G4String psName;
0266   G4PSNofCollision* scorer0 = new G4PSNofCollision(psName = "Collisions");
0267   MFDet->RegisterPrimitive(scorer0);
0268 
0269   G4PSNofCollision* scorer1 = new G4PSNofCollision(psName = "CollWeight");
0270   scorer1->Weighted(true);
0271   MFDet->RegisterPrimitive(scorer1);
0272 
0273   G4PSPopulation* scorer2 = new G4PSPopulation(psName = "Population");
0274   MFDet->RegisterPrimitive(scorer2);
0275 
0276   G4PSTrackCounter* scorer3 = new G4PSTrackCounter(psName = "TrackEnter", fCurrent_In);
0277   MFDet->RegisterPrimitive(scorer3);
0278 
0279   G4PSTrackLength* scorer4 = new G4PSTrackLength(psName = "SL");
0280   MFDet->RegisterPrimitive(scorer4);
0281 
0282   G4PSTrackLength* scorer5 = new G4PSTrackLength(psName = "SLW");
0283   scorer5->Weighted(true);
0284   MFDet->RegisterPrimitive(scorer5);
0285 
0286   G4PSTrackLength* scorer6 = new G4PSTrackLength(psName = "SLWE");
0287   scorer6->Weighted(true);
0288   scorer6->MultiplyKineticEnergy(true);
0289   MFDet->RegisterPrimitive(scorer6);
0290 
0291   G4PSTrackLength* scorer7 = new G4PSTrackLength(psName = "SLW_V");
0292   scorer7->Weighted(true);
0293   scorer7->DivideByVelocity(true);
0294   MFDet->RegisterPrimitive(scorer7);
0295 
0296   G4PSTrackLength* scorer8 = new G4PSTrackLength(psName = "SLWE_V");
0297   scorer8->Weighted(true);
0298   scorer8->MultiplyKineticEnergy(true);
0299   scorer8->DivideByVelocity(true);
0300   MFDet->RegisterPrimitive(scorer8);
0301 }
0302 
0303 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......