Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LogicalCrystalVolume Class Reference

#include <G4LogicalCrystalVolume.hh>

Inheritance diagram for G4LogicalCrystalVolume:
Collaboration diagram for G4LogicalCrystalVolume:

Public Member Functions

 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.)
 
 ~G4LogicalCrystalVolume ()
 
G4bool IsExtended () const
 
void SetMillerOrientation (G4int h, G4int k, G4int l, G4double rot=0.)
 
const G4ThreeVectorRotateToLattice (G4ThreeVector &dir) const
 
const G4ThreeVectorRotateToSolid (G4ThreeVector &dir) const
 
const G4CrystalExtensionGetCrystal () const
 
const G4ThreeVectorGetBasis (G4int i) const
 
void SetVerbose (G4int aInt)
 
- Public Member Functions inherited from G4LogicalVolume
 G4LogicalVolume (G4VSolid *pSolid, G4Material *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=0, G4VSensitiveDetector *pSDetector=0, G4UserLimits *pULimits=0, G4bool optimise=true)
 
virtual ~G4LogicalVolume ()
 
const G4StringGetName () const
 
void SetName (const G4String &pName)
 
G4int GetNoDaughters () const
 
G4VPhysicalVolumeGetDaughter (const G4int i) const
 
void AddDaughter (G4VPhysicalVolume *p)
 
G4bool IsDaughter (const G4VPhysicalVolume *p) const
 
G4bool IsAncestor (const G4VPhysicalVolume *p) const
 
void RemoveDaughter (const G4VPhysicalVolume *p)
 
void ClearDaughters ()
 
G4int TotalVolumeEntities () const
 
EVolume CharacteriseDaughters () const
 
G4VSolidGetSolid () const
 
void SetSolid (G4VSolid *pSolid)
 
G4MaterialGetMaterial () const
 
void SetMaterial (G4Material *pMaterial)
 
void UpdateMaterial (G4Material *pMaterial)
 
G4double GetMass (G4bool forced=false, G4bool propagate=true, G4Material *parMaterial=0)
 
void ResetMass ()
 
G4FieldManagerGetFieldManager () const
 
void SetFieldManager (G4FieldManager *pFieldMgr, G4bool forceToAllDaughters)
 
G4VSensitiveDetectorGetSensitiveDetector () const
 
void SetSensitiveDetector (G4VSensitiveDetector *pSDetector)
 
G4UserLimitsGetUserLimits () const
 
void SetUserLimits (G4UserLimits *pULimits)
 
G4SmartVoxelHeaderGetVoxelHeader () const
 
void SetVoxelHeader (G4SmartVoxelHeader *pVoxel)
 
G4double GetSmartless () const
 
void SetSmartless (G4double s)
 
G4bool IsToOptimise () const
 
void SetOptimisation (G4bool optim)
 
G4bool IsRootRegion () const
 
void SetRegionRootFlag (G4bool rreg)
 
G4bool IsRegion () const
 
void SetRegion (G4Region *reg)
 
G4RegionGetRegion () const
 
void PropagateRegion ()
 
const G4MaterialCutsCoupleGetMaterialCutsCouple () const
 
void SetMaterialCutsCouple (G4MaterialCutsCouple *cuts)
 
G4bool operator== (const G4LogicalVolume &lv) const
 
const G4VisAttributesGetVisAttributes () const
 
void SetVisAttributes (const G4VisAttributes *pVA)
 
void SetVisAttributes (const G4VisAttributes &VA)
 
G4FastSimulationManagerGetFastSimulationManager () const
 
void SetBiasWeight (G4double w)
 
G4double GetBiasWeight () const
 
 G4LogicalVolume (__void__ &)
 
G4FieldManagerGetMasterFieldManager () const
 
G4VSensitiveDetectorGetMasterSensitiveDetector () const
 
G4VSolidGetMasterSolid () const
 
G4int GetInstanceID () const
 
void Lock ()
 
void InitialiseWorker (G4LogicalVolume *ptrMasterObject, G4VSolid *pSolid, G4VSensitiveDetector *pSDetector)
 
void TerminateWorker (G4LogicalVolume *ptrMasterObject)
 
void AssignFieldManager (G4FieldManager *fldMgr)
 

