Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
XrayFluoMercuryDetectorConstruction.hh
Go to the documentation of this file.
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 //
27 // $Id: XrayFluoMercuryDetectorConstruction.hh
28 // GEANT4 tag $Name:XrayFluo-V05-02-06
29 //
30 // Author: Alfonso Mantero (Alfonso.Mantero@ge.infn.it)
31 //
32 // History:
33 // -----------
34 //
35 // 08 Sep 2003 Alfonso Mantero created
36 //
37 // -------------------------------------------------------------------
38 
39 #ifndef XrayFluoMercuryDetectorConstruction_hh
40 #define XrayFluoMercuryDetectorConstruction_hh 1
41 
42 #include "globals.hh"
44 #include "G4RotationMatrix.hh"
45 #include "G4Cache.hh"
48 #include "XrayFluoSD.hh"
49 
50 class G4Box;
51 class G4Tubs;
52 class G4Sphere;
53 class G4LogicalVolume;
54 class G4VPhysicalVolume;
55 class G4Material;
58 
59 //class XrayFluoSD;
60 //class XrayFluoVDetectorType;
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63 
65 {
66 public:
67 
68 
70 
71 public:
72 
74 
75  void ConstructSDandField();
76 
77  void UpdateGeometry();
78 
79 
80  void SetMercuryMaterial(G4String newMaterial);
81 
82  void SetDetectorType(G4String type);
83 
85 
86  inline void SetOribitHeight(G4double size)
87  {distDe = size;};
88  inline void SetLatitude(G4double lat)
89  {ThetaHPGe = 180.* CLHEP::deg - lat;};
90 
92 
94 
95 
96  G4double GetWorldSizeZ() const {return WorldSizeZ;};
97  G4double GetWorldSizeXY() const {return WorldSizeXY;};
98 
99  G4double GetDeviceThickness() const {return DeviceThickness;};
100  G4double GetDeviceSizeX() const {return DeviceSizeX;};
101  G4double GetDeviceSizeY() const {return DeviceSizeY;};
102  G4double GetPixelSizeXY() const {return PixelSizeXY;};
103  G4double GetContactSizeXY() const {return ContactSizeXY;};
104 
105  G4int GetNbOfPixels() const {return NbOfPixels;}; //mandatory for XrayFluoSD
106  G4int GetNbOfPixelRows() const {return NbOfPixelRows;};
107  G4int GetNbOfPixelColumns() const {return NbOfPixelColumns;};
108 
109  G4Material* GetOhmicPosMaterial() const {return OhmicPosMaterial;};
110  G4double GetOhmicPosThickness() const {return OhmicPosThickness;};
111 
112  G4Material* GetOhmicNegMaterial() const {return OhmicNegMaterial;};
113  G4double GetOhmicNegThickness() const {return OhmicNegThickness;};
114 
115  const G4VPhysicalVolume* GetphysiWorld() const {return physiWorld;};
116  const G4VPhysicalVolume* GetHPGe() const {return physiHPGe;};
117  const G4VPhysicalVolume* GetMercury() const {return physiMercury;};
118 
119  const G4VPhysicalVolume* GetphysiPixel() const {return physiPixel;};
120  const G4VPhysicalVolume* GetOhmicPos() const {return physiOhmicPos;};
121  const G4VPhysicalVolume* GetOhmicNeg() const {return physiOhmicNeg;};
122 
123  const G4VPhysicalVolume* GetOptic() const {return physiOptic;};
124 
125 private:
126 
128 
129  static XrayFluoMercuryDetectorConstruction* instance;
130 
131  XrayFluoVDetectorType* detectorType;
132 
133  G4bool mercuryGranularity;
134 
135  G4double DeviceSizeX;
136  G4double DeviceSizeY;
137  G4double DeviceThickness;
138 
139  G4Box* solidWorld; //pointer to the solid World
140  G4LogicalVolume* logicWorld; //pointer to the logical World
141  G4VPhysicalVolume* physiWorld; //pointer to the physical World
142 
143  G4Box* solidHPGe; //pointer to the solid Sensor
144  G4LogicalVolume* logicHPGe; //pointer to the logical Sensor
145  G4VPhysicalVolume* physiHPGe; //pointer to the physical Sensor
146 
147  G4Box* solidScreen; //pointer to the solid Screen
148  G4LogicalVolume* logicScreen; //pointer to the logical Screen
149  G4VPhysicalVolume* physiScreen; //pointer to the physical Screen
150 
151  G4Sphere* solidMercury; //pointer to the solid Mercury
152  G4LogicalVolume* logicMercury; //pointer to the logical Mercury
153  G4VPhysicalVolume* physiMercury; //pointer to the physical Mercury
154 
155  // G4Tubs* solidDia1; //pointer to the solid Diaphragm
156  // G4LogicalVolume* logicDia1; //pointer to the logical Diaphragm
157  // G4VPhysicalVolume* physiDia1; //pointer to the physical Diaphragm
158 
159  // G4Tubs* solidDia3; //pointer to the solid Diaphragm
160  // G4LogicalVolume* logicDia3; //pointer to the logical Diaphragm
161  // G4VPhysicalVolume* physiDia3; //pointer to the physical Diaphragm
162 
163  G4Box* solidOhmicPos;
164  G4LogicalVolume* logicOhmicPos;
165  G4VPhysicalVolume* physiOhmicPos;
166 
167  G4Box* solidOhmicNeg;
168  G4LogicalVolume* logicOhmicNeg;
169  G4VPhysicalVolume* physiOhmicNeg;
170 
171  G4Box* solidPixel;
172  G4LogicalVolume* logicPixel;
173  G4VPhysicalVolume* physiPixel;
174 
175  G4Tubs* solidOptic;
176  G4LogicalVolume* logicOptic;
177  G4VPhysicalVolume* physiOptic;
178 
179  G4LogicalVolume* logicGrain;
180 
181  //materials management
182  XrayFluoNistMaterials* materials;
183 
184  G4Material* screenMaterial;
185  G4Material* OhmicPosMaterial;
186  G4Material* OhmicNegMaterial;
187  G4Material* pixelMaterial;
188  G4Material* mercuryMaterial;
189  // G4Material* Dia3Material;
190  G4Material* defaultMaterial;
191 
192  //apparate parameters
193 
194  G4double OhmicPosThickness;
195  G4double OhmicNegThickness;
196 
197  G4double opticDia;
198  G4double opticThickness;
199 
200  G4double screenSizeXY;
201  G4double screenThickness;
202 
203  G4int PixelCopyNb;
204  G4int grainCopyNb;
205  G4int NbOfPixels;
206  G4int NbOfPixelRows;
207  G4int NbOfPixelColumns;
208  G4double PixelThickness;
209  G4double PixelSizeXY;
210  G4double ContactSizeXY;
211 
212  G4double opticAperture;
213 
214  G4double mercuryDia;
215  G4double sunDia;
216 
217  G4double mercurySunDistance;
218 
219  G4double ThetaHPGe;
220 
221  G4double distDe;
222  G4double distScreen;
223  G4double distOptic;
224 
225 
226  G4double PhiHPGe;
227 
228  G4RotationMatrix zRotPhiHPGe;
229 
230  G4double WorldSizeXY;
231  G4double WorldSizeZ;
232 
233 
234  XrayFluoMercuryDetectorMessenger* detectorMessenger; //pointer to the Messenger
235 
236  G4Cache<XrayFluoSD*> HPGeSD; //pointer to the sensitive detector
237 
238 public:
239 
240  G4Material* GetMercuryMaterial() const {return mercuryMaterial;};
241  G4Material* GetPixelMaterial() const {return pixelMaterial;};
242 
243  G4double GetMercuryDia() const {return mercuryDia;};
244  G4double GetSunDia() const {return sunDia;};
245 
246  //Inclinaton of the orbit respect Mercury respect the equator (latitude)
247 
248  G4double GetOrbitInclination() const {return 180 * CLHEP::deg - ThetaHPGe;};
249  G4double GetOrbitDistance() const {return distDe;};
250  G4double GetOpticAperture() const {return opticAperture;};
251 
252 
253 
254 private:
255 
256  void DefineDefaultMaterials();
257  G4VPhysicalVolume* ConstructApparate();
258 
259  //calculates some quantities used to construct geometry
260  void ComputeApparateParameters();
261 
262 };
263 
264 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
265 
266 inline void XrayFluoMercuryDetectorConstruction::ComputeApparateParameters()
267 {
268  // Compute derived parameters of the apparate
269 
270  DeviceThickness = PixelThickness+OhmicNegThickness+OhmicPosThickness;
271 
273 
274  DeviceSizeY =(NbOfPixelRows * std::max(ContactSizeXY,PixelSizeXY));
275  DeviceSizeX =(NbOfPixelColumns * std::max(ContactSizeXY,PixelSizeXY));
276 
277  screenSizeXY = opticDia;
278 
279  G4cout << "DeviceSizeX(cm): "<< DeviceSizeX/CLHEP::cm <<G4endl;
280  G4cout << "DeviceSizeY(cm): "<< DeviceSizeY/CLHEP::cm << G4endl;
281 
282  //*********************************************************************
283  //** Astronomical distances reduce by a factor 10^-7 due to G4 Bug **
284  //*********************************************************************
285 
286  WorldSizeZ = 2 * mercurySunDistance ;
287  WorldSizeXY = (2 * distDe) + 2000 * CLHEP::km ;
288  //WorldSizeZ = WorldSizeXY;
289 }
290 
291 #endif
static constexpr double cm
Definition: SystemOfUnits.h:99
Definition: G4Box.hh:64
Definition: G4Tubs.hh:85
int G4int
Definition: G4Types.hh:78
static XrayFluoMercuryDetectorConstruction * GetInstance()
static constexpr double km
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
static constexpr double deg
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76