Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:20:44

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 biasing/B03/src/B03ImportanceDetectorConstruction.cc
0027 /// \brief Implementation of the B03ImportanceDetectorConstruction class
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 // For Primitive Scorers
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0055 
0056 B03ImportanceDetectorConstruction::B03ImportanceDetectorConstruction(G4String worldName)
0057   : G4VUserParallelWorld(worldName), fLogicalVolumeVector()
0058 {
0059   //  Construct();
0060 }
0061 
0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0063 
0064 B03ImportanceDetectorConstruction::~B03ImportanceDetectorConstruction()
0065 {
0066   fLogicalVolumeVector.clear();
0067 }
0068 
0069 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0070 
0071 void B03ImportanceDetectorConstruction::Construct()
0072 {
0073   G4cout << " constructing parallel world " << G4endl;
0074 
0075   G4Material* dummyMat = 0;
0076 
0077   // GetWorld methods create a clone of the mass world to the parallel world (!)
0078   //  via the transportation manager
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   //  fPVolumeStore.AddPVolume(G4GeometryCell(*pWorldVolume, 0));
0087   fPVolumeStore.AddPVolume(G4GeometryCell(*fGhostWorld, 0));
0088 
0089   // creating 18 slobs of 10 cm thicknes
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   // logical parallel cells
0101 
0102   G4LogicalVolume* aShield_log_imp = new G4LogicalVolume(aShield, dummyMat, "aShield_log_imp");
0103   fLogicalVolumeVector.push_back(aShield_log_imp);
0104 
0105   // physical parallel cells
0106 
0107   G4int i = 1;
0108   G4double startz = -85 * cm;
0109   //  for (i=1; i<=18; ++i) {
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     //                        0);
0119     G4GeometryCell cell(*pvol, i);
0120     //    G4GeometryCell cell(*pvol, 0);
0121     fPVolumeStore.AddPVolume(cell);
0122   }
0123 
0124   // filling the rest of the world volumr behind the concrete with
0125   // another slob which should get the same importance value as the
0126   // last slob
0127   innerRadiusShield = 0 * cm;
0128   //  outerRadiusShield = 110*cm; exceeds world volume!!!!
0129   outerRadiusShield = 100 * cm;
0130   //  heightShield       = 10*cm;
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   //  G4double pos_z = 100*cm;
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   //                      0);
0151   G4GeometryCell cell(*pvol, 19);
0152   //  G4GeometryCell cell(*pvol, 0);
0153   fPVolumeStore.AddPVolume(cell);
0154 
0155   SetSensitive();
0156 }
0157 
0158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0159 
0160 const G4VPhysicalVolume&
0161 B03ImportanceDetectorConstruction::GetPhysicalVolumeByName(const G4String& name) const
0162 {
0163   return *fPVolumeStore.GetPVolume(name);
0164 }
0165 
0166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0167 
0168 G4String B03ImportanceDetectorConstruction::ListPhysNamesAsG4String()
0169 {
0170   G4String names(fPVolumeStore.GetPNames());
0171   return names;
0172 }
0173 
0174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0206 
0207 G4VPhysicalVolume& B03ImportanceDetectorConstruction::GetWorldVolumeAddress() const
0208 {
0209   return *fGhostWorld;
0210 }
0211 
0212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0213 
0214 G4VPhysicalVolume* B03ImportanceDetectorConstruction::GetWorldVolume()
0215 {
0216   return fGhostWorld;
0217 }
0218 
0219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0220 
0221 void B03ImportanceDetectorConstruction::SetSensitive()
0222 {
0223   //  -------------------------------------------------
0224   //   The collection names of defined Primitives are
0225   //   0       ConcreteSD/Collisions
0226   //   1       ConcreteSD/CollWeight
0227   //   2       ConcreteSD/Population
0228   //   3       ConcreteSD/TrackEnter
0229   //   4       ConcreteSD/SL
0230   //   5       ConcreteSD/SLW
0231   //   6       ConcreteSD/SLWE
0232   //   7       ConcreteSD/SLW_V
0233   //   8       ConcreteSD/SLWE_V
0234   //  -------------------------------------------------
0235 
0236   // now moved to ConstructSD()
0237 }
0238 
0239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0240 void B03ImportanceDetectorConstruction::ConstructSD()
0241 {
0242   G4SDManager* SDman = G4SDManager::GetSDMpointer();
0243   //
0244   // Sensitive Detector Name
0245   G4String concreteSDname = "ConcreteSD";
0246 
0247   //------------------------
0248   // MultiFunctionalDetector
0249   //------------------------
0250   //
0251   // Define MultiFunctionalDetector with name.
0252   G4MultiFunctionalDetector* MFDet = new G4MultiFunctionalDetector(concreteSDname);
0253   SDman->AddNewDetector(MFDet);  // Register SD to SDManager
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     //      (*it)->SetSensitiveDetector(MFDet);
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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......