Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/exoticphysics/channeling/ch1/src/DetectorConstruction.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 /// \file DetectorConstruction.cc
 27 /// \brief Implementation of the DetectorConstruction class
 28 //
 29 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 30 
 31 #include "DetectorConstruction.hh"
 32 
 33 #include "G4RunManager.hh"
 34 #include "G4NistManager.hh"
 35 #include "G4Box.hh"
 36 #include "G4LogicalVolume.hh"
 37 #include "G4RegionStore.hh"
 38 #include "G4VisAttributes.hh"
 39 #include "PrimaryGeneratorAction.hh"
 40 #include <CLHEP/Units/SystemOfUnits.h>
 41 
 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 43 
 44 DetectorConstruction::DetectorConstruction()
 45 {}
 46 
 47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 48 
 49 G4VPhysicalVolume* DetectorConstruction::Construct()
 50 {
 51     //Check overlap option
 52     G4bool checkOverlaps = true;
 53 
 54     //Materials
 55     G4NistManager* nist = G4NistManager::Instance();
 56     G4Material* world_mat = nist->FindOrBuildMaterial("G4_Galactic");
 57     G4Material* silicon = nist->FindOrBuildMaterial("G4_Si");
 58 
 59     //World
 60     G4Box* solidWorld = new G4Box("World", 0.2*CLHEP::m, 0.2*CLHEP::m, 10.*CLHEP::m);
 61     G4LogicalVolume* logicWorld = new G4LogicalVolume(solidWorld, world_mat, "World");  
 62     G4VPhysicalVolume* physWorld = new G4PVPlacement
 63                                     (0,                    // no rotation
 64                                     G4ThreeVector(),       // centre position
 65                                     logicWorld,            // its logical volume
 66                                     "World",               // its name
 67                                     0,                     // its mother volume
 68                                     false,                 // no boolean operation
 69                                     0,                     // copy number
 70                                     checkOverlaps);        // overlaps checking
 71     logicWorld->SetVisAttributes(G4VisAttributes::GetInvisible());
 72 
 73 
 74     // --------------- Crystal ------------------------------------
 75 
 76     /*parameters of the following experiment:
 77       Only channeling: A. Mazzolari et al. Phys. Rev. Lett. 112, 135503 (2014)
 78       Radition:        L. Bandiera et al. Phys. Rev. Lett. 115, 025504 (2015)
 79       Published experimental validation of G4ChannelingFastSimModel (only channeling):
 80       A. Sytov et al. Journal of the Korean Physical Society 83, 132–139 (2023)
 81     */
 82 
 83     //Select crystal material
 84     fCrystalMaterial = nist->FindOrBuildMaterial("G4_Si");
 85 
 86     //Setting crystal rotation angle (also the angle of crystal planes vs the beam)
 87     //Crystal rotation angle (also the angle of crystal planes vs the beam)
 88     G4double angleX = 0.*1e-6; //rad
 89     G4RotationMatrix* crystalRotationMatrix = new G4RotationMatrix;
 90     crystalRotationMatrix->rotateY(-angleX);
 91 
 92     //Crystal bending angle
 93     fBendingAngle = 0.905*CLHEP::mrad;
 94 
 95     //setting crystal dimensions:
 96     G4ThreeVector crystalSize = G4ThreeVector(20.*CLHEP::mm,
 97                                               20.*CLHEP::mm,
 98                                               0.0305*CLHEP::mm);
 99 
100     //Setting crystal position
101     G4ThreeVector posCrystal = G4ThreeVector(0., 0., crystalSize.z()/2.);
102 
103     //crystal volume
104     G4Box* solidCrystal = new G4Box("Crystal",
105                                     crystalSize.x()/2,
106                                     crystalSize.y()/2,
107                                     crystalSize.z()/2.);
108     
109     fLogicCrystal = new G4LogicalVolume(solidCrystal,
110                                        fCrystalMaterial,
111                                        "Crystal");
112     new G4PVPlacement(crystalRotationMatrix,
113                       posCrystal,
114                       fLogicCrystal,
115                       "Crystal",
116                       logicWorld,
117                       false,
118                       0,
119                       checkOverlaps);
120     
121     //crystal region (necessary for the FastSim model)
122     G4Region* regionCh = new G4Region("Crystal");
123     regionCh->AddRootLogicalVolume(fLogicCrystal);
124 
125     //visualization attributes
126     G4VisAttributes* crystalVisAttribute =
127         new G4VisAttributes(G4Colour(1., 0., 0.));
128     crystalVisAttribute->SetForceSolid(true);
129     fLogicCrystal->SetVisAttributes(crystalVisAttribute);
130         
131     //print crystal info
132     G4cout << "Crystal size: " << crystalSize.x()/CLHEP::mm
133            << " " << crystalSize.y()/CLHEP::mm
134            << " " << crystalSize.z()/CLHEP::mm << " mm3" << G4endl;
135     G4cout << "Crystal bending angle: " << fBendingAngle << " rad" << G4endl;
136     G4cout << "Crystal angleX: " << angleX << " rad" << G4endl;
137 
138     // --------------- Detector -----------------------------------
139     //Setting detector position
140     G4ThreeVector posDetector = G4ThreeVector(0, 0, 5973*CLHEP::mm);
141 
142     //particle detector volume
143     G4Box* detector = new G4Box("Detector",
144                                 10*CLHEP::cm/2,
145                                 10*CLHEP::cm/2,
146                                 0.3*CLHEP::mm/2);
147     
148     G4LogicalVolume* logicDetector = new G4LogicalVolume(detector,
149                                                          silicon,
150                                                          "Detector");
151     new G4PVPlacement(0, 
152                       posDetector, 
153                       logicDetector,
154                       "Detector", 
155                       logicWorld, 
156                       false, 
157                       0, 
158                       checkOverlaps);
159 
160     //visualization attributes
161     G4VisAttributes* detectorVisAttribute =
162         new G4VisAttributes(G4Colour(0., 0., 1));
163     detectorVisAttribute->SetForceSolid(true);
164     logicDetector->SetVisAttributes(detectorVisAttribute);
165     
166     //always return the physical World
167     return physWorld;
168 }
169 
170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171 
172 void DetectorConstruction::ConstructSDandField()
173 {
174     // --------------- fast simulation ----------------------------
175     //extract the region of the crystal from the store
176     G4RegionStore* regionStore = G4RegionStore::GetInstance();
177     G4Region* regionCh = regionStore->GetRegion("Crystal");
178 
179     ///create the channeling model for this region
180     G4ChannelingFastSimModel* channelingModel =
181         new G4ChannelingFastSimModel("ChannelingModel", regionCh);
182 
183     ///Crystal planes or axes considered
184     ///Use brackets (...) for planes and <...> for axes
185     G4String lattice = "(111)";
186 
187     ///activate the channeling model
188     channelingModel->Input(fCrystalMaterial, lattice);
189     ///setting bending angle of the crystal planes (default is 0)
190     channelingModel->GetCrystalData()->SetBendingAngle(fBendingAngle,fLogicCrystal);
191 
192     /*
193     activate radiation model (do it only when you want to take into account the
194     radiation production in an oriented crystal; it reduces simulation speed.)
195     */
196     G4bool activateRadiationModel = true;
197     if (activateRadiationModel)
198     {
199         channelingModel->RadiationModelActivate();
200         G4cout << "Radiation model activated" << G4endl;
201     }
202 }
203 
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205 
206