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