Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LogicalCrystalVolume.cc
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// Implementation of G4LogicalCrystalVolume
27//
28// 21-04-16, created by E.Bagli
29// --------------------------------------------------------------------
30
32#include "G4ExtendedMaterial.hh"
33#include "G4CrystalExtension.hh"
35
36std::vector<G4LogicalVolume*> G4LogicalCrystalVolume::fLCVvec;
37
38// --------------------------------------------------------------------
39
42 const G4String& name, G4FieldManager* pFieldMgr,
43 G4VSensitiveDetector* pSDetector,
44 G4UserLimits* pULimits, G4bool optimise,
45 G4int h, G4int k, G4int l, G4double rot)
46: G4LogicalVolume(pSolid,pMaterial,name,pFieldMgr,pSDetector,pULimits,optimise)
47{
48 SetMillerOrientation(h, k, l, rot);
49 fLCVvec.push_back(this);
50}
51
52// --------------------------------------------------------------------
53
55{
56 fLCVvec.erase( std::remove(fLCVvec.begin(),fLCVvec.end(), this ),
57 fLCVvec.end() );
58}
59
60// --------------------------------------------------------------------
61
63{
64 return std::find(fLCVvec.begin(), fLCVvec.end(), aLV) != fLCVvec.end();
65}
66
67// --------------------------------------------------------------------
68
70{
71 return dynamic_cast<G4CrystalExtension*>(
72 dynamic_cast<G4ExtendedMaterial*>(GetMaterial() )
73 ->RetrieveExtension("crystal"));
74}
75
76// --------------------------------------------------------------------
77
79{
80 return GetCrystal()->GetUnitCell()->GetBasis(i);
81}
82
83// --------------------------------------------------------------------
84
86 G4int k,
87 G4int l,
88 G4double rot)
89{
90 // Align Miller normal vector (hkl) with +Z axis, and rotation about axis
91
92 if (verboseLevel)
93 {
94 G4cout << "G4LatticePhysical::SetMillerOrientation(" << h << " "
95 << k << " " << l << ", " << rot/CLHEP::deg << " deg)" << G4endl;
96 }
97
98 hMiller = h;
99 kMiller = k;
100 lMiller = l;
101 fRot = rot;
102
103 G4ThreeVector norm = (h*GetBasis(0)+k*GetBasis(1)+l*GetBasis(2)).unit();
104
105 if (verboseLevel>1) G4cout << " norm = " << norm << G4endl;
106
107 // Aligns geometry +Z axis with lattice (hkl) normal
109 fOrient.rotateZ(rot).rotateY(norm.theta()).rotateZ(norm.phi());
110 fInverse = fOrient.inverse();
111
112 if (verboseLevel>1) G4cout << " fOrient = " << fOrient << G4endl;
113
114 // FIXME: Is this equivalent to (phi,theta,rot) Euler angles???
115}
116
117// --------------------------------------------------------------------
118
119// Rotate input vector between lattice and solid orientations
120
121const G4ThreeVector&
123{
124 return dir.transform(fOrient);
125}
126
127const G4ThreeVector&
129{
130 return dir.transform(fInverse);
131}
132
133// --------------------------------------------------------------------
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double phi() const
double theta() const
Hep3Vector & transform(const HepRotation &)
Definition: ThreeVectorR.cc:20
HepRotation inverse() const
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:366
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87
HepRotation & rotateY(double delta)
Definition: Rotation.cc:74
G4CrystalUnitCell * GetUnitCell() const
const G4ThreeVector & GetBasis(G4int idx) const
G4LogicalCrystalVolume(G4VSolid *pSolid, G4ExtendedMaterial *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=nullptr, G4VSensitiveDetector *pSDetector=nullptr, G4UserLimits *pULimits=nullptr, G4bool optimise=true, G4int h=0, G4int k=0, G4int l=0, G4double rot=0.0)
const G4ThreeVector & RotateToSolid(G4ThreeVector &dir) const
static G4bool IsLattice(G4LogicalVolume *aLV)
const G4ThreeVector & GetBasis(G4int i) const
const G4ThreeVector & RotateToLattice(G4ThreeVector &dir) const
const G4CrystalExtension * GetCrystal() const
void SetMillerOrientation(G4int h, G4int k, G4int l, G4double rot=0.0)
G4Material * GetMaterial() const