Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/parallel/ThreadsafeScorers/include/TSDetectorConstruction.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 /// \file parallel/ThreadsafeScorers/include/TSDetectorConstruction.hh
 27 /// \brief Definition of the TSDetectorConstruction class
 28 //
 29 //
 30 //
 31 //
 32 /// Construction of a target material (default = boron) surrounded by a
 33 ///     casing material (default = water) and a vacuum world (default =
 34 ///     target and casing fill world). The target + casing is brick
 35 ///     geometry with fTargetSections defining the number of divisions
 36 ///     in each dimension. The end sections in each dimension
 37 ///     is set to the casing. So a fTargetSections = G4ThreeVector(3, 3, 3)
 38 ///     would be one section of boron and 8 sections of water.
 39 /// The idea behind this geometry is just to create a simple geometry that
 40 ///     scatters and produces a lot neutrons with a minimal number of sections
 41 ///     (i.e. coarse meshing) such that the contention in operating on
 42 ///     the atomic hits maps is higher and round-off errors in the
 43 ///     thread-local hits maps are detectable (printed out in TSRunAction)
 44 ///     from the sheer number of floating point sum operations.
 45 /// Two scorers are implemented: EnergyDeposit and Number of steps
 46 ///     The energy deposit is to (possibly) show the round-off error seen
 47 ///     with thread-local hits maps. The # of steps scorer is to verify
 48 ///     the thread-safe and thread-local hits maps provide the same results.
 49 //
 50 //
 51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 53 
 54 #ifndef tsdetectorconstruction_hh
 55 #define tsdetectorconstruction_hh 1
 56 
 57 #include "G4ThreeVector.hh"
 58 #include "G4VUserDetectorConstruction.hh"
 59 #include "globals.hh"
 60 
 61 #include <map>
 62 #include <set>
 63 
 64 class G4Box;
 65 class G4Tubs;
 66 class G4Sphere;
 67 class G4LogicalVolume;
 68 class G4VPhysicalVolume;
 69 class G4Material;
 70 
 71 class TSDetectorConstruction : public G4VUserDetectorConstruction
 72 {
 73   public:
 74     typedef std::map<G4String, G4Material*> MaterialCollection_t;
 75     typedef std::set<G4LogicalVolume*> ScoringVolumes_t;
 76 
 77   public:
 78     TSDetectorConstruction();
 79     virtual ~TSDetectorConstruction();
 80 
 81     static TSDetectorConstruction* Instance();
 82 
 83   public:
 84     G4VPhysicalVolume* Construct();
 85     inline const G4ThreeVector& GetWorldDimensions() const { return fWorldDim; }
 86     inline const ScoringVolumes_t& GetScoringVolumes() const { return fScoringVolumes; }
 87     inline const G4String& GetMFDName() const { return fMfdName; }
 88     inline G4int GetTotalTargets() const
 89     {
 90       return fTargetSections.x() * fTargetSections.y() * fTargetSections.z();
 91     }
 92 
 93   protected:
 94     virtual MaterialCollection_t ConstructMaterials();
 95     virtual G4VPhysicalVolume* ConstructWorld(const MaterialCollection_t&);
 96     virtual void ConstructSDandField();
 97 
 98   private:
 99     static TSDetectorConstruction* fgInstance;
100     G4VPhysicalVolume* fWorldPhys;
101     ScoringVolumes_t fScoringVolumes;
102     G4String fWorldMaterialName;
103     G4String fTargetMaterialName;
104     G4String fCasingMaterialName;
105     G4ThreeVector fWorldDim;
106     G4ThreeVector fTargetDim;
107     G4ThreeVector fTargetSections;
108     G4String fMfdName;
109 };
110 
111 #endif
112