Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LogicalCrystalVolume Class Reference

#include <G4LogicalCrystalVolume.hh>

+ Inheritance diagram for G4LogicalCrystalVolume:

Public Member Functions

 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)
 
 ~G4LogicalCrystalVolume () override
 
G4bool IsExtended () const override
 
void SetMillerOrientation (G4int h, G4int k, G4int l, G4double rot=0.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=nullptr, G4VSensitiveDetector *pSDetector=nullptr, G4UserLimits *pULimits=nullptr, G4bool optimise=true)
 
virtual ~G4LogicalVolume ()
 
 G4LogicalVolume (const G4LogicalVolume &)=delete
 
G4LogicalVolumeoperator= (const G4LogicalVolume &)=delete
 
const G4StringGetName () const
 
void SetName (const G4String &pName)
 
std::size_t GetNoDaughters () const
 
G4VPhysicalVolumeGetDaughter (const std::size_t 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
 
EVolume DeduceDaughtersType () 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=nullptr)
 
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)
 
G4bool ChangeDaughtersType (EVolume atype)
 

Static Public Member Functions

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

Detailed Description

Definition at line 46 of file G4LogicalCrystalVolume.hh.

Constructor & Destructor Documentation

◆ G4LogicalCrystalVolume()

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

Definition at line 40 of file G4LogicalCrystalVolume.cc.

46: G4LogicalVolume(pSolid,pMaterial,name,pFieldMgr,pSDetector,pULimits,optimise)
47{
48 SetMillerOrientation(h, k, l, rot);
49 fLCVvec.push_back(this);
50}
void SetMillerOrientation(G4int h, G4int k, G4int l, G4double rot=0.0)
G4LogicalVolume(G4VSolid *pSolid, G4Material *pMaterial, const G4String &name, G4FieldManager *pFieldMgr=nullptr, G4VSensitiveDetector *pSDetector=nullptr, G4UserLimits *pULimits=nullptr, G4bool optimise=true)

◆ ~G4LogicalCrystalVolume()

G4LogicalCrystalVolume::~G4LogicalCrystalVolume ( )
override

Definition at line 54 of file G4LogicalCrystalVolume.cc.

55{
56 fLCVvec.erase( std::remove(fLCVvec.begin(),fLCVvec.end(), this ),
57 fLCVvec.end() );
58}

Member Function Documentation

◆ GetBasis()

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

Definition at line 78 of file G4LogicalCrystalVolume.cc.

79{
80 return GetCrystal()->GetUnitCell()->GetBasis(i);
81}
G4CrystalUnitCell * GetUnitCell() const
const G4ThreeVector & GetBasis(G4int idx) const
const G4CrystalExtension * GetCrystal() const

Referenced by SetMillerOrientation().

◆ GetCrystal()

const G4CrystalExtension * G4LogicalCrystalVolume::GetCrystal ( ) const

Definition at line 69 of file G4LogicalCrystalVolume.cc.

70{
71 return dynamic_cast<G4CrystalExtension*>(
72 dynamic_cast<G4ExtendedMaterial*>(GetMaterial() )
73 ->RetrieveExtension("crystal"));
74}
G4VMaterialExtension * RetrieveExtension(const G4String &name)
G4Material * GetMaterial() const

Referenced by GetBasis().

◆ IsExtended()

G4bool G4LogicalCrystalVolume::IsExtended ( ) const
inlineoverridevirtual

Reimplemented from G4LogicalVolume.

Definition at line 64 of file G4LogicalCrystalVolume.hh.

64{ return true; }

◆ IsLattice()

G4bool G4LogicalCrystalVolume::IsLattice ( G4LogicalVolume * aLV)
static

Definition at line 62 of file G4LogicalCrystalVolume.cc.

63{
64 return std::find(fLCVvec.begin(), fLCVvec.end(), aLV) != fLCVvec.end();
65}

Referenced by G4Channeling::GetMeanFreePath(), and G4Channeling::PostStepDoIt().

◆ RotateToLattice()

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

Definition at line 122 of file G4LogicalCrystalVolume.cc.

123{
124 return dir.transform(fOrient);
125}
Hep3Vector & transform(const HepRotation &)

◆ RotateToSolid()

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

Definition at line 128 of file G4LogicalCrystalVolume.cc.

129{
130 return dir.transform(fInverse);
131}

Referenced by G4Channeling::PostStepDoIt().

◆ SetMillerOrientation()

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

Definition at line 85 of file G4LogicalCrystalVolume.cc.

89{
90 // Align Miller normal vector (hkl) with +Z axis, and rotation about axis
91
92 if (verboseLevel != 0)
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}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double phi() const
double theta() const
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
const G4ThreeVector & GetBasis(G4int i) const

Referenced by G4LogicalCrystalVolume().

◆ SetVerbose()

void G4LogicalCrystalVolume::SetVerbose ( G4int aInt)
inline

Definition at line 81 of file G4LogicalCrystalVolume.hh.

81{ verboseLevel = aInt; }

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