Public Member Functions | Private Member Functions | Private Attributes | List of all members
genie::ZExpAxialFormFactorModel Class Reference

Concrete implementation of the AxialFormFactorModelI interface. Computes the axial form factor using the model-independent z-expansion technique. More...

#include <ZExpAxialFormFactorModel.h>

Inheritance diagram for genie::ZExpAxialFormFactorModel:
genie::AxialFormFactorModelI genie::Algorithm

Public Member Functions

 ZExpAxialFormFactorModel ()
 
 ZExpAxialFormFactorModel (string config)
 
virtual ~ZExpAxialFormFactorModel ()
 
double FA (const Interaction *interaction) const
 Compute the axial form factor. More...
 
void Configure (const Registry &config)
 
void Configure (string param_set)
 
- Public Member Functions inherited from genie::AxialFormFactorModelI
virtual ~AxialFormFactorModelI ()
 
- Public Member Functions inherited from genie::Algorithm
virtual ~Algorithm ()
 
virtual void FindConfig (void)
 
virtual const RegistryGetConfig (void) const
 
RegistryGetOwnedConfig (void)
 
virtual const AlgIdId (void) const
 Get algorithm ID. More...
 
virtual AlgStatus_t GetStatus (void) const
 Get algorithm status. More...
 
virtual bool AllowReconfig (void) const
 
virtual AlgCmp_t Compare (const Algorithm *alg) const
 Compare with input algorithm. More...
 
virtual void SetId (const AlgId &id)
 Set algorithm ID. More...
 
virtual void SetId (string name, string config)
 
const AlgorithmSubAlg (const RgKey &registry_key) const
 
void AdoptConfig (void)
 
void AdoptSubstructure (void)
 
virtual void Print (ostream &stream) const
 Print algorithm info. More...
 

Private Member Functions

double CalculateZ (double q2) const
 
void FixCoeffs (void)
 
void FixA0 (void)
 
void FixQ4Limit (void)
 
void LoadConfig (void)
 

Private Attributes

bool fQ4limit
 
int fKmax
 
double fT0
 
double fTcut
 
double fFA0
 
double * fZ_An
 

Additional Inherited Members

- Static Public Member Functions inherited from genie::Algorithm
static string BuildParamVectKey (const std::string &comm_name, unsigned int i)
 
static string BuildParamVectSizeKey (const std::string &comm_name)
 
- Protected Member Functions inherited from genie::AxialFormFactorModelI
 AxialFormFactorModelI ()
 
 AxialFormFactorModelI (string name)
 
 AxialFormFactorModelI (string name, string config)
 
- Protected Member Functions inherited from genie::Algorithm
 Algorithm ()
 
 Algorithm (string name)
 
 Algorithm (string name, string config)
 
void Initialize (void)
 
void DeleteConfig (void)
 
void DeleteSubstructure (void)
 
RegistryExtractLocalConfig (const Registry &in) const
 
RegistryExtractLowerConfig (const Registry &in, const string &alg_key) const
 Split an incoming configuration Registry into a block valid for the sub-algo identified by alg_key. More...
 
template<class T >
bool GetParam (const RgKey &name, T &p, bool is_top_call=true) const
 
template<class T >
bool GetParamDef (const RgKey &name, T &p, const T &def) const
 
template<class T >
int GetParamVect (const std::string &comm_name, std::vector< T > &v, bool is_top_call=true) const
 Handle to load vectors of parameters. More...
 
int GetParamVectKeys (const std::string &comm_name, std::vector< RgKey > &k, bool is_top_call=true) const
 
int AddTopRegistry (Registry *rp, bool owns=true)
 add registry with top priority, also update ownership More...
 
int AddLowRegistry (Registry *rp, bool owns=true)
 add registry with lowest priority, also update ownership More...
 
int MergeTopRegistry (const Registry &r)
 
int AddTopRegisties (const vector< Registry * > &rs, bool owns=false)
 Add registries with top priority, also udated Ownerships. More...
 
- Protected Attributes inherited from genie::Algorithm
bool fAllowReconfig
 
bool fOwnsSubstruc
 true if it owns its substructure (sub-algs,...) More...
 
AlgId fID
 algorithm name and configuration set More...
 
vector< Registry * > fConfVect
 
vector< boolfOwnerships
 ownership for every registry in fConfVect More...
 
AlgStatus_t fStatus
 algorithm execution status More...
 
AlgMapfOwnedSubAlgMp
 local pool for owned sub-algs (taken out of the factory pool) More...
 

Detailed Description

Concrete implementation of the AxialFormFactorModelI interface. Computes the axial form factor using the model-independent z-expansion technique.

Hill et al. arXiv:1008.4619 DOI: 10.1103/PhysRevD.82.113005

Author
Aaron Meyer <asmeyer2012 uchicago.edu>
      based off DipoleELFormFactorsModel by
      Costas Andreopoulos <constantinos.andreopoulos \at cern.ch>
      STFC, Rutherford Appleton Laboratory

August 16, 2013

Copyright (c) 2003-2020, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 35 of file ZExpAxialFormFactorModel.h.

Constructor & Destructor Documentation

ZExpAxialFormFactorModel::ZExpAxialFormFactorModel ( )

Definition at line 31 of file ZExpAxialFormFactorModel.cxx.

31  :
32 AxialFormFactorModelI("genie::ZExpAxialFormFactorModel")
33 {
34 
35 }
ZExpAxialFormFactorModel::ZExpAxialFormFactorModel ( string  config)

Definition at line 37 of file ZExpAxialFormFactorModel.cxx.

37  :
38 AxialFormFactorModelI("genie::ZExpAxialFormFactorModel", config)
39 {
40 
41 }
static Config * config
Definition: config.cpp:1054
ZExpAxialFormFactorModel::~ZExpAxialFormFactorModel ( )
virtual

Definition at line 43 of file ZExpAxialFormFactorModel.cxx.

44 {
45  delete[] fZ_An;
46 }

Member Function Documentation

double ZExpAxialFormFactorModel::CalculateZ ( double  q2) const
private

Definition at line 67 of file ZExpAxialFormFactorModel.cxx.

68 {
69 
70  // calculate z expansion parameter
71  double znum = TMath::Sqrt(fTcut - q2) - TMath::Sqrt(fTcut - fT0);
72  double zden = TMath::Sqrt(fTcut - q2) + TMath::Sqrt(fTcut - fT0);
73 
74  return znum/zden;
75 }
void ZExpAxialFormFactorModel::Configure ( const Registry config)
virtual

Configure the algorithm with an external registry The registry is merged with the top level registry if it is owned, Otherwise a copy of it is added with the highest priority

Reimplemented from genie::Algorithm.

Definition at line 205 of file ZExpAxialFormFactorModel.cxx.

206 {
207  Algorithm::Configure(config);
208  this->LoadConfig();
209 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
void ZExpAxialFormFactorModel::Configure ( string  config)
virtual

Configure the algorithm from the AlgoConfigPool based on param_set string given in input An algorithm contains a vector of registries coming from different xml configuration files, which are loaded according a very precise prioriy This methods will load a number registries in order of priority: 1) "Tunable" parameter set from CommonParametes. This is loaded with the highest prioriry and it is designed to be used for tuning procedure Usage not expected from the user. 2) For every string defined in "CommonParame" the corresponding parameter set will be loaded from CommonParameter.xml 3) parameter set specified by the config string and defined in the xml file of the algorithm 4) if config is not "Default" also the Default parameter set from the same xml file will be loaded Effectively this avoids the repetion of a parameter when it is not changed in the requested configuration

Reimplemented from genie::Algorithm.

Definition at line 211 of file ZExpAxialFormFactorModel.cxx.

212 {
213  Algorithm::Configure(param_set);
214  this->LoadConfig();
215 }
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
double ZExpAxialFormFactorModel::FA ( const Interaction interaction) const
virtual

Compute the axial form factor.

Implements genie::AxialFormFactorModelI.

Definition at line 48 of file ZExpAxialFormFactorModel.cxx.

49 {
50  // calculate and return FA
51  double q2 = interaction->KinePtr()->q2();
52  double zparam = this->CalculateZ(q2);
53  if (zparam != zparam) // checks for nan
54  {
55  LOG("ZExpAxialFormFactorModel",pWARN) << "Undefined expansion parameter";
56  return 0.;
57  }
58  double fa = 0.;
59  for (int ki=0;ki<=fKmax+(fQ4limit ? 4 : 0);ki++)
60  {
61  fa = fa + TMath::Power(zparam,ki) * fZ_An[ki];
62  }
63 
64  return fa;
65 }
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
double q2(bool selected=false) const
Definition: Kinematics.cxx:141
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pWARN
Definition: Messenger.h:60
void ZExpAxialFormFactorModel::FixA0 ( void  )
private

Definition at line 86 of file ZExpAxialFormFactorModel.cxx.

87 {
88  // Function to fix form factor such that FA(q2=0) = gA
89  // For T0 = 0, this will set A0 = gA
90  double zparam = this->CalculateZ(0.);
91  double fa = 0.;
92  for (int ki=1;ki<=fKmax;ki++)
93  {
94  fa = fa + TMath::Power(zparam,ki) * fZ_An[ki];
95  }
96  fZ_An[0] = fFA0 - fa;
97 
98 }
void ZExpAxialFormFactorModel::FixCoeffs ( void  )
private

Definition at line 77 of file ZExpAxialFormFactorModel.cxx.

78 {
79  //if (fKmax < 1 ) { fKmax = 1; }
80  //else if (fKmax > 10) { fKmax = 10; }
81 
82  if (fQ4limit) this->FixQ4Limit();
83  else this->FixA0();
84 }
void ZExpAxialFormFactorModel::FixQ4Limit ( void  )
private

Definition at line 100 of file ZExpAxialFormFactorModel.cxx.

101 {
102  // fixes parameters such that the model limits to 1/Q^4 at large t
103  // -- requires at least 5 parameters to do so --
104  // 4 parameters for Q^4 behavior, 1 for normalization to FA(q2=0)=gA
105 
106  // will use A_0 and A_Kmax through A_Kmax-3 to do the fitting
107  // calculate some useful numbers (redundancy for readability)
108  double kp4 = (double)fKmax+4;
109  double kp3 = (double)fKmax+3;
110  double kp2 = (double)fKmax+2;
111  double kp1 = (double)fKmax+1;
112  double kp0 = (double)fKmax ;
113  //double km5 = (double)fKmax-5;
114  double z0 = this->CalculateZ(0.);
115  double zkp4 = TMath::Power(z0,(int)kp4);
116  double zkp3 = TMath::Power(z0,(int)kp3);
117  double zkp2 = TMath::Power(z0,(int)kp2);
118  double zkp1 = TMath::Power(z0,(int)kp1);
119 
120  // denominator
121  double denom = \
122  6. - kp4*kp3*kp2*zkp1 + 3.*kp4*kp3*kp1*zkp2 \
123  - 3.*kp4*kp2*kp1*zkp3 + kp3*kp2*kp1*zkp4;
124 
125  // extra parameters (effectively constants)
126  // number refers to the number of derivatives
127  double b0 = 0.;
128  for (int ki = 1;ki <= fKmax;ki++)
129  {
130  b0 = b0 + fZ_An[ki];
131  }
132 
133  double b0z = -fFA0;
134  for (int ki = 1;ki <= fKmax;ki++)
135  {
136  b0z = b0z + fZ_An[ki]*TMath::Power(z0,ki);
137  }
138 
139  double b1 = 0.;
140  for (int ki = 1;ki <= fKmax;ki++)
141  {
142  b1 = b1 + ki*fZ_An[ki];
143  }
144 
145  double b2 = 0.;
146  for (int ki = 1;ki <= fKmax;ki++)
147  {
148  b2 = b2 + ki*(ki-1)*fZ_An[ki];
149  }
150 
151  double b3 = 0.;
152  for (int ki = 1;ki <= fKmax;ki++)
153  {
154  b3 = b3 + ki*(ki-1)*(ki-2)*fZ_An[ki];
155  }
156 
157  // Assign new parameters
158  fZ_An[(int)kp4] = (1./denom) * \
159  ( (b0-b0z)*kp3*kp2*kp1 \
160  + b3*( -1. + .5*kp3*kp2*zkp1 - kp3*kp1*zkp2 \
161  +.5*kp2*kp1*zkp3 ) \
162  + b2*( 3.*kp1 - kp3*kp2*kp1*zkp1 \
163  +kp3*kp1*(2*fKmax+1)*zkp2 - kp2*kp1*kp0*zkp3 ) \
164  + b1*( -3.*kp2*kp1 + .5*kp3*kp2*kp2*kp1*zkp1 \
165  -kp3*kp2*kp1*kp0*zkp2 + .5*kp2*kp1*kp1*kp0*zkp3) );
166 
167  fZ_An[(int)kp3] = (1./denom) * \
168  ( -3.*(b0-b0z)*kp4*kp2*kp1 \
169  + b3*( 3. - kp4*kp2*zkp1 + (3./2.)*kp4*kp1*zkp2 \
170  -.5*kp2*kp1*zkp4 ) \
171  + b2*( -3.*(3*fKmax+4) + kp4*kp2*(2*fKmax+3)*zkp1 \
172  -3.*kp4*kp1*kp1*zkp2 + kp2*kp1*kp0*zkp4 ) \
173  + b1*( 3.*kp1*(3*fKmax+8) - kp4*kp3*kp2*kp1*zkp1 \
174  +(3./2.)*kp4*kp3*kp1*kp0*zkp2 - .5*kp2*kp1*kp1*kp0*zkp4) );
175 
176  fZ_An[(int)kp2] = (1./denom) * \
177  ( 3.*(b0-b0z)*kp4*kp3*kp1 \
178  + b3*( -3. + .5*kp4*kp3*zkp1 - (3./2.)*kp4*kp1*zkp3 \
179  +kp3*kp1*zkp4 ) \
180  + b2*( 3.*(3*fKmax+5) - kp4*kp3*kp2*zkp1 \
181  +3.*kp4*kp1*kp1*zkp3 - kp3*kp1*(2*fKmax+1)*zkp4) \
182  + b1*( -3.*kp3*(3*fKmax+4) + .5*kp4*kp3*kp3*kp2*zkp1 \
183  -(3./2.)*kp4*kp3*kp1*kp0*zkp3 + kp3*kp2*kp1*kp0*zkp4) );
184 
185  fZ_An[(int)kp1] = (1./denom) * \
186  ( -(b0-b0z)*kp4*kp3*kp2 \
187  + b3*( 1. - .5*kp4*kp3*zkp2 + kp4*kp2*zkp3 \
188  -.5*kp3*kp2*zkp4 ) \
189  + b2*( -3.*kp2 + kp4*kp3*kp2*zkp2 \
190  -kp4*kp2*(2*fKmax+3)*zkp3 + kp3*kp2*kp1*zkp4) \
191  + b1*( 3.*kp3*kp2 - .5*kp4*kp3*kp3*kp2*zkp2 \
192  +kp4*kp3*kp2*kp1*zkp3 - .5*kp3*kp2*kp2*kp1*zkp4) );
193 
194  fZ_An[0] = (1./denom) * \
195  ( -6.*b0z \
196  + b0*( kp4*kp3*kp2*zkp1 - 3.*kp4*kp3*kp1*zkp2 \
197  +3.*kp4*kp2*kp1*zkp3 - kp3*kp2*kp1*zkp4 ) \
198  + b3*( -zkp1 + 3.*zkp2 - 3.*zkp3 + zkp4 ) \
199  + b2*( 3.*kp2*zkp1 - 3.*(3*fKmax+5)*zkp2 \
200  +3.*(3*fKmax+4)*zkp3 - 3.*kp1*zkp4 ) \
201  + b1*( -3.*kp3*kp2*zkp1 + 3.*kp3*(3*fKmax+4)*zkp2 \
202  -3.*kp1*(3*fKmax+8)*zkp3 + 3.*kp2*kp1*zkp4 ) );
203 }
void ZExpAxialFormFactorModel::LoadConfig ( void  )
private

Definition at line 217 of file ZExpAxialFormFactorModel.cxx.

218 {
219 // get config options from the configuration registry or set defaults
220 // from the global parameter list
221 
222  GetParam( "QEL-Q4limit", fQ4limit ) ;
223  GetParam( "QEL-Kmax", fKmax ) ;
224 
225  GetParam( "QEL-T0", fT0 ) ;
226  GetParam( "QEL-T0", fT0 ) ;
227  GetParam( "QEL-Tcut", fTcut ) ;
228 
229  GetParam( "QEL-FA0", fFA0 ) ;
230  assert(fKmax > 0);
231 
232  // z expansion coefficients
233  if (fQ4limit) fZ_An = new double [fKmax+5];
234  else fZ_An = new double [fKmax+1];
235 
236  // load the user-defined coefficient values
237  // -- A0 and An for n<fKmax are calculated from other means
238  for (int ip=1;ip<fKmax+1;ip++) {
239  ostringstream alg_key;
240  alg_key << "QEL-Z_A" << ip;
241  GetParam( alg_key.str(), fZ_An[ip] ) ;
242  }
243 
244  this->FixCoeffs();
245 }
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const

Member Data Documentation

double genie::ZExpAxialFormFactorModel::fFA0
private

Definition at line 62 of file ZExpAxialFormFactorModel.h.

int genie::ZExpAxialFormFactorModel::fKmax
private

Definition at line 59 of file ZExpAxialFormFactorModel.h.

bool genie::ZExpAxialFormFactorModel::fQ4limit
private

Definition at line 58 of file ZExpAxialFormFactorModel.h.

double genie::ZExpAxialFormFactorModel::fT0
private

Definition at line 60 of file ZExpAxialFormFactorModel.h.

double genie::ZExpAxialFormFactorModel::fTcut
private

Definition at line 61 of file ZExpAxialFormFactorModel.h.

double* genie::ZExpAxialFormFactorModel::fZ_An
private

Definition at line 64 of file ZExpAxialFormFactorModel.h.


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