Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | Static Private Attributes | List of all members
ExGflash2ParallelWorld Class Reference

#include <Doxymodules_parameterisations.h>

Inheritance diagram for ExGflash2ParallelWorld:
G4VUserParallelWorld

Public Member Functions

 ExGflash2ParallelWorld (G4String aWorldName)
 
 ~ExGflash2ParallelWorld () override
 
void Construct () final
 
void ConstructSD () final
 

Private Attributes

G4RegionfRegion {nullptr}
 

Static Private Attributes

static G4ThreadLocal GFlashShowerModelfFastShowerModel = nullptr
 
static G4ThreadLocal GFlashHomoShowerParameterisationfParameterisation = nullptr
 
static G4ThreadLocal GFlashParticleBoundsfParticleBounds = nullptr
 
static G4ThreadLocal GFlashHitMakerfHitMaker = nullptr
 

Detailed Description

Definition at line 123 of file Doxymodules_parameterisations.h.

Constructor & Destructor Documentation

◆ ExGflash2ParallelWorld()

ExGflash2ParallelWorld::ExGflash2ParallelWorld ( G4String  aWorldName)

Definition at line 60 of file ExGflash2ParallelWorld.cc.

61 : G4VUserParallelWorld(aWorldName)
62{
63 G4cout << "ExGflash2ParallelWorld::Parralel world constructor" << G4endl;
64}

◆ ~ExGflash2ParallelWorld()

ExGflash2ParallelWorld::~ExGflash2ParallelWorld ( )
override

Definition at line 68 of file ExGflash2ParallelWorld.cc.

69{
70 delete fFastShowerModel;
71 delete fParameterisation;
72 delete fParticleBounds;
73 delete fHitMaker;
74}
static G4ThreadLocal GFlashHomoShowerParameterisation * fParameterisation
static G4ThreadLocal GFlashParticleBounds * fParticleBounds
static G4ThreadLocal GFlashHitMaker * fHitMaker
static G4ThreadLocal GFlashShowerModel * fFastShowerModel

Member Function Documentation

◆ Construct()

void ExGflash2ParallelWorld::Construct ( )
final

Definition at line 78 of file ExGflash2ParallelWorld.cc.

79{
80 // In parallel world material does not matter
81 G4Material* dummy = nullptr;
82
83 // Build parallel/ghost geometry:
84 auto ghostLogicalVolume = GetWorld()->GetLogicalVolume();
85
86 // Use part of the Ex1GflashDetectorConstruction (without individual crystals)
87
88 G4int nbOfCrystals = 10; // this are the crystals PER ROW in this example
89 // cube of 10 x 10 crystals
90 // don't change it @the moment, since
91 // the readout in event action assumes this
92 // dimensions and is not automatically adapted
93 // in this version of the example :-(
94 // Simplified `CMS-like` PbWO4 crystal calorimeter
95 G4double calo_xside = 31 * cm;
96 G4double calo_yside = 31 * cm;
97 G4double calo_zside = 24 * cm;
98
99 G4double crystalWidth = 3 * cm;
100 G4double crystalLength = 24 * cm;
101
102 calo_xside = (crystalWidth * nbOfCrystals) + 1 * cm;
103 calo_yside = (crystalWidth * nbOfCrystals) + 1 * cm;
104 calo_zside = crystalLength;
105
106 auto calo_box = new G4Box("CMS calorimeter", // its name
107 calo_xside / 2., // size
108 calo_yside / 2., calo_zside / 2.);
109 auto calo_log = new G4LogicalVolume(calo_box, // its solid
110 dummy, // its material
111 "calo log", // its name
112 nullptr, // opt: fieldManager
113 nullptr, // opt: SensitiveDetector
114 nullptr); // opt: UserLimit
115
116 G4double xpos = 0.0;
117 G4double ypos = 0.0;
118 G4double zpos = 100.0 * cm;
119 new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), calo_log, "calorimeter",
120 ghostLogicalVolume, false, 1);
121
122 auto caloVisAtt = new G4VisAttributes(G4Colour(1.0, 1.0, 1.0));
123 calo_log->SetVisAttributes(caloVisAtt);
124
125 // define the fParameterisation region
126 fRegion = new G4Region("crystals");
127 calo_log->SetRegion(fRegion);
128 fRegion->AddRootLogicalVolume(calo_log);
129}

◆ ConstructSD()

void ExGflash2ParallelWorld::ConstructSD ( )
final

Definition at line 133 of file ExGflash2ParallelWorld.cc.

134{
135 // Get nist material manager
136 G4NistManager* nistManager = G4NistManager::Instance();
137 G4Material* pbWO4 = nistManager->FindOrBuildMaterial("G4_PbWO4");
138 // -- fast simulation models:
139 // **********************************************
140 // * Initializing shower modell
141 // ***********************************************
142 G4cout << "Creating shower parameterization models" << G4endl;
143 fFastShowerModel = new GFlashShowerModel("fFastShowerModel", fRegion);
145 fFastShowerModel->SetParameterisation(*fParameterisation);
146 // Energy Cuts to kill particles:
148 fFastShowerModel->SetParticleBounds(*fParticleBounds);
149 // Makes the EnergieSpots
151 fFastShowerModel->SetHitMaker(*fHitMaker);
152 G4cout << "end shower parameterization." << G4endl;
153 // **********************************************
154}

Member Data Documentation

◆ fRegion

G4Region* ExGflash2ParallelWorld::fRegion {nullptr}
private

Definition at line 55 of file ExGflash2ParallelWorld.hh.

55{nullptr};

◆ fFastShowerModel

G4ThreadLocal GFlashShowerModel* ExGflash2ParallelWorld::fFastShowerModel = nullptr
inlinestaticprivate

Definition at line 57 of file ExGflash2ParallelWorld.hh.

◆ fParameterisation

G4ThreadLocal GFlashHomoShowerParameterisation* ExGflash2ParallelWorld::fParameterisation = nullptr
inlinestaticprivate

Definition at line 58 of file ExGflash2ParallelWorld.hh.

◆ fParticleBounds

G4ThreadLocal GFlashParticleBounds* ExGflash2ParallelWorld::fParticleBounds = nullptr
inlinestaticprivate

Definition at line 59 of file ExGflash2ParallelWorld.hh.

◆ fHitMaker

G4ThreadLocal GFlashHitMaker* ExGflash2ParallelWorld::fHitMaker = nullptr
inlinestaticprivate

Definition at line 60 of file ExGflash2ParallelWorld.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration