Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/xray_TESdetector/src/XrayTESdetDetectorConstruction.cc

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 //
 27 /// \file XrayTESdetDetectorConstruction.cc
 28 /// \brief Implementation of the XrayTESdetDetectorConstruction class
 29 //
 30 // Authors: P.Dondero (paolo.dondero@cern.ch), R.Stanzani (ronny.stanzani@cern.ch)
 31 //
 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 33 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 34 
 35 
 36 #include "XrayTESdetDetectorConstruction.hh"
 37 #include "XrayTESdetDetectorMessenger.hh"
 38 #include "XrayTESdetDetParameterisation.hh"
 39 #include "G4Material.hh"
 40 #include "G4Sphere.hh"
 41 #include "G4Box.hh"
 42 #include "G4Tubs.hh"
 43 #include "G4LogicalVolume.hh"
 44 #include "G4ThreeVector.hh"
 45 #include "G4PVPlacement.hh"
 46 #include "G4PVParameterised.hh"
 47 #include "globals.hh"
 48 #include "G4SDManager.hh"
 49 #include "G4GeometryTolerance.hh"
 50 #include "G4GeometryManager.hh"
 51 #include "G4UserLimits.hh"
 52 #include "G4VisAttributes.hh"
 53 #include "G4Colour.hh"
 54 #include "G4ios.hh"
 55 #include "G4SubtractionSolid.hh"
 56 #include "G4Cons.hh"
 57 #include "G4RotationMatrix.hh"
 58 #include "G4PVReplica.hh"
 59 #include "G4NistManager.hh"
 60 #include "G4Cons.hh"
 61 #include "G4Polyhedra.hh"
 62 #include "G4PVParameterised.hh"
 63 #include "G4Trap.hh"
 64 #include "G4Trd.hh"
 65 #include "G4Hype.hh"
 66 #include "G4ExtrudedSolid.hh"
 67 #include "G4SystemOfUnits.hh"
 68 #include "G4PhysicalConstants.hh"
 69 
 70 // #include <stdio.h>
 71 #include <cmath>
 72 #include <iostream>
 73 #include <fstream>
 74 #include <iomanip>
 75 #include <string>
 76 #include <vector>
 77 #include <algorithm>
 78 #include <sstream>
 79 
 80 #include "G4GDMLParser.hh"
 81 #include "G4PhysicalVolumeStore.hh"
 82 #include "G4ProductionCuts.hh"
 83 #include "G4ProductionCutsTable.hh"
 84 
 85 
 86 XrayTESdetDetectorConstruction::XrayTESdetDetectorConstruction()
 87 : fStepLimitPolyfilta(nullptr),fStepLimitAlfilta(nullptr),fStepLimitmembrane(nullptr),fExperimentalHall_phys(nullptr),fDetectorMessenger(nullptr)
 88 {
 89   fReadFile = "";
 90   fDetectorMessenger = new XrayTESdetDetectorMessenger(this);
 91 }
 92 
 93 
 94 XrayTESdetDetectorConstruction::~XrayTESdetDetectorConstruction()
 95 {
 96   delete fStepLimitAlfilta;
 97   delete fStepLimitPolyfilta;
 98   delete fStepLimitmembrane;
 99   delete fDetectorMessenger;
100 }
101 
102 
103 G4VPhysicalVolume* XrayTESdetDetectorConstruction::Construct()
104 {
105   return ConstructDetector();
106 }
107 
108 
109 G4VPhysicalVolume* XrayTESdetDetectorConstruction::ConstructDetector()
110 {
111   G4VPhysicalVolume* fWorld_phys;
112   auto *InnerRegion = new G4Region("InnerRegion");
113 
114   if (fReadFile == "")
115   {
116     G4cout << "Build a new geometry" << G4endl;
117 
118     // -------------------------------------------------------------------------------------------
119     // --------------- Define materials ----------------------------------------------------------
120     // -------------------------------------------------------------------------------------------
121     G4String symbol;
122     G4int ncomponents;
123     G4double density;
124     G4double fractionmass;
125 
126     auto *man = G4NistManager::Instance();
127     auto *Si = man->FindOrBuildMaterial("G4_Si");
128     auto *Cu = man->FindOrBuildMaterial("G4_Cu");
129     auto *Bi = man->FindOrBuildMaterial("G4_Bi");
130     auto *Al = man->FindOrBuildMaterial("G4_Al");
131     auto *Ni = man->FindOrBuildMaterial("G4_Ni");
132     auto *Au = man->FindOrBuildMaterial("G4_Au");
133     auto *N = man->FindOrBuildMaterial("G4_N");
134     auto *Fe = man->FindOrBuildMaterial("G4_Fe");
135     auto *C = man->FindOrBuildMaterial("G4_C");
136     auto *Mn = man->FindOrBuildMaterial("G4_Mn");
137     auto *Ph = man->FindOrBuildMaterial("G4_P");
138     auto *S = man->FindOrBuildMaterial("G4_S");
139     auto *Cr = man->FindOrBuildMaterial("G4_Cr");
140     auto *Mo = man->FindOrBuildMaterial("G4_Mo");
141     auto *Ti = man->FindOrBuildMaterial("G4_Ti");
142 
143     // Define the material of the ring
144     auto *RingAl = new G4Material("RingAl", density=0.736*g/cm3, ncomponents=1);
145     RingAl->AddMaterial(Al, fractionmass=1);
146 
147     // Define cryoperm
148     auto *Cryoperm = new G4Material("Cryoperm", density=8.7*g/cm3, ncomponents=4);
149     Cryoperm->AddMaterial(Ni, fractionmass=0.775);
150     Cryoperm->AddMaterial(Cu, fractionmass=0.045);
151     Cryoperm->AddMaterial(Mo, fractionmass=0.025);
152     Cryoperm->AddMaterial(Fe, fractionmass=0.155);
153 
154     // Define vacuum
155     auto *vacuum = man->FindOrBuildMaterial("G4_Galactic");
156 
157     // Define the Stainless Steel
158     G4double SS_density = 8.03*g/cm3;
159     auto *SS = new G4Material("Stainless Steel",density=SS_density,ncomponents=8);
160     SS->AddMaterial(N, fractionmass=0.01);
161     SS->AddMaterial(Ni, fractionmass=0.1);
162     SS->AddMaterial(Cr, fractionmass=0.185);
163     SS->AddMaterial(S, fractionmass=0.03);
164     SS->AddMaterial(Ph, fractionmass=0.045);
165     SS->AddMaterial(Mn, fractionmass=0.02);
166     SS->AddMaterial(C, fractionmass=0.08);
167     SS->AddMaterial(Fe, fractionmass=0.53);
168 
169     // Si3N4
170     G4double A;
171     G4double Z;
172     auto *elSi = new G4Element("Silicon", "Si", Z=14., A=28.0855*g/mole);
173     auto *elN = new G4Element("Nitrogen","N", Z=7., A=14.00674*g/mole);
174     G4double Si3N4density = 3.44*g/cm3;
175     auto *Si3N4 = new G4Material("Si3N4", Si3N4density, 2);
176     Si3N4->AddElement (elSi, 3);
177     Si3N4->AddElement (elN, 4);
178 
179     // -------------------------------------------------------------------------------------------
180     // ----------------------------------- Define World volume -----------------------------------
181     G4double expHall_x = 3*m;
182     G4double expHall_y = 3*m;
183     G4double expHall_z = 3*m;
184 
185     auto *experimentalHall_box = new G4Box("experimentalHall_box",expHall_x,expHall_y,expHall_z);
186     fExperimentalHall_log = new G4LogicalVolume(experimentalHall_box,vacuum,"experimentalHall_log",0,0,0);// vacuum ok, its the mother volume
187     fExperimentalHall_phys = new G4PVPlacement(nullptr,G4ThreeVector(), fExperimentalHall_log,"expHall",0,false,0);
188 
189     // -------------------------------------------------------------------------------------------
190     // ------------------------------------ Define TES array -------------------------------------
191     //absorber: 3 µm Bi
192     G4double Biabsorberthickness = 3*um;
193     G4double Gridthickness = 200*um;
194     G4double membrane_thickness = 0.75*um;
195     G4double wafer_thickness = 0.25*mm;
196 
197     //height of the entire volume to be replicated
198     G4double xifu_z = Biabsorberthickness+Gridthickness+membrane_thickness;
199 
200     //size of the mother volume to be replicated
201     G4double pxl_pitch = 0.814*mm;//V2 res90
202     G4double pxl_gap = 11*um;
203     G4double pxl_size = pxl_pitch-pxl_gap;
204     G4double grid_size = pxl_pitch;
205 
206     G4double DetPos_x = 0.0*mm;
207     G4double DetPos_y = 0.0*mm;
208     G4double DetPos_z = -xifu_z; //So that the absorbers surface is at z=0
209 
210     //mother volume of the detector has hexagonal shape
211     G4double phiStart = 0;
212     G4double phiTotal = 0;
213     G4int numSide = 6;
214     G4int numZPlanes = 2;
215     G4double rInner[2] = {0,0};
216 
217     //end of parameters definitions
218     auto *element_box = new G4Box("element_box",grid_size*0.5,grid_size*0.5,xifu_z*0.5);
219     fElement_log = new G4LogicalVolume(element_box,vacuum,"element_log");
220 
221     //definition of the mother volume containing the replicated elements
222     G4String pName = "Detector_hex";
223     G4double zPlane[2] = {-xifu_z*0.5,xifu_z*0.5};
224 
225     G4double rOuter[2] = {8.5*mm,8.5*mm};//apothem
226 
227     G4Polyhedra* Detector_hex = new G4Polyhedra(pName,
228     phiStart,
229     phiTotal,
230     numSide,
231     numZPlanes,
232     zPlane,
233     rInner,
234     rOuter);
235 
236     fDetector_log = new G4LogicalVolume(Detector_hex,vacuum,"Detector_log");//vacuum ok, it's the mother volume
237     fDetector_phys = new G4PVPlacement(nullptr, G4ThreeVector(DetPos_x,DetPos_y,DetPos_z+0.0001*um),fDetector_log,"Detector",fExperimentalHall_log, false, 0);
238 
239     //Bi absorber 3 um
240     auto *Bipxl_box = new G4Box("Bipxl_box",pxl_size*0.5,pxl_size*0.5,Biabsorberthickness*0.5);
241     fBipxl_log = new G4LogicalVolume(Bipxl_box,Bi,"Bipxl_log");
242     fBipxl_phys = new G4PVPlacement(nullptr,G4ThreeVector(0., 0., xifu_z*0.5-Biabsorberthickness*0.5), fBipxl_log, "Bipxl", fElement_log, false, 0);
243 
244     //SiN membrane sustaing the absorbers
245     auto *membranepxl_box = new G4Box("membranepxl_box",pxl_size*0.5,pxl_size*0.5,membrane_thickness*0.5);
246     fMembranepxl_log = new G4LogicalVolume(membranepxl_box,Si3N4,"membranepxl_log");
247     fMembranepxl_phys = new G4PVPlacement(nullptr,G4ThreeVector(0., 0., xifu_z*0.5-Biabsorberthickness-membrane_thickness*0.5), fMembranepxl_log, "membrane", fElement_log, false, 0);
248 
249     //step limiter: the membrane is really thin, thinner than the typical interaction step. This forces the particles to perform at least 5 interactions inside
250     G4double maxStepmembrane = 0.2*membrane_thickness;
251     fStepLimitmembrane = new G4UserLimits(maxStepmembrane);
252     fMembranepxl_log->SetUserLimits(fStepLimitmembrane);
253 
254     // Grid element
255     auto *gridbox = new G4Box("grid_box",grid_size*0.5,grid_size*0.5,Gridthickness*0.5);
256     G4double gridbeam_size = 165*um;
257     G4double gridhole_size = grid_size-gridbeam_size;
258     G4double gridhole_thickness = Gridthickness+1*mm;
259     auto *gridhole_box = new G4Box("gridhole_box",gridhole_size*0.5,gridhole_size*0.5,gridhole_thickness*0.5);
260     G4SubtractionSolid* gridsubtraction = new G4SubtractionSolid("gridsubtraction",gridbox, gridhole_box);
261     fGridpiece_log = new G4LogicalVolume(gridsubtraction,Si,"gridpiece_log");
262     fGridpiece_phys = new G4PVPlacement(nullptr,G4ThreeVector(0., 0., xifu_z*0.5-Biabsorberthickness-membrane_thickness-Gridthickness*0.5), fGridpiece_log, "gridpiece", fElement_log, false, 0);
263 
264     //parametrization of element_log
265     auto *DetParam = new XrayTESdetDetParameterisation(
266                 317,//V1 and V2 res80 and res90 // to be commented for overlap check
267                 //2,  // To enable volumes overlap check the number of replicas must be reduced. Keep it commented if not needed
268                 -grid_size*38,
269                 grid_size,
270                 grid_size*0.5,
271                 grid_size );
272 
273     fPhysiDet = new G4PVParameterised(
274                   "Det",             // their name
275                   fElement_log,       // their logical volume
276                   fDetector_log,      // Mother logical volume
277                   kXAxis,            // Are placed along this axis
278                   317,               //V1 and V2 res80 and res90
279                   //2,               // overlap check. Same as before.
280                   DetParam);         // The parametrisation
281 
282     // ------------ Define support wafer ----------------------------------------------
283     G4double rOutersupp[2] = {39.4*mm,39.4*0.5*mm};
284     G4double zwaferPlane[2] = {-0.5*wafer_thickness,0.5*wafer_thickness};
285 
286     auto *wafer_hex = new G4Polyhedra("wafer_hex",
287     phiStart,
288     phiTotal,
289     numSide,
290     numZPlanes,
291     zwaferPlane,
292     rOuter,
293     rOutersupp);
294 
295     fWafer_log = new G4LogicalVolume(wafer_hex,Si,"wafer_log");
296     G4double wafer_z = -xifu_z*0.5-Biabsorberthickness-Gridthickness*0.5;
297     fWafer_phys = new G4PVPlacement(nullptr, G4ThreeVector(DetPos_x,DetPos_y,wafer_z), fWafer_log,"wafer",fExperimentalHall_log, false, 0);
298 
299     //------------------------------------------------------------------------------------------------------
300     //----------------------------------   BSC electron detector     ---------------------------------------
301     // This fake surface placed just above the absorbers gathers information on the electrons crossing it
302     // it can be used to check if the same particle crossed it twice, to identify backscattering events
303     //------------------------------------------------------------------------------------------------------
304     G4double BSC_z = 0.05*mm;  // half thickness
305 
306     G4String pBSCName = "BSCdetector";
307     G4double zPlaneBSC[2] = {-BSC_z,BSC_z};
308     G4double rOuterBSC[2] = {rOuter[0]+2*mm,rOuter[1]+2*mm};
309 
310     auto *BSC_hex = new G4Polyhedra(pBSCName,
311       phiStart,
312       phiTotal,
313       numSide,
314       numZPlanes,
315       zPlaneBSC,
316       rInner,
317       rOuterBSC);
318 
319     fBSC_log = new G4LogicalVolume(BSC_hex,vacuum,"BSC_log");
320 
321     G4double BSCPos_x = 0.0*mm;
322     G4double BSCPos_y = 0.0*mm;
323     G4double BSCPos_z = wafer_z + wafer_thickness+BSC_z*0.5;
324 
325     fBSC_phys = new G4PVPlacement(nullptr, G4ThreeVector(BSCPos_x,BSCPos_y,BSCPos_z),fBSC_log,"BSC",fExperimentalHall_log, false, 0);
326 
327     // --------------------------------------------------------------------------------------------------------------------------------------------
328     // -------------------------------------------------------- Define ACD inner array ------------------------------------------------------------
329     // Anti Coincidence Detectors are usually thicker than the main ones. It can generate a huge output file if low cuts are used
330     // defining a inner part of the ACD and assigning it to a lower precision region, thus allowing only the creation of high energy secondaries
331     // can help speeding up the simulation and reducing the output size. In this example only the external 25 um belong to the higher precision region
332     // --------------------------------------------------------------------------------------------------------------------------------------------
333     G4double ACDthickness = 200*um;
334     G4double ACD_distance = 0.75*mm; // distance from the lower edge of the main detector
335     G4double ACDarrayPos_z = -1*(ACD_distance+Biabsorberthickness/2+ACDthickness/2);
336     G4double ACDpxlGAP = 50*um;
337 
338     G4double pxlheight = 8.5*mm;
339     G4double pxlwideside = pxlheight*1.191*mm;
340     G4double pxlshortside = pxlwideside*0.5;
341 
342     G4double inACDthickness = ACDthickness-50*um;
343 
344     G4double inpxlheight = pxlheight-50*um;
345     G4double inpxlwideside = pxlwideside-50*um;
346     G4double inpxlshortside = inpxlwideside*0.5;
347 
348     auto *inACDpxl = new G4Trap("inACDpxl",
349       inACDthickness,
350       inpxlheight,
351       inpxlwideside,
352       inpxlshortside);
353 
354     fInACDpxl_log = new G4LogicalVolume(inACDpxl,Si,"inACDpxl_log");
355     G4RotationMatrix rotm = G4RotationMatrix();
356 
357     //pxl1
358     G4double ACDarrayPos1_x = (pxlshortside*1.5+ACDpxlGAP)*0.5;
359     G4double ACDarrayPos1_y = (pxlheight+ACDpxlGAP)*0.5;
360     fInACDpxl_phys = new G4PVPlacement(nullptr, G4ThreeVector(ACDarrayPos1_x,ACDarrayPos1_y,ACDarrayPos_z),fInACDpxl_log,"pxl1",fExperimentalHall_log, false, 0);
361 
362     //pxl2
363     rotm.rotateZ(180*deg);
364     G4double ACDarrayPos2_x = -ACDarrayPos1_x;
365     G4double ACDarrayPos2_y = -ACDarrayPos1_y;
366     G4ThreeVector position2 = G4ThreeVector(ACDarrayPos2_x,ACDarrayPos2_y,ACDarrayPos_z);
367     G4Transform3D transform2 = G4Transform3D(rotm,position2);
368     fInACDpxl_phys = new G4PVPlacement(transform2,fInACDpxl_log,"pxl2",fExperimentalHall_log, false, 0);
369 
370     //pxl3
371     rotm.rotateY(180*deg);
372     G4double ACDarrayPos3_x = ACDarrayPos1_x;
373     G4double ACDarrayPos3_y = -ACDarrayPos1_y;
374     G4ThreeVector position3 = G4ThreeVector(ACDarrayPos3_x,ACDarrayPos3_y,ACDarrayPos_z);
375     G4Transform3D transform3 = G4Transform3D(rotm,position3);
376     fInACDpxl_phys = new G4PVPlacement(transform3,fInACDpxl_log,"pxl3",fExperimentalHall_log, false, 0);
377 
378     //pxl4
379     rotm.rotateZ(180*deg);
380     G4double ACDarrayPos4_x = -ACDarrayPos1_x;
381     G4double ACDarrayPos4_y = ACDarrayPos1_y;
382     G4ThreeVector position4 = G4ThreeVector(ACDarrayPos4_x,ACDarrayPos4_y,ACDarrayPos_z);
383     G4Transform3D transform4 = G4Transform3D(rotm,position4);
384     fInACDpxl_phys = new G4PVPlacement(transform4,fInACDpxl_log,"pxl4",fExperimentalHall_log, false, 0);
385 
386     // ----------------------------------------------------------------------------------
387     // -------------------- Define ACD outer layer array --------------------------------
388     auto *ACDpxl = new G4Trap("ACDpxl",
389       ACDthickness,
390       pxlheight,
391       pxlwideside,
392       pxlshortside);
393     G4SubtractionSolid* ACDsub = new G4SubtractionSolid("ACDsub",ACDpxl, inACDpxl);
394     fACDpxl_log = new G4LogicalVolume(ACDsub,Si,"ACDpxl_log");
395 
396     //pxl1
397     fACDpxl_phys = new G4PVPlacement(nullptr, G4ThreeVector(ACDarrayPos1_x,ACDarrayPos1_y,ACDarrayPos_z),fACDpxl_log,"pxl1a",fExperimentalHall_log, false, 0);
398 
399     //pxl2
400     fACDpxl_phys = new G4PVPlacement(transform2,fACDpxl_log,"pxl2a",fExperimentalHall_log, false, 0);
401 
402     //pxl3
403     fACDpxl_phys = new G4PVPlacement(transform3,fACDpxl_log,"pxl3a",fExperimentalHall_log, false, 0);
404 
405     //pxl4
406     fACDpxl_phys = new G4PVPlacement(transform4,fACDpxl_log,"pxl4a",fExperimentalHall_log, false, 0);
407 
408     //------------------------------------------------------------------------------------------------------
409     //-----------------------------------------SUPPORTS------------------------------------------------//
410 
411     //---------------------   ACD plate     --------------
412     G4double ACplate_thickness = 3.5*mm;  // support thickness
413 
414     G4String ACplateName = "ACplateName";
415     G4double ACplatezPlane[2] = {-0.5*ACplate_thickness,0.5*ACplate_thickness};
416     G4double ACplaterInner[2] = {pxlwideside,pxlwideside};//apothem
417     G4double ACplaterOuter[2] = {43*mm,43*mm};
418 
419     auto *ACplate_hex = new G4Polyhedra(ACplateName,
420     phiStart,
421     phiTotal,
422     numSide,
423     numZPlanes,
424     ACplatezPlane,
425     ACplaterInner,
426     ACplaterOuter );
427 
428     fACplate_log = new G4LogicalVolume(ACplate_hex,Cu,"ACplate_log");//tbc
429 
430     G4double ACplatePos_x = DetPos_x;
431     G4double ACplatePos_y = DetPos_y;
432     G4double ACplatePos_z = -2.25*mm;
433 
434     fACplate_phys = new G4PVPlacement(nullptr, G4ThreeVector(ACplatePos_x,ACplatePos_y,ACplatePos_z),fACplate_log,"ACDplate",fExperimentalHall_log, false, 0);
435 
436     //------------------------------------------------------------------------------------------------------
437     //----------------------------------------------   cage   ----------------------------------------------
438     // ALL THE FOLLOWING SUPPORTING STRUCTURES DEFINE A COMPLEX SHAPE THAT WAS SIMPLIFIED AND ADAPTED FROM AN ENGINEERING CAD MODEL
439     // ALL THE NUMBERS RELATIVE TO SIZES AND POSITIONS HAVE BEEN EXTRACTED FROM THE CAD FILE, THAT IS WHY THE NUMBERS ARE SO ODD
440     // THIS IS MEANT JUST TO GIVE AN EXAMPLE OF THE SUPPORTING STRUCTURES PRESENT IN AN ACTUAL CRYOSTAT
441 
442     //------------------supporting columns (full)
443     G4double hightOfcagecolumn = 60*mm;
444 
445     auto *cagecolumn = new G4Tubs("cagecolumn", 0*mm, 3.99*mm, hightOfcagecolumn*0.5, 0*deg, 360*deg);
446     G4LogicalVolume* cagecolumn_log = new G4LogicalVolume(cagecolumn,Cu,"cagecolumn_log",0,0,0);
447 
448     fCagecolumn_phys = new G4PVPlacement(
449       0,                                        // no rotation
450       G4ThreeVector(45.75*mm,0*mm,-36.5027*mm), // translation position
451       cagecolumn_log,                           // its logical volume
452       "cagecolumn1",                             // its name
453       fExperimentalHall_log,                      // its mother (logical) volume
454       false,                                    // no boolean operations
455       0);
456 
457     fCagecolumn2_phys = new G4PVPlacement(
458       0,
459       G4ThreeVector(-45.75*mm,0*mm,-36.5027*mm),
460       cagecolumn_log,
461       "cagecolumn2",
462       fExperimentalHall_log,
463       false,
464       0);
465 
466     fCagecolumn3_phys = new G4PVPlacement(
467       0,
468       G4ThreeVector(22.875*mm,39.621*mm,-36.5027*mm),
469       cagecolumn_log,
470       "cagecolumn3",
471       fExperimentalHall_log,
472       false,
473       0);
474 
475     fCagecolumn4_phys = new G4PVPlacement(
476       0,
477       G4ThreeVector(22.875*mm,-39.621*mm,-36.5027*mm),
478       cagecolumn_log,
479       "cagecolumn4",
480       fExperimentalHall_log,
481       false,
482       0);
483 
484     fCagecolumn5_phys = new G4PVPlacement(
485       0,
486       G4ThreeVector(-22.875*mm,-39.621*mm,-36.5027*mm),
487       cagecolumn_log,
488       "cagecolumn5",
489       fExperimentalHall_log,
490       false,
491       0);
492 
493     fCagecolumn6_phys = new G4PVPlacement(
494       0,
495       G4ThreeVector(-22.875*mm,39.621*mm,-36.5027*mm),
496       cagecolumn_log,
497       "cagecolumn6",
498       fExperimentalHall_log,
499       false,
500       0);
501 
502     //----------------------------------supporting columns (empty)
503     auto *emptycagecolumn = new G4Tubs("emptycagecolumn", 1.25*mm, 3*mm, hightOfcagecolumn*0.5, 0*deg, 360*deg);
504     G4LogicalVolume* emptycagecolumn_log = new G4LogicalVolume(emptycagecolumn,Cu,"emptycagecolumn_log",0,0,0);
505 
506     fEmptycagecolumn_phys = new G4PVPlacement(
507       0,
508       G4ThreeVector(22.875*mm,0*mm,-36.5027*mm),
509       emptycagecolumn_log,
510       "emptycagecolumn1",
511       fExperimentalHall_log,
512       false,
513       0);
514 
515     fEmptycagecolumn2_phys = new G4PVPlacement(
516       0,
517       G4ThreeVector(-11.4375*mm,-19.8103311*mm,-36.5027*mm),
518       emptycagecolumn_log,
519       "emptycagecolumn2",
520       fExperimentalHall_log,
521       false,
522       0);
523 
524     fEmptycagecolumn3_phys = new G4PVPlacement(
525       0,
526       G4ThreeVector(-11.4375*mm,19.810331*mm,-36.5027*mm),
527       emptycagecolumn_log,
528       "emptycagecolumn3",
529       fExperimentalHall_log,
530       false,
531       0);
532 
533     //-----------------------walls
534     G4double cagewall_x = 1*mm;
535     G4double cagewall_y = 32.62*mm;
536     G4double cagewall_z = 65*mm;
537 
538     auto *cagewall_box = new G4Box("_box",cagewall_x*0.5,cagewall_y*0.5,cagewall_z*0.5);
539 
540     fCagewall_log = new G4LogicalVolume(cagewall_box,Cu,"cagewall_log",0,0,0);  //tbc
541     fCagewall_phys = new G4PVPlacement(
542       0,
543       G4ThreeVector(22.875*mm,19.3103311*mm,-36.5027*mm),
544       fCagewall_log,
545       "cagewall1",
546       fExperimentalHall_log,
547       false,
548       0);
549 
550     ///and the other 5 pieces
551     fCagewall2_phys = new G4PVPlacement(
552       0,
553       G4ThreeVector(22.875*mm,-19.3103311*mm,-36.5027*mm),
554       fCagewall_log,
555       "cagewall2",
556       fExperimentalHall_log,
557       false,
558       0);
559 
560     G4RotationMatrix wallrotm = G4RotationMatrix();
561 
562     wallrotm.rotateZ(-60*deg);
563     G4ThreeVector wallposition3 = G4ThreeVector(5.2857373*mm,29.4654967*mm,-36.5027*mm);
564     G4Transform3D walltransform3 = G4Transform3D(wallrotm,wallposition3);
565     fCagewall3_phys = new G4PVPlacement(walltransform3,fCagewall_log,"cagewall3",fExperimentalHall_log, true, 0);
566 
567     G4ThreeVector wallposition4 = G4ThreeVector(-28.1607373*mm,10.1551656*mm,-36.5027*mm);
568     G4Transform3D walltransform4 = G4Transform3D(wallrotm,wallposition4);
569     fCagewall4_phys = new G4PVPlacement(walltransform4,fCagewall_log,"cagewall4",fExperimentalHall_log, true, 0);
570 
571     wallrotm.rotateZ(120*deg);
572     G4ThreeVector wallposition5 = G4ThreeVector(-28.1607373*mm,-10.1551656*mm,-36.5027*mm);
573     G4Transform3D walltransform5 = G4Transform3D(wallrotm,wallposition5);
574     fCagewall5_phys = new G4PVPlacement(walltransform5,fCagewall_log,"cagewall5",fExperimentalHall_log, true, 0);
575 
576     G4ThreeVector wallposition6 = G4ThreeVector(5.2857373*mm,-29.4654967*mm,-36.5027*mm);
577     G4Transform3D walltransform6 = G4Transform3D(wallrotm,wallposition6);
578     fCagewall6_phys = new G4PVPlacement(walltransform6,fCagewall_log,"cagewall6",fExperimentalHall_log, true, 0);
579 
580     //-----------------------miniwalls
581     G4double minicagewall_x = 15.875*mm;
582     G4double minicagewall_y = 1*mm;
583 
584     auto *minicagewall_box = new G4Box("minicagewall_box",minicagewall_x*0.5,minicagewall_y*0.5,cagewall_z*0.5);
585     fMinicagewall_log = new G4LogicalVolume(minicagewall_box,Cu,"minicagewall_log",0,0,0);//tbc
586     fMinicagewall_phys = new G4PVPlacement(
587       0,
588       G4ThreeVector(33.8125*mm,0*mm,-36.5027*mm),
589       fMinicagewall_log,
590       "minicagewall1",
591       fExperimentalHall_log,
592       false,
593       0);
594 
595     G4RotationMatrix wallrotm2 = G4RotationMatrix();
596 
597     wallrotm2.rotateZ(-60*deg);
598     G4ThreeVector miniwallposition2 = G4ThreeVector(-16.90625*mm,29.282484*mm,-36.5027*mm);
599     G4Transform3D miniwalltransform2 = G4Transform3D(wallrotm2,miniwallposition2);
600     fMinicagewall2_phys = new G4PVPlacement(miniwalltransform2,fMinicagewall_log,"minicagewall2",fExperimentalHall_log, true, 0);
601 
602     wallrotm2.rotateZ(-60*deg);
603     G4ThreeVector miniwallposition3 = G4ThreeVector(-16.90625*mm,-29.282484*mm,-36.5027*mm);
604     G4Transform3D miniwalltransform3 = G4Transform3D(wallrotm2,miniwallposition3);
605     fMinicagewall3_phys = new G4PVPlacement(miniwalltransform3,fMinicagewall_log,"minicagewall3",fExperimentalHall_log, true, 0);
606 
607     //--------------extwalls
608     G4double extcagewall_x = 1.85*mm;
609     G4double extcagewall_y = 37.75*mm;
610     G4double extcagewall_z = 65*mm;
611 
612     auto *extcagewall_box = new G4Box("extcagewall_box",extcagewall_x*0.5,extcagewall_y*0.5,extcagewall_z*0.5);
613     fExtcagewall_log = new G4LogicalVolume(extcagewall_box,Cu,"extcagewall_log",0,0,0);
614     G4RotationMatrix wallrotm3 = G4RotationMatrix();
615     wallrotm3.rotateZ(-30*deg);
616     G4ThreeVector extwallposition = G4ThreeVector(-35.0724512*mm,20.2490891*mm,-36.5027*mm);
617     G4Transform3D extwalltransform = G4Transform3D(wallrotm3,extwallposition);
618     fExtcagewall_phys = new G4PVPlacement(extwalltransform,fExtcagewall_log,"extcagewall1",fExperimentalHall_log, true, 0);
619 
620     G4ThreeVector extwallposition2 = G4ThreeVector(35.0724512*mm,-20.2490891*mm,-36.5027*mm);
621     G4Transform3D extwalltransform2 = G4Transform3D(wallrotm3,extwallposition2);
622     fExtcagewall2_phys = new G4PVPlacement(extwalltransform2,fExtcagewall_log,"extcagewall2",fExperimentalHall_log, true, 0);
623 
624     //--------------extboards
625     G4double extboard_x = 0.5*mm;
626     G4double extboard_y = 36.5*mm;
627     G4double extboard_z = 64.5*mm;
628 
629     auto *extboard_box = new G4Box("extboard_box",extboard_x*0.5,extboard_y*0.5,extboard_z*0.5);
630     fExtboard_log = new G4LogicalVolume(extboard_box,Si,"extboard_log",0,0,0);
631     G4ThreeVector extboardposition = G4ThreeVector(37.932*mm,-21.9*mm,-36.253*mm);
632     G4Transform3D extboardtransform = G4Transform3D(wallrotm3,extboardposition);
633     fExtboard_phys = new G4PVPlacement(extboardtransform,fExtboard_log,"extboard1",fExperimentalHall_log, true, 0);
634 
635     G4ThreeVector extboardposition2 = G4ThreeVector(-37.932*mm,21.9*mm,-36.253*mm);
636     G4Transform3D extboardtransform2 = G4Transform3D(wallrotm3,extboardposition2);
637     fExtboard2_phys = new G4PVPlacement(extboardtransform2,fExtboard_log,"extboard2",fExperimentalHall_log, true, 0);
638 
639     //-------------------------------------------------------------------------------------------------
640     // --- END OF THE SUPPORTING STRUCTURES. HERE: STAGES OF THE CRYOSTAT -----------------------------
641 
642     //-------------------------------------------------------------------------------------------------
643     //----------------------------------------   First thermal shield  -----------------------------------------
644 
645     auto *FirstShield_botplate = new G4Tubs("Tube_FirstShield_botplate", 0.*mm, 80*mm, 1*mm*0.5, 0*deg, 360*deg);
646     G4LogicalVolume* FirstShield_botplate_log = new G4LogicalVolume(FirstShield_botplate,Cryoperm,"FirstShield_botplate_log",0,0,0);
647     fFirstShield_botplate_phys = new G4PVPlacement(
648       0,
649       G4ThreeVector(0,0,-101*mm),
650       FirstShield_botplate_log,
651       "FirstShield_botplate",
652       fExperimentalHall_log,
653       false,
654       0);
655 
656     //------------------------
657     auto *FirstShield_side = new G4Tubs("Tube_FirstShield_side", 79*mm, 80*mm, 145*mm*0.5, 0*deg, 360*deg);
658     G4LogicalVolume* FirstShield_side_log = new G4LogicalVolume(FirstShield_side,Cryoperm,"FirstShield_side_log",0,0,0);
659     fFirstShield_side_phys = new G4PVPlacement(
660       0,
661       G4ThreeVector(0,0,-28*mm),
662       FirstShield_side_log,
663       "FirstShield_side",
664       fExperimentalHall_log,
665       false,
666       0);
667 
668     //---------------------------
669     auto *FirstShield_topplate = new G4Tubs("Tube_FirstShield_topplate", 31.85*mm, 80*mm, 1*mm*0.5, 0*deg, 360*deg);
670     G4LogicalVolume* FirstShield_topplate_log = new G4LogicalVolume(FirstShield_topplate,Cryoperm,"FirstShield_topplate_log",0,0,0);
671     fFirstShield_topplate_phys = new G4PVPlacement(
672       0,
673       G4ThreeVector(0,0,45*mm),
674       FirstShield_topplate_log,
675       "FirstShield_topplate",
676       fExperimentalHall_log,
677       false,
678       0);
679 
680     //------------------------------
681     auto *FirstShield_topcyl = new G4Tubs("Tube_FirstShield_topcyl", 31.85*mm, 32.85*mm, 74*mm*0.5, 0*deg, 360*deg);
682     G4LogicalVolume* FirstShield_topcyl_log = new G4LogicalVolume(FirstShield_topcyl,Cryoperm,"FirstShield_topcyl_log",0,0,0);
683     fFirstShield_topcyl_phys = new G4PVPlacement(
684       0,
685       G4ThreeVector(0,0,82.5*mm),
686       FirstShield_topcyl_log,
687       "FirstShield_topcyl",
688       fExperimentalHall_log,
689       false,
690       0);
691 
692     //-------------------------------------------------------------------------------------------------
693     //------------------------------------------   Ring  ----------------------------------------------
694 
695     // This solid has a complex shape, filled with holes and supporting beams. Since it is not in direct sight of
696     // the detector, it can be simplified as a ring-like shape.
697 
698     // The material has been defined as a custom Al with density equal to the average density in the original solid
699     G4double hightOfRing = 85*mm;
700     auto *Ring = new G4Tubs("Ring_tubs", 80*mm, 97.1*mm, hightOfRing*0.5, 0*deg, 360*deg);
701 
702     // Three boxes are used for a boolean subtraction, cutting 3 sections of the ring
703     G4double cuts_x = 38.781347139*mm;
704     G4double cuts_y = 116.620103593*mm;
705     G4double cuts_z = 85.1*mm;
706 
707     auto *cuts_box = new G4Box("cuts_box",cuts_x*0.5,cuts_y*0.5,cuts_z*0.5);
708 
709     G4ThreeVector cutspos = G4ThreeVector(-105.055643631*mm,0,0);
710     G4SubtractionSolid* Ringsub = new G4SubtractionSolid("Ringsub",Ring,cuts_box,0,cutspos);
711 
712     auto *cutsrotm = new G4RotationMatrix();
713     cutsrotm->rotateZ(60*deg);
714     G4ThreeVector cutspos2 = G4ThreeVector(52.527821885*mm,-90.980856208*mm,0);
715     G4SubtractionSolid* Ringsub2 = new G4SubtractionSolid("Ringsub2",Ringsub,cuts_box,cutsrotm,cutspos2);
716     cutsrotm->rotateZ(60*deg);
717     G4ThreeVector cutspos3 = G4ThreeVector(52.527821885*mm,90.980856208*mm,0);
718     auto *Ringsub3 = new G4SubtractionSolid("Ringsub3",Ringsub2,cuts_box,cutsrotm,cutspos3);
719     auto *Ring_log = new G4LogicalVolume(Ringsub3,RingAl,"Ring_log",0,0,0);//tbc
720     fRing_phys = new G4PVPlacement(
721       0,
722       G4ThreeVector(0.,0.,-22*mm),
723       Ring_log,
724       "Ring",
725       fExperimentalHall_log,
726       false,
727       0);
728 
729     //-------------------------------------------------------------------------------------------------
730     //----------------------------------------   Second thermal shield  -----------------------------------------
731     // This is another irregular and complex shape derived from a CAD engineering model that was simplified for Geant4 implementation
732 
733     // Define the overall shape and bottom plate defining the base points to extrude
734     std::vector<G4TwoVector> newShape(6);
735     newShape[0] = G4TwoVector(58.936*mm, 154.115*mm);
736     newShape[1] = G4TwoVector(104*mm, 128.098*mm);
737     newShape[2] = G4TwoVector(104*mm, -128.098*mm);
738     newShape[3] = G4TwoVector(58.936*mm, -154.115*mm);
739     newShape[4] = G4TwoVector(-162.936*mm, -26.018*mm);
740     newShape[5] = G4TwoVector(-162.936*mm, 26.018*mm);
741 
742     // ------------------------------------------- the lower plate  -------------------------------------------
743     // Extrusion of the solid defined by the 6 points above
744     auto *MyShape = new G4ExtrudedSolid("MyShape", newShape, 2.65*0.5*mm, G4TwoVector(0, 0), 1.0, G4TwoVector(0, 0), 1.0);
745 
746     // Inserting a hole in the plate
747     auto *MyHole = new G4Tubs("MyHole", 0, 19*mm, 3*mm*0.5, 0*deg, 360*deg);
748     G4ThreeVector holepos = G4ThreeVector(0,0,0);
749     auto *SecondShieldbotPlate = new G4SubtractionSolid("SecondShieldbotPlate_sub",MyShape,MyHole,0,holepos);
750     auto *SecondShieldbotPlate_log = new G4LogicalVolume(SecondShieldbotPlate,Al,"SecondShieldbotPlate_log",0,0,0);
751     fSecondShieldbotPlate_phys = new G4PVPlacement(
752       0,
753       G4ThreeVector(0,0,-110.178*mm),
754       SecondShieldbotPlate_log,
755       "SecondShieldbotPlate",
756       fExperimentalHall_log,
757       false,
758       0);
759 
760     //---------------------------------------------------------------------------------------------------
761     //------------------------------------------ the central part ---------------------------------------
762 
763     // Definitino of the central block
764     auto *MyShape2 = new G4ExtrudedSolid("MyShape2", newShape, 215.35*0.5*mm, G4TwoVector(0, 0), 1.0, G4TwoVector(0, 0), 1.0);
765 
766     // The central part needs to be empty inside, so we define 6 new points for the subtraction of the central area
767     std::vector<G4TwoVector> newHole(6);
768     newHole[0] = G4TwoVector(58.936*mm, 149.762*mm);
769     newHole[1] = G4TwoVector(100.23*mm, 125.921*mm);
770     newHole[2] = G4TwoVector(100.23*mm, -125.921*mm);
771     newHole[3] = G4TwoVector(58.936*mm, -149.762*mm);
772     newHole[4] = G4TwoVector(-159.166*mm, -23.841*mm);
773     newHole[5] = G4TwoVector(-159.166*mm, 23.841*mm);
774 
775     // Extrusion of the subtraction solid
776     auto *MyInnerShape = new G4ExtrudedSolid("MyInnerShape", newHole, 216*0.5*mm, G4TwoVector(0, 0), 1.0, G4TwoVector(0, 0), 1.0);
777 
778     // creation of the central part: subtract the  subtraction solid from the central block
779     auto *SecondShieldside = new G4SubtractionSolid("SecondShieldside_sub",MyShape2,MyInnerShape,0,holepos);
780     auto *SecondShieldside_log = new G4LogicalVolume(SecondShieldside,Al,"SecondShieldside_log",0,0,0);
781     fSecondShieldside_phys = new G4PVPlacement(
782       0,
783       G4ThreeVector(0,0,-1.178*mm),
784       SecondShieldside_log,
785       "SecondShieldside",
786       fExperimentalHall_log,
787       false,
788       0);
789 
790     //-------------------------------------- top plate cover ---------------------------------------
791     G4ExtrudedSolid *MyShape3 = new G4ExtrudedSolid("MyShape3", newShape, 10.87*0.5*mm, G4TwoVector(0, 0), 1.0, G4TwoVector(0, 0), 1.0);
792     auto *MyHole2 = new G4Tubs("MyHole2", 0, 36.85*mm, 11*mm*0.5, 0*deg, 360*deg);
793     auto *SecondShieldtopPlate = new G4SubtractionSolid("SecondShieldtopPlate_sub",MyShape3,MyHole2,0,holepos);
794     auto *SecondShieldtopPlate_log = new G4LogicalVolume(SecondShieldtopPlate,Al,"SecondShieldtopPlate_log",0,0,0);
795     fSecondShieldtopPlate_phys = new G4PVPlacement(
796       0,
797       G4ThreeVector(0,0,111.932*mm),
798       SecondShieldtopPlate_log,
799       "SecondShieldtopPlate",
800       fExperimentalHall_log,
801       false,
802       0);
803 
804     //--------------------------------------
805     auto *SecondShieldshield_topcyl = new G4Tubs("SecondShieldshield", 36.85*mm, 39*mm, 12*mm*0.5, 0*deg, 360*deg);
806     auto *SecondShieldshield_topcyl_log = new G4LogicalVolume(SecondShieldshield_topcyl,Al,"SecondShieldshield_topcyl_log",0,0,0);
807     fSecondShieldshield_topcyl_phys = new G4PVPlacement(
808       0,
809       G4ThreeVector(0,0,123.367*mm),
810       SecondShieldshield_topcyl_log,
811       "SecondShieldshield_topcyl",
812       fExperimentalHall_log,
813       false,
814       0);
815 
816     //------------------------------------------------------------------------------------------------------
817     //-------------------------------------    EXTERNAL ZONE     -------------------------------------------
818 
819     //------------------------------------------------------------------------------------------------------
820     //---------------------------------------- Third thermal shield ---------------------------------------------------
821     G4double TopConeDOWNradius = 174.6*mm;
822     G4double TopConeUPradius = 93*mm;
823     G4double TopConeheight = 40*mm;
824     G4double SideHeight = 280*mm;
825     G4double ConeDOWNradius = 80*mm;
826     G4double Conethick = 0.4*mm;
827     G4double ConeUPradius = 174.6*mm;
828     G4double Coneheight = 40*mm;
829 
830     //---------------
831     auto *ThirdShieldside = new G4Tubs("ThirdShieldside_tubs", TopConeDOWNradius, TopConeDOWNradius+Conethick, SideHeight*0.5, 0*deg, 360*deg);
832     G4LogicalVolume* ThirdShieldside_log = new G4LogicalVolume(ThirdShieldside,Al,"ThirdShieldside_log",0,0,0);
833     fThirdShieldside_phys = new G4PVPlacement(
834       0,
835       G4ThreeVector(0,0,0),
836       ThirdShieldside_log,
837       "ThirdShieldside",
838       fExperimentalHall_log,
839       false,
840       0);
841 
842     //-------------------
843     auto *ThirdShieldtopCone = new G4Cons("ThirdShieldtopCone_cons", TopConeDOWNradius,Conethick+TopConeDOWNradius,TopConeUPradius,TopConeUPradius+Conethick, TopConeheight*0.5, 0*deg, 360*deg);
844     auto *ThirdShieldtopCone_log = new G4LogicalVolume(ThirdShieldtopCone,Al,"ThirdShieldtopCone_log",0,0,0);
845     fThirdShieldtopCone_phys = new G4PVPlacement(
846       0,
847       G4ThreeVector(0.,0.,(SideHeight+TopConeheight)*0.5),
848       ThirdShieldtopCone_log,
849       "ThirdShieldtopCone",
850       fExperimentalHall_log,
851       false,
852       0);
853 
854     //-------------------
855     auto *ThirdShieldbotCone = new G4Cons("ThirdShieldbotCone_cons", ConeDOWNradius,Conethick+ConeDOWNradius,ConeUPradius,ConeUPradius+Conethick, Coneheight*0.5, 0*deg, 360*deg);
856     auto *ThirdShieldbotCone_log = new G4LogicalVolume(ThirdShieldbotCone,Al,"ThirdShieldbotCone_log",0,0,0);
857     fThirdShieldbotCone_phys = new G4PVPlacement(
858       0,
859       G4ThreeVector(0.,0.,(-SideHeight-Coneheight)*0.5),
860       ThirdShieldbotCone_log,
861       "ThirdShieldbotCone",
862       fExperimentalHall_log,
863       false,
864       0);
865 
866     //-------------------
867     auto *ThirdShieldbotPlate = new G4Tubs("ThirdShieldbotPlate_tubs", 0*mm, Conethick+ConeDOWNradius,Conethick*0.5, 0*deg, 360*deg);
868     auto *ThirdShieldbotPlate_log = new G4LogicalVolume(ThirdShieldbotPlate,Al,"ThirdShieldbotPlate_log",0,0,0);
869     fThirdShieldbotPlate_phys = new G4PVPlacement(
870       0,
871       G4ThreeVector(0,0,(-SideHeight-Conethick)*0.5-Coneheight),
872       ThirdShieldbotPlate_log,
873       "ThirdShieldbotPlate",
874       fExperimentalHall_log,
875       false,
876       0);
877 
878     //----------------------------------------------------------------------------------------------------------------
879     //-------------------------------------------------- External layer of the Cryostat --------------------------------------------------
880     //----------------------------------------------------------------------------------------------------------------
881     G4double sideradius = 200*mm;
882     G4double midsidethick = 8.8*mm;
883     G4double midsideheigth = 120*mm;
884     auto *CryostatmidSide = new G4Tubs("CryostatmidSide_tubs", sideradius, sideradius+midsidethick, midsideheigth*0.5, 0*deg, 360*deg);
885     auto *CryostatmidSide_log = new G4LogicalVolume(CryostatmidSide,Al,"CryostatmidSide_log",0,0,0);
886     fCryostatmidSide_phys = new G4PVPlacement(
887       0,
888       G4ThreeVector(0,0,0),
889       CryostatmidSide_log,
890       "CryostatmidSide",
891       fExperimentalHall_log,
892       false,
893       0);
894 
895     //---------------
896     G4double topsidethick = 6.5*mm;
897     G4double topsideheigth = 85*mm;
898     auto *CryostattopSide = new G4Tubs("CryostattopSide_tubs", sideradius, sideradius+topsidethick, topsideheigth*0.5, 0*deg, 360*deg);
899     auto *CryostattopSide_log = new G4LogicalVolume(CryostattopSide,Al,"CryostattopSide_log",0,0,0);
900     fCryostattopSide_phys = new G4PVPlacement(
901       0,
902       G4ThreeVector(0,0,midsideheigth*0.5+topsideheigth*0.5),
903       CryostattopSide_log,
904       "CryostattopSide",
905       fExperimentalHall_log,
906       false,
907       0);
908 
909     //-------------------
910     ConeUPradius = 97*mm;
911     Coneheight = 75*mm;
912     auto *CryostattopCone = new G4Cons("CryostattopCone_cons", sideradius,topsidethick+sideradius,ConeUPradius,ConeUPradius+topsidethick, Coneheight*0.5, 0*deg, 360*deg);
913     auto *CryostattopCone_log = new G4LogicalVolume(CryostattopCone,Al,"CryostattopCone_log",0,0,0);
914     fCryostattopCone_phys = new G4PVPlacement(
915       0,
916       G4ThreeVector(0.,0.,midsideheigth*0.5+topsideheigth+Coneheight*0.5),  // translation position
917       CryostattopCone_log,
918       "CryostattopCone",
919       fExperimentalHall_log,
920       false,
921       0);
922 
923     //---------------
924     G4double botsideheigth = 85*mm;
925     auto *CryostatbotSide = new G4Tubs("CryostatbotSide_tubs", sideradius, sideradius+topsidethick, botsideheigth*0.5, 0*deg, 360*deg);
926     auto *CryostatbotSide_log = new G4LogicalVolume(CryostatbotSide,Al,"CryostatbotSide_log",0,0,0);
927     fCryostatbotSide_phys = new G4PVPlacement(
928       0,
929       G4ThreeVector(0,0,-midsideheigth*0.5-botsideheigth*0.5),
930       CryostatbotSide_log,
931       "CryostatbotSide",
932       fExperimentalHall_log,
933       false,
934       0);
935 
936     //-------------------
937     ConeDOWNradius = 120*mm;
938     auto *CryostatbotCone = new G4Cons("CryostatbotCone_cons", ConeDOWNradius,topsidethick+ConeDOWNradius,sideradius,sideradius+topsidethick, Coneheight*0.5, 0*deg, 360*deg);
939     auto *CryostatbotCone_log = new G4LogicalVolume(CryostatbotCone,Al,"CryostatbotCone_log",0,0,0);
940     fCryostatbotCone_phys = new G4PVPlacement(
941       0,
942       G4ThreeVector(0.,0.,-midsideheigth*0.5-botsideheigth-Coneheight*0.5), // translation position
943       CryostatbotCone_log,
944       "CryostatbotCone",
945       fExperimentalHall_log,
946       false,
947       0);
948 
949     //---------------------
950     auto *CryostatbotPlate = new G4Tubs("CryostatbotPlate_tubs", 0*mm, ConeDOWNradius, midsidethick*0.5, 0*deg, 360*deg);
951     auto *CryostatbotPlate_log = new G4LogicalVolume(CryostatbotPlate,Al,"CryostatbotPlate_log",0,0,0);
952     fCryostatbotPlate_phys = new G4PVPlacement(
953       0,
954       G4ThreeVector(0,0,-midsideheigth*0.5-botsideheigth-Coneheight+midsidethick*0.5),  // translation position
955       CryostatbotPlate_log,
956       "CryostatbotPlate",
957       fExperimentalHall_log,
958       false,
959       0);
960 
961     //------------------------------------------------------------------------------------------------------
962     //-------------------------------------  radiation filter   --------------------------------------------
963     //------------------------------------------------------------------------------------------------------
964     G4double halfhightOfAlFilter = 20*0.5*nm;
965     G4double mesh_thickness = 80*um;
966     G4double filtercarrierthickness = 5*mm;
967     G4double outerRadiusOfFilter = 128*0.5*mm;
968     G4double ZposFilter = 180*mm;
969 
970     auto *AlFilter = new G4Tubs("Filter", 0, outerRadiusOfFilter, halfhightOfAlFilter, 0.*deg, 360.*deg);
971     auto *AlFilter_log = new G4LogicalVolume(AlFilter,Al,"AlFilter_log",0,0,0);
972     fAlFilter_phys = new G4PVPlacement(nullptr,G4ThreeVector(0.,0.,ZposFilter), AlFilter_log, "AlFilter" ,fExperimentalHall_log,false,0);
973 
974     //------step limiter filter
975     AlFilter_log->SetUserLimits(fStepLimitAlfilta);
976 
977     //----- mesh: Nb
978     auto *Mesh = new G4Tubs("meshcage", 0, outerRadiusOfFilter+3.3*mm, mesh_thickness*0.5 , 0.*deg, 360.*deg);
979     fMesh_log = new G4LogicalVolume(Mesh,vacuum,"Mesh_log");   //vacuum, it's the mother volume and not the grid
980     G4double meshPos_z = ZposFilter-halfhightOfAlFilter-mesh_thickness*0.5;
981     fMesh_phys = new G4PVPlacement(nullptr, G4ThreeVector(DetPos_x,DetPos_y,meshPos_z),fMesh_log,"meshphys",fExperimentalHall_log, false, 0);
982 
983     //////// mesh pixels (grid support)
984     G4double half_x1 = 4827/2.*um;
985     G4double half_x2 = 4827/2.*um;
986     G4double half_x3 = 4880/2.*um;
987     G4double half_x4 = 4880/2.*um;
988     G4double half_y1 = 80/2.*um;
989     G4double half_y2 = 80/2.*um;
990     G4double dZ = 26.50/2.*um;
991 
992     G4double L1 = half_x1;
993     G4double L3 = dZ;
994 
995     // create the trapezoids for the element of the grid
996     auto *trapezoid_A = new G4Trap("Trapezoid_A",dZ,0,0,half_y1,half_x1,half_x2,0,half_y2,half_x3,half_x4,0);
997     auto *trapezoid_B = new G4Trap("Trapezoid_B",dZ,0,0,half_y1,half_x1,half_x2,0,half_y2,half_x3,half_x4,0);
998     auto *trapezoid_C = new G4Trap("Trapezoid_C",dZ,0,0,half_y1,half_x1,half_x2,0,half_y2,half_x3,half_x4,0);
999     auto *trapezoid_D = new G4Trap("Trapezoid_D",dZ,0,0,half_y1,half_x1,half_x2,0,half_y2,half_x3,half_x4,0);
1000 
1001     // create their rotation matrices and placement vectors
1002     auto *rotmat_A = new G4RotationMatrix();
1003     auto *rotmat_B = new G4RotationMatrix();
1004     auto *rotmat_C = new G4RotationMatrix();
1005     auto *rotmat_D = new G4RotationMatrix();
1006     rotmat_A->rotateX(90.*deg);
1007     rotmat_B->rotateX(90.*deg);
1008     rotmat_B->rotateY(-90.*deg);
1009     rotmat_C->rotateX(270.*deg);
1010     rotmat_D->rotateX(90.*deg);
1011     rotmat_D->rotateY(90.*deg);
1012     G4ThreeVector pos_trap_A = G4ThreeVector(0,L1+L3,meshPos_z);
1013     G4ThreeVector pos_trap_B = G4ThreeVector(L1+L3,0,meshPos_z);
1014     G4ThreeVector pos_trap_C = G4ThreeVector(0,-L1-L3,meshPos_z);
1015     G4ThreeVector pos_trap_D = G4ThreeVector(-L1-L3,0,meshPos_z);
1016 
1017     //Test: create their logical and physical volumes
1018     fTrapezoid_A_log = new G4LogicalVolume(trapezoid_A,SS,"trapezoid_A_log");
1019     fTrapezoid_B_log = new G4LogicalVolume(trapezoid_B,SS,"trapezoid_B_log");
1020     fTrapezoid_C_log = new G4LogicalVolume(trapezoid_C,SS,"trapezoid_C_log");
1021     fTrapezoid_D_log = new G4LogicalVolume(trapezoid_D,SS,"trapezoid_D_log");
1022 
1023     G4ThreeVector pos_munion_base = G4ThreeVector(0,0,0);
1024     auto *rotmat_final = new G4RotationMatrix();
1025     rotmat_final->rotateX(0.*deg);
1026     G4Transform3D transf_base = G4Transform3D(*rotmat_final, pos_munion_base);
1027 
1028     G4int index = 0;
1029     G4double coordx = 0;
1030     G4double coordy = 0;
1031     const G4int num = 13;
1032     G4double distance = 0;
1033 
1034     G4ThreeVector pos_base = G4ThreeVector(0,0,0);
1035     G4ThreeVector posA = G4ThreeVector(0,0,0); // the z coord is 0 because it will refer to Mesh_log
1036     G4ThreeVector posB = G4ThreeVector(0,0,0);
1037     G4ThreeVector posC = G4ThreeVector(0,0,0);
1038     G4ThreeVector posD = G4ThreeVector(0,0,0);
1039 
1040     std::vector<G4VPhysicalVolume*> phys_array;
1041     for (int i=-num; i<=num; i++)
1042     {
1043       for (int j=-num; j<=num; j++)
1044       {
1045         std::string sA = "trapezoidA_";
1046         std::string sB = "trapezoidB_";
1047         std::string sC = "trapezoidC_";
1048         std::string sD = "trapezoidD_";
1049 
1050         std::string index_str = std::to_string(index);
1051         sA += index_str;
1052         sB += index_str;
1053         sC += index_str;
1054         sD += index_str;
1055 
1056         coordx = i*half_x3*2;
1057         coordy = j*half_x3*2;
1058 
1059         distance = std::sqrt(coordx*coordx + coordy*coordy);
1060 
1061         if (distance < outerRadiusOfFilter)
1062         {
1063           pos_base = G4ThreeVector(coordx,coordy,meshPos_z);
1064           posA = G4ThreeVector(coordx,coordy+L1+L3,0); // the z coord is 0 because it will refer to fMesh_log
1065           posB = G4ThreeVector(coordx+L1+L3,coordy,0);
1066           posC = G4ThreeVector(coordx,coordy-L1-L3,0);
1067           posD = G4ThreeVector(coordx-L1-L3,coordy,0);
1068 
1069           fTrapezoid_A_phys = new G4PVPlacement(rotmat_A,posA,fTrapezoid_A_log,sA,fMesh_log, false, 0);
1070           fTrapezoid_B_phys = new G4PVPlacement(rotmat_B,posB,fTrapezoid_B_log,sB,fMesh_log, false, 0);
1071           fTrapezoid_C_phys = new G4PVPlacement(rotmat_C,posC,fTrapezoid_C_log,sC,fMesh_log, false, 0);
1072           fTrapezoid_D_phys = new G4PVPlacement(rotmat_D,posD,fTrapezoid_D_log,sD,fMesh_log, false, 0);
1073 
1074           phys_array.push_back(fTrapezoid_A_phys);
1075           phys_array.push_back(fTrapezoid_B_phys);
1076           phys_array.push_back(fTrapezoid_C_phys);
1077           phys_array.push_back(fTrapezoid_D_phys);
1078           index += 1;
1079         }
1080       }
1081     }
1082 
1083     auto *filtercarrier = new G4Tubs("filtercarrier_tubs", outerRadiusOfFilter,outerRadiusOfFilter+3.85*mm ,filtercarrierthickness*0.5,0.*deg,360.*deg);
1084     auto *filtercarrier_log = new G4LogicalVolume(filtercarrier,Al,"filtercarrier_log",0,0,0);
1085     fFiltercarrier_phys = new G4PVPlacement(nullptr,G4ThreeVector(0.,0.,meshPos_z+mesh_thickness*0.5+filtercarrierthickness*0.5),filtercarrier_log,"filtercarrier",fExperimentalHall_log,false,0);
1086 
1087     //------------------------------------------------------------------------------------------------------
1088     //----------------------------------------   aperture cylinder   ---------------------------------------
1089     //------------------------------------------------------------------------------------------------------
1090     //300K piece
1091     //G4double ACinnerDOWNradius = 106*mm;
1092     G4double ACinnerUPradius = 67.8644*mm;
1093     G4double ACthick = 1*mm;
1094     G4double ACheight = midsideheigth*0.5+topsideheigth+Coneheight-ZposFilter-filtercarrierthickness;
1095 
1096     //ZposFilter
1097     auto *AC = new G4Cons("AC_cons",ACinnerUPradius-ACthick,ACinnerUPradius,ConeUPradius-ACthick,ConeUPradius,ACheight*0.5,0*deg,360*deg);
1098     auto *AC_log = new G4LogicalVolume(AC,Ti,"AC_log",0,0,0);
1099     fAC_phys = new G4PVPlacement(
1100       0,
1101       G4ThreeVector(0.,0.,midsideheigth*0.5+topsideheigth+Coneheight-ACheight*0.5),
1102       AC_log,
1103       "AC",
1104       fExperimentalHall_log,
1105       false,
1106       0);
1107 
1108     //-----inner coating
1109     auto *ACIC = new G4Cons("ACIC_cons",ACinnerUPradius-ACthick-5*um,ACinnerUPradius-ACthick,ConeUPradius-ACthick-5*um,ConeUPradius-ACthick,ACheight*0.5,0*deg,360*deg);
1110     auto *ACIC_log = new G4LogicalVolume(ACIC,Au,"ACIC_log",0,0,0);
1111     fACIC_phys = new G4PVPlacement(
1112       0,
1113       G4ThreeVector(0.,0.,midsideheigth*0.5+topsideheigth+Coneheight-ACheight*0.5),
1114       ACIC_log,
1115       "ACIC",
1116       fExperimentalHall_log,
1117       false,
1118       0);
1119 
1120     //------------------------------------------------------------------------------------------------------
1121     //------------------------------    The Al sphere representing the spacecraft     ----------------------
1122     //------------------------------------------------------------------------------------------------------
1123     auto *sphere = new G4Sphere("sphere_sphere",26*cm,27*cm,0.,twopi,0.,pi);
1124     fSphere_log = new G4LogicalVolume(sphere,Al,"sphere_log");
1125     fSphere_phys = new G4PVPlacement(nullptr,G4ThreeVector(0.,0.,0*mm),fSphere_log,"sphere",fExperimentalHall_log,false,0);
1126 
1127     //------------------------------------------------------------------------------------------------------
1128     //---------------------    regions assignment: closer to the detector, lower the Ecut     --------------
1129     //------------------------------------------------------------------------------------------------------
1130     // everything directly seen by the detector
1131     InnerRegion->AddRootLogicalVolume(fElement_log);
1132     InnerRegion->AddRootLogicalVolume(fDetector_log);
1133     InnerRegion->AddRootLogicalVolume(fBipxl_log);
1134     InnerRegion->AddRootLogicalVolume(fMembranepxl_log);
1135     InnerRegion->AddRootLogicalVolume(fGridpiece_log);
1136     InnerRegion->AddRootLogicalVolume(fWafer_log);
1137     InnerRegion->AddRootLogicalVolume(fACDpxl_log);
1138     InnerRegion->AddRootLogicalVolume(AlFilter_log);  // it was in the intermediate region
1139     InnerRegion->AddRootLogicalVolume(ACIC_log);
1140 
1141     fWorld_phys = fExperimentalHall_phys;
1142   }
1143   else
1144   {
1145     G4cout << "Build geometry from GDML file: " << fReadFile << G4endl;
1146     G4VPhysicalVolume* fWorldPhysVol;
1147     fParser.Read(fReadFile);
1148     fWorldPhysVol = fParser.GetWorldVolume();
1149     fWorld_log =  fWorldPhysVol->GetLogicalVolume();
1150     G4int num_volumes = 9;
1151     auto *inner_volumes = new G4String[num_volumes];
1152     inner_volumes[0] = "element_log"; inner_volumes[1] = "Detector_log"; inner_volumes[2] = "Bipxl_log"; inner_volumes[3] = "membranepxl_log";
1153     inner_volumes[4] = "gridpiece_log"; inner_volumes[5] = "wafer_log"; inner_volumes[6] = "ACDpxl_log"; inner_volumes[7] = "AlFilter_log";
1154     inner_volumes[8] = "ACIC_log";
1155 
1156     G4cout << "Adding volumes to the InnerRegion" << G4endl;
1157     // Add volumes to the InnerRegion
1158     for (G4int j = 0; j <= num_volumes; j++)
1159     {
1160       for (G4int i = 0; i < (G4int)fWorld_log->GetNoDaughters(); i++)
1161       {
1162         if (fWorld_log->GetDaughter(i)->GetLogicalVolume()->GetName() == inner_volumes[j])
1163         {
1164           InnerRegion->AddRootLogicalVolume(fWorld_log->GetDaughter(i)->GetLogicalVolume());
1165           G4cout << "Added: " << fWorld_log->GetDaughter(i)->GetLogicalVolume()->GetName() << G4endl;
1166         }
1167       }
1168     }
1169 
1170     // Read volume names
1171     auto *pvs = G4PhysicalVolumeStore::GetInstance();
1172     if (pvs == nullptr)
1173     {
1174       G4cout << "PhysicalVolumeStore not accessible" << G4endl;
1175     }
1176     else
1177     {
1178       G4cout << "Volumes imported in the PhysicalVolumeStore. " << G4endl;
1179       G4int length = 0;
1180       length = pvs->size();
1181       for (G4int i= 0; i<length; i++)
1182       {
1183         G4String volName = ((*pvs)[i])->GetName();
1184       }
1185     }
1186     fWorld_phys = fWorldPhysVol;
1187   }
1188 
1189   // Set production cuts
1190   G4ProductionCutsTable::GetProductionCutsTable()->SetEnergyRange(200*eV, 1000*GeV);
1191   auto* InnerCuts = new G4ProductionCuts;
1192   InnerCuts->SetProductionCut(1*mm);
1193   InnerRegion->SetProductionCuts(InnerCuts);
1194 
1195   return fWorld_phys;
1196 }
1197 
1198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1199 
1200 void XrayTESdetDetectorConstruction::SetMaxStep(G4double maxStep)
1201 {
1202   if ((fStepLimitAlfilta != nullptr)&&(maxStep>0.)){ fStepLimitAlfilta->SetMaxAllowedStep(maxStep); }
1203   if ((fStepLimitPolyfilta  != nullptr)&&(maxStep>0.)){ fStepLimitPolyfilta->SetMaxAllowedStep(maxStep); }
1204   if ((fStepLimitmembrane  != nullptr)&&(maxStep>0.)){ fStepLimitmembrane->SetMaxAllowedStep(maxStep); }
1205 }
1206 
1207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1208 
1209 void XrayTESdetDetectorConstruction::SetReadFile(G4String input_geometry)
1210 {
1211   fReadFile = std::move(input_geometry);
1212 }
1213 
1214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1215