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

#include <G4LivermoreGammaConversionModel.hh>

Inheritance diagram for G4LivermoreGammaConversionModel:
Collaboration diagram for G4LivermoreGammaConversionModel:

Public Member Functions

 G4LivermoreGammaConversionModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermoreConversion")
 
virtual ~G4LivermoreGammaConversionModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 41 of file G4LivermoreGammaConversionModel.hh.

Constructor & Destructor Documentation

G4LivermoreGammaConversionModel::G4LivermoreGammaConversionModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermoreConversion" 
)

Definition at line 47 of file G4LivermoreGammaConversionModel.cc.

48 :G4VEmModel(nam),isInitialised(false),smallEnergy(2.*MeV)
49 {
50  fParticleChange = nullptr;
51 
52  lowEnergyLimit = 2.0*electron_mass_c2;
53 
54  verboseLevel= 0;
55  // Verbosity scale for debugging purposes:
56  // 0 = nothing
57  // 1 = calculation of cross sections, file openings...
58  // 2 = entering in methods
59 
60  if(verboseLevel > 0)
61  {
62  G4cout << "G4LivermoreGammaConversionModel is constructed " << G4endl;
63  }
64 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static constexpr double electron_mass_c2
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4LivermoreGammaConversionModel::~G4LivermoreGammaConversionModel ( )
virtual

Definition at line 68 of file G4LivermoreGammaConversionModel.cc.

69 {
70  if(IsMaster()) {
71  for(G4int i=0; i<maxZ; ++i) {
72  if(data[i]) {
73  delete data[i];
74  data[i] = 0;
75  }
76  }
77  }
78 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
int G4int
Definition: G4Types.hh:78
const XML_Char const XML_Char * data
Definition: expat.h:268

Member Function Documentation

G4double G4LivermoreGammaConversionModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 212 of file G4LivermoreGammaConversionModel.cc.

216 {
217  if (verboseLevel > 1)
218  {
219  G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModel"
220  << G4endl;
221  }
222 
223  if (GammaEnergy < lowEnergyLimit) { return 0.0; }
224 
225  G4double xs = 0.0;
226 
227  G4int intZ=G4int(Z);
228 
229  if(intZ < 1 || intZ > maxZ) { return xs; }
230 
231  G4LPhysicsFreeVector* pv = data[intZ];
232 
233  // if element was not initialised
234  // do initialisation safely for MT mode
235  if(!pv)
236  {
237  InitialiseForElement(0, intZ);
238  pv = data[intZ];
239  if(!pv) { return xs; }
240  }
241  // x-section is taken from the table
242  xs = pv->Value(GammaEnergy);
243 
244  if(verboseLevel > 0)
245  {
246  G4int n = pv->GetVectorLength() - 1;
247  G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
248  << GammaEnergy/MeV << G4endl;
249  G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
250  G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
251  G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
252  G4cout << "*********************************************************" << G4endl;
253  }
254 
255  return xs;
256 
257 }
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
const XML_Char const XML_Char * data
Definition: expat.h:268
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4LivermoreGammaConversionModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 82 of file G4LivermoreGammaConversionModel.cc.

85 {
86  if (verboseLevel > 1)
87  {
88  G4cout << "Calling Initialise() of G4LivermoreGammaConversionModel."
89  << G4endl
90  << "Energy range: "
91  << LowEnergyLimit() / MeV << " MeV - "
92  << HighEnergyLimit() / GeV << " GeV"
93  << G4endl;
94  }
95 
96  if(IsMaster())
97  {
98  // Initialise element selector
99  InitialiseElementSelectors(particle, cuts);
100 
101  // Access to elements
102  char* path = getenv("G4LEDATA");
103 
104  G4ProductionCutsTable* theCoupleTable =
106 
107  G4int numOfCouples = theCoupleTable->GetTableSize();
108 
109  for(G4int i=0; i<numOfCouples; ++i)
110  {
111  const G4Material* material =
112  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
113  const G4ElementVector* theElementVector = material->GetElementVector();
114  G4int nelm = material->GetNumberOfElements();
115 
116  for (G4int j=0; j<nelm; ++j)
117  {
118  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
119  if(Z < 1) { Z = 1; }
120  else if(Z > maxZ) { Z = maxZ; }
121  if(!data[Z]) { ReadData(Z, path); }
122  }
123  }
124  }
125  if(isInitialised) { return; }
126  fParticleChange = GetParticleChangeForGamma();
127  isInitialised = true;
128 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
std::vector< G4Element * > G4ElementVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
const XML_Char const XML_Char * data
Definition: expat.h:268
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const

Here is the call graph for this function:

void G4LivermoreGammaConversionModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 471 of file G4LivermoreGammaConversionModel.cc.

474 {
475  G4AutoLock l(&LivermoreGammaConversionModelMutex);
476  // G4cout << "G4LivermoreGammaConversionModel::InitialiseForElement Z= "
477  // << Z << G4endl;
478  if(!data[Z]) { ReadData(Z); }
479  l.unlock();
480 }
const XML_Char const XML_Char * data
Definition: expat.h:268

Here is the call graph for this function:

void G4LivermoreGammaConversionModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 132 of file G4LivermoreGammaConversionModel.cc.

134 {
136 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809

Here is the call graph for this function:

G4double G4LivermoreGammaConversionModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 141 of file G4LivermoreGammaConversionModel.cc.

144 {
145  return lowEnergyLimit;
146 }
void G4LivermoreGammaConversionModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 261 of file G4LivermoreGammaConversionModel.cc.

266 {
267 
268 // The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
269 // cross sections with Coulomb correction. A modified version of the random
270 // number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
271 
272 // Note 1 : Effects due to the breakdown of the Born approximation at low
273 // energy are ignored.
274 // Note 2 : The differential cross section implicitly takes account of
275 // pair creation in both nuclear and atomic electron fields. However triplet
276 // prodution is not generated.
277 
278  if (verboseLevel > 1) {
279  G4cout << "Calling SampleSecondaries() of G4LivermoreGammaConversionModel"
280  << G4endl;
281  }
282 
283  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
284  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
285 
286  G4double epsilon ;
287  G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
288 
289  // Do it fast if photon energy < 2. MeV
290  if (photonEnergy < smallEnergy )
291  {
292  epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
293  }
294  else
295  {
296  // Select randomly one element in the current material
297 
298  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
299  const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
300 
301  if (element == 0)
302  {
303  G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - element = 0"
304  << G4endl;
305  return;
306  }
307  G4IonisParamElm* ionisation = element->GetIonisation();
308  if (ionisation == 0)
309  {
310  G4cout << "G4LivermoreGammaConversionModel::SampleSecondaries - ionisation = 0"
311  << G4endl;
312  return;
313  }
314 
315  // Extract Coulomb factor for this Element
316  G4double fZ = 8. * (ionisation->GetlogZ3());
317  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
318 
319  // Limits of the screening variable
320  G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
321  G4double screenMax = G4Exp ((42.24 - fZ)/8.368) - 0.952 ;
322  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
323 
324  // Limits of the energy sampling
325  G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
326  G4double epsilonMin = std::max(epsilon0Local,epsilon1);
327  G4double epsilonRange = 0.5 - epsilonMin ;
328 
329  // Sample the energy rate of the created electron (or positron)
330  G4double screen;
331  G4double gReject ;
332 
333  G4double f10 = ScreenFunction1(screenMin) - fZ;
334  G4double f20 = ScreenFunction2(screenMin) - fZ;
335  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
336  G4double normF2 = std::max(1.5 * f20,0.);
337 
338  do
339  {
340  if (normF1 / (normF1 + normF2) > G4UniformRand() )
341  {
342  epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.333333) ;
343  screen = screenFactor / (epsilon * (1. - epsilon));
344  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
345  }
346  else
347  {
348  epsilon = epsilonMin + epsilonRange * G4UniformRand();
349  screen = screenFactor / (epsilon * (1 - epsilon));
350  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
351  }
352  } while ( gReject < G4UniformRand() );
353 
354  } // End of epsilon sampling
355 
356  // Fix charges randomly
357 
358  G4double electronTotEnergy;
359  G4double positronTotEnergy;
360 
361  if (G4UniformRand() > 0.5)
362  {
363  electronTotEnergy = (1. - epsilon) * photonEnergy;
364  positronTotEnergy = epsilon * photonEnergy;
365  }
366  else
367  {
368  positronTotEnergy = (1. - epsilon) * photonEnergy;
369  electronTotEnergy = epsilon * photonEnergy;
370  }
371 
372  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
373  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
374  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
375 
376  G4double u;
377  const G4double a1 = 0.625;
378  G4double a2 = 3. * a1;
379  // G4double d = 27. ;
380 
381  // if (9. / (9. + d) > G4UniformRand())
382  if (0.25 > G4UniformRand())
383  {
384  u = - G4Log(G4UniformRand() * G4UniformRand()) / a1 ;
385  }
386  else
387  {
388  u = - G4Log(G4UniformRand() * G4UniformRand()) / a2 ;
389  }
390 
391  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
392  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
393  G4double phi = twopi * G4UniformRand();
394 
395  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
396  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
397 
398 
399  // Kinematics of the created pair:
400  // the electron and positron are assumed to have a symetric angular
401  // distribution with respect to the Z axis along the parent photon
402 
403  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
404 
405  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
406  electronDirection.rotateUz(photonDirection);
407 
409  electronDirection,
410  electronKineEnergy);
411 
412  // The e+ is always created
413  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
414 
415  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
416  positronDirection.rotateUz(photonDirection);
417 
418  // Create G4DynamicParticle object for the particle2
420  positronDirection,
421  positronKineEnergy);
422  // Fill output vector
423  fvect->push_back(particle1);
424  fvect->push_back(particle2);
425 
426  // kill incident photon
427  fParticleChange->SetProposedKineticEnergy(0.);
428  fParticleChange->ProposeTrackStatus(fStopAndKill);
429 
430 }
G4double GetKineticEnergy() const
G4double GetfCoulomb() const
Definition: G4Element.hh:191
G4ParticleDefinition * GetDefinition() const
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4double GetlogZ3() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4Positron * Positron()
Definition: G4Positron.cc:94
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:199
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double GetZ3() const
double epsilon(double density, double temperature)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541

Here is the call graph for this function:


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