Static Public Member Functions

static G4bool IsLattice (G4LogicalVolume *aLV)
 
- Static Public Member Functions inherited from G4LogicalVolume
static const G4LVManagerGetSubInstanceManager ()
 
static G4VSolidGetSolid (G4LVData &instLVdata)
 
static void SetSolid (G4LVData &instLVdata, G4VSolid *pSolid)
 

Detailed Description

Definition at line 55 of file G4LogicalCrystalVolume.hh.

Constructor & Destructor Documentation

G4LogicalCrystalVolume::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. 
)

Definition at line 46 of file G4LogicalCrystalVolume.cc.

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 }
G4LogicalVolume(G4VSolid *pSolid, G4Material *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=0, G4VSensitiveDetector *pSDetector=0, G4UserLimits *pULimits=0, G4bool optimise=true)
void SetMillerOrientation(G4int h, G4int k, G4int l, G4double rot=0.)

Here is the call graph for this function:

G4LogicalCrystalVolume::~G4LogicalCrystalVolume ( )

Definition at line 60 of file G4LogicalCrystalVolume.cc.

61 {
62  fLCVvec.erase( std::remove(fLCVvec.begin(),fLCVvec.end(), this ),
63  fLCVvec.end() );
64 }

Member Function Documentation

const G4ThreeVector & G4LogicalCrystalVolume::GetBasis ( G4int  i) const

Definition at line 83 of file G4LogicalCrystalVolume.cc.

84 {
85  return GetCrystal()->GetUnitCell()->GetBasis(i);
86 }
G4CrystalUnitCell * GetUnitCell() const
const G4ThreeVector & GetBasis(G4int idx) const
const G4CrystalExtension * GetCrystal() const

Here is the call graph for this function:

Here is the caller graph for this function:

const G4CrystalExtension * G4LogicalCrystalVolume::GetCrystal ( ) const

Definition at line 75 of file G4LogicalCrystalVolume.cc.

76 {
77  return dynamic_cast<G4CrystalExtension*>(dynamic_cast<G4ExtendedMaterial*>(GetMaterial())
78  ->RetrieveExtension("crystal"));
79 }
G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4LogicalCrystalVolume::IsExtended ( ) const
inlinevirtual

Reimplemented from G4LogicalVolume.

Definition at line 73 of file G4LogicalCrystalVolume.hh.

73 { return true; }
G4bool G4LogicalCrystalVolume::IsLattice ( G4LogicalVolume aLV)
static

Definition at line 68 of file G4LogicalCrystalVolume.cc.

69 {
70  return std::find(fLCVvec.begin(), fLCVvec.end(), aLV) != fLCVvec.end();
71 }
const G4ThreeVector & G4LogicalCrystalVolume::RotateToLattice ( G4ThreeVector dir) const

Definition at line 127 of file G4LogicalCrystalVolume.cc.

128 {
129  return dir.transform(fOrient);
130 }
Hep3Vector & transform(const HepRotation &)
Definition: ThreeVectorR.cc:24

Here is the call graph for this function:

const G4ThreeVector & G4LogicalCrystalVolume::RotateToSolid ( G4ThreeVector dir) const

Definition at line 133 of file G4LogicalCrystalVolume.cc.

134 {
135  return dir.transform(fInverse);
136 }
Hep3Vector & transform(const HepRotation &)
Definition: ThreeVectorR.cc:24

Here is the call graph for this function:

void G4LogicalCrystalVolume::SetMillerOrientation ( G4int  h,
G4int  k,
G4int  l,
G4double  rot = 0. 
)

Definition at line 90 of file G4LogicalCrystalVolume.cc.

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 }
HepRotation & rotateY(double delta)
Definition: Rotation.cc:79
HepRotation inverse() const
static DLL_API const HepRotation IDENTITY
Definition: Rotation.h:369
G4GLOB_DLL std::ostream G4cout
double phi() const
double theta() const
static constexpr double deg
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:92
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetBasis(G4int i) const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4LogicalCrystalVolume::SetVerbose ( G4int  aInt)
inline

Definition at line 90 of file G4LogicalCrystalVolume.hh.

90 { verboseLevel = aInt; }

The documentation for this class was generated from the following files: