Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/coherent_elastic/src/G4NeutrinoElectronNcModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

Diff markup

Differences between /processes/hadronic/models/coherent_elastic/src/G4NeutrinoElectronNcModel.cc (Version 11.3.0) and /processes/hadronic/models/coherent_elastic/src/G4NeutrinoElectronNcModel.cc (Version 10.5)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
 26 //                                                 26 //
 27 // Geant4 Header : G4NeutrinoElectronNcModel       27 // Geant4 Header : G4NeutrinoElectronNcModel
 28 //                                                 28 //
 29 // Author : V.Grichine 6.4.17                      29 // Author : V.Grichine 6.4.17
 30 //                                                 30 //  
 31                                                    31 
 32 #include "G4NeutrinoElectronNcModel.hh"            32 #include "G4NeutrinoElectronNcModel.hh"
 33 #include "G4SystemOfUnits.hh"                      33 #include "G4SystemOfUnits.hh"
 34 #include "G4ParticleTable.hh"                      34 #include "G4ParticleTable.hh"
 35 #include "G4ParticleDefinition.hh"                 35 #include "G4ParticleDefinition.hh"
 36 #include "G4IonTable.hh"                           36 #include "G4IonTable.hh"
 37 #include "Randomize.hh"                            37 #include "Randomize.hh"
 38 #include "G4Electron.hh"                           38 #include "G4Electron.hh"
 39 #include "G4HadronicParameters.hh"                 39 #include "G4HadronicParameters.hh"
 40 #include "G4PhysicsModelCatalog.hh"            << 
 41                                                    40 
 42 using namespace std;                               41 using namespace std;
 43 using namespace CLHEP;                             42 using namespace CLHEP;
 44                                                    43 
 45 G4NeutrinoElectronNcModel::G4NeutrinoElectronN     44 G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(const G4String& name) 
 46   : G4HadronElastic(name)                          45   : G4HadronElastic(name)
 47 {                                                  46 {
 48   secID = G4PhysicsModelCatalog::GetModelID( " << 
 49                                                << 
 50   SetMinEnergy( 0.0*GeV );                         47   SetMinEnergy( 0.0*GeV );
 51   SetMaxEnergy( G4HadronicParameters::Instance     48   SetMaxEnergy( G4HadronicParameters::Instance()->GetMaxEnergy() );
 52   SetLowestEnergyLimit(1.e-6*eV);                  49   SetLowestEnergyLimit(1.e-6*eV);  
 53                                                    50 
 54   theElectron = G4Electron::Electron();            51   theElectron = G4Electron::Electron();
 55   // PDG2016: sin^2 theta Weinberg                 52   // PDG2016: sin^2 theta Weinberg
 56                                                    53 
 57   fSin2tW = 0.23129; // 0.2312;                    54   fSin2tW = 0.23129; // 0.2312;
 58                                                    55 
 59   fCutEnergy = 0.; // default value                56   fCutEnergy = 0.; // default value
                                                   >>  57 
 60 }                                                  58 }
 61                                                    59 
 62                                                    60 
 63 G4NeutrinoElectronNcModel::~G4NeutrinoElectron     61 G4NeutrinoElectronNcModel::~G4NeutrinoElectronNcModel()
 64 {}                                                 62 {}
 65                                                    63 
                                                   >>  64 
 66 void G4NeutrinoElectronNcModel::ModelDescripti     65 void G4NeutrinoElectronNcModel::ModelDescription(std::ostream& outFile) const
 67 {                                                  66 {
 68   outFile << "G4NeutrinoElectronNcModel is a n <<  67 
 69     << "model which uses the standard model \n <<  68     outFile << "G4NeutrinoElectronNcModel is a neutrino-electron (neutral current) elastic scattering\n"
 70     << "transfer parameterization.  The model  <<  69             << "model which uses the standard model \n"
                                                   >>  70             << "transfer parameterization.  The model is fully relativistic\n";
                                                   >>  71 
 71 }                                                  72 }
 72                                                    73 
 73 //////////////////////////////////////////////     74 /////////////////////////////////////////////////////////
 74                                                    75 
 75 G4bool G4NeutrinoElectronNcModel::IsApplicable <<  76 G4bool G4NeutrinoElectronNcModel::IsApplicable(const G4HadProjectile & aTrack, 
                                                   >>  77               G4Nucleus & targetNucleus)
 76 {                                                  78 {
 77   G4bool result  = false;                          79   G4bool result  = false;
 78   G4String pName = aTrack.GetDefinition()->Get     80   G4String pName = aTrack.GetDefinition()->GetParticleName();
 79   G4double minEnergy = 0.;                     <<  81   G4double minEnergy = 0., energy = aTrack.GetTotalEnergy();
 80   G4double energy = aTrack.GetTotalEnergy();   << 
 81                                                    82 
 82   if( fCutEnergy > 0. ) // min detected recoil     83   if( fCutEnergy > 0. ) // min detected recoil electron energy
 83   {                                                84   {
 84     minEnergy = 0.5*(fCutEnergy+sqrt(fCutEnerg     85     minEnergy = 0.5*(fCutEnergy+sqrt(fCutEnergy*(fCutEnergy+2.*electron_mass_c2)));
 85   }                                                86   }
 86   if( ( pName == "nu_e"   || pName == "anti_nu     87   if( ( pName == "nu_e"   || pName == "anti_nu_e"   || 
 87         pName == "nu_mu"  || pName == "anti_nu     88         pName == "nu_mu"  || pName == "anti_nu_nu"  || 
 88         pName == "nu_tau" || pName == "anti_nu     89         pName == "nu_tau" || pName == "anti_nu_tau"   ) &&
 89         energy > minEnergy                         90         energy > minEnergy                                 )
 90   {                                                91   {
 91     result = true;                                 92     result = true;
 92   }                                                93   }
                                                   >>  94   G4int Z = targetNucleus.GetZ_asInt();
                                                   >>  95         Z *= 1;
                                                   >>  96 
 93   return result;                                   97   return result;
 94 }                                                  98 }
 95                                                    99 
 96 //////////////////////////////////////////////    100 ////////////////////////////////////////////////
 97 //                                                101 //
 98 //                                                102 //
 99                                                   103 
100 G4HadFinalState* G4NeutrinoElectronNcModel::Ap    104 G4HadFinalState* G4NeutrinoElectronNcModel::ApplyYourself(
101      const G4HadProjectile& aTrack, G4Nucleus& << 105      const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
102 {                                                 106 {
103   theParticleChange.Clear();                      107   theParticleChange.Clear();
104                                                   108 
105   const G4HadProjectile* aParticle = &aTrack;     109   const G4HadProjectile* aParticle = &aTrack;
106   G4double nuTkin = aParticle->GetKineticEnerg    110   G4double nuTkin = aParticle->GetKineticEnergy();
107                                                   111 
108   if( nuTkin <= LowestEnergyLimit() )             112   if( nuTkin <= LowestEnergyLimit() ) 
109   {                                               113   {
110     theParticleChange.SetEnergyChange(nuTkin);    114     theParticleChange.SetEnergyChange(nuTkin);
111     theParticleChange.SetMomentumChange(aTrack    115     theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
112     return &theParticleChange;                    116     return &theParticleChange;
113   }                                               117   }
114   // sample and make final state in lab frame     118   // sample and make final state in lab frame
115                                                   119 
116   G4double eTkin = SampleElectronTkin( aPartic    120   G4double eTkin = SampleElectronTkin( aParticle );
117                                                   121 
118   if( eTkin > fCutEnergy )                        122   if( eTkin > fCutEnergy )
119   {                                               123   {
120     G4double ePlab = sqrt( eTkin*(eTkin + 2.*e    124     G4double ePlab = sqrt( eTkin*(eTkin + 2.*electron_mass_c2) );
121                                                   125 
122     G4double cost2  = eTkin*(nuTkin + electron    126     G4double cost2  = eTkin*(nuTkin + electron_mass_c2)*(nuTkin + electron_mass_c2);
123              cost2 /= nuTkin*nuTkin*(eTkin + 2    127              cost2 /= nuTkin*nuTkin*(eTkin + 2.*electron_mass_c2);
124                                                   128 
125     if( cost2 > 1. ) cost2 = 1.;                  129     if( cost2 > 1. ) cost2 = 1.;
126     if( cost2 < 0. ) cost2 = 0.;                  130     if( cost2 < 0. ) cost2 = 0.;
127                                                   131 
128     G4double cost = sqrt(cost2);                  132     G4double cost = sqrt(cost2);
129     G4double sint = std::sqrt( (1.0 - cost)*(1    133     G4double sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
130     G4double phi  = G4UniformRand()*CLHEP::two    134     G4double phi  = G4UniformRand()*CLHEP::twopi;
131                                                   135 
132     G4ThreeVector eP( sint*std::cos(phi), sint    136     G4ThreeVector eP( sint*std::cos(phi), sint*std::sin(phi), cost );
133     eP *= ePlab;                                  137     eP *= ePlab;
134     G4LorentzVector lvt2( eP, eTkin + electron    138     G4LorentzVector lvt2( eP, eTkin + electron_mass_c2 );
135     G4DynamicParticle * aSec = new G4DynamicPa    139     G4DynamicParticle * aSec = new G4DynamicParticle( theElectron, lvt2 );
136     theParticleChange.AddSecondary( aSec, secI << 140     theParticleChange.AddSecondary( aSec );
137                                                   141 
138     G4LorentzVector lvp1 = aParticle->Get4Mome    142     G4LorentzVector lvp1 = aParticle->Get4Momentum();
139     G4LorentzVector lvt1(0.,0.,0.,electron_mas    143     G4LorentzVector lvt1(0.,0.,0.,electron_mass_c2);
140     G4LorentzVector lvsum = lvp1+lvt1;            144     G4LorentzVector lvsum = lvp1+lvt1;
141                                                   145 
142     G4LorentzVector lvp2 = lvsum-lvt2;            146     G4LorentzVector lvp2 = lvsum-lvt2;
143     G4double nuTkin2 = lvp2.e()-aParticle->Get    147     G4double nuTkin2 = lvp2.e()-aParticle->GetDefinition()->GetPDGMass();
144     theParticleChange.SetEnergyChange(nuTkin2)    148     theParticleChange.SetEnergyChange(nuTkin2);
145     theParticleChange.SetMomentumChange(lvp2.v    149     theParticleChange.SetMomentumChange(lvp2.vect().unit());
146   }                                               150   }
147   else if( eTkin > 0.0 )                          151   else if( eTkin > 0.0 ) 
148   {                                               152   {
149     theParticleChange.SetLocalEnergyDeposit( e    153     theParticleChange.SetLocalEnergyDeposit( eTkin );
150     nuTkin -= eTkin;                              154     nuTkin -= eTkin;
151                                                   155 
152     if( nuTkin > 0. )                             156     if( nuTkin > 0. )
153     {                                             157     {
154       theParticleChange.SetEnergyChange( nuTki    158       theParticleChange.SetEnergyChange( nuTkin );
155       theParticleChange.SetMomentumChange( aTr    159       theParticleChange.SetMomentumChange( aTrack.Get4Momentum().vect().unit() );
156     }                                             160     }
157   }                                               161   }
158   else                                            162   else 
159   {                                               163   {
160     theParticleChange.SetEnergyChange( nuTkin     164     theParticleChange.SetEnergyChange( nuTkin );
161     theParticleChange.SetMomentumChange( aTrac    165     theParticleChange.SetMomentumChange( aTrack.Get4Momentum().vect().unit() );
162   }                                               166   }
                                                   >> 167   G4int Z = targetNucleus.GetZ_asInt();
                                                   >> 168         Z *= 1;
                                                   >> 169  
163   return &theParticleChange;                      170   return &theParticleChange;
164 }                                                 171 }
165                                                   172 
166 //////////////////////////////////////////////    173 //////////////////////////////////////////////////////
167 //                                                174 //
168 // sample recoil electron energy in lab frame     175 // sample recoil electron energy in lab frame
169                                                   176 
170 G4double G4NeutrinoElectronNcModel::SampleElec    177 G4double G4NeutrinoElectronNcModel::SampleElectronTkin(const G4HadProjectile* aParticle)
171 {                                                 178 {
172   G4double result = 0., xi, cofL, cofR, cofL2,    179   G4double result = 0., xi, cofL, cofR, cofL2, cofR2, cofLR;
173                                                   180 
174   G4double energy = aParticle->GetTotalEnergy(    181   G4double energy = aParticle->GetTotalEnergy();
175   if( energy == 0.) return result; // vmg: < t    182   if( energy == 0.) return result; // vmg: < th?? as in xsc 
176                                                   183 
177   G4String pName  = aParticle->GetDefinition()    184   G4String pName  = aParticle->GetDefinition()->GetParticleName();
178                                                   185 
179   if( pName == "nu_e")                            186   if( pName == "nu_e")
180   {                                               187   {
181     cofL = 0.5 + fSin2tW;                         188     cofL = 0.5 + fSin2tW;
182     cofR = fSin2tW;                               189     cofR = fSin2tW;
183   }                                               190   }
184   else if( pName == "anti_nu_e")                  191   else if( pName == "anti_nu_e")
185   {                                               192   {
186     cofL = fSin2tW;                               193     cofL = fSin2tW;
187     cofR = 0.5 + fSin2tW;                         194     cofR = 0.5 + fSin2tW;
188   }                                               195   }
189   else if( pName == "nu_mu")                      196   else if( pName == "nu_mu")
190   {                                               197   {
191     cofL = -0.5 + fSin2tW;                        198     cofL = -0.5 + fSin2tW;
192     cofR = fSin2tW;                               199     cofR = fSin2tW;
193   }                                               200   }
194   else if( pName == "anti_nu_mu")                 201   else if( pName == "anti_nu_mu")
195   {                                               202   {
196     cofL = fSin2tW;                               203     cofL = fSin2tW;
197     cofR = -0.5 + fSin2tW;                        204     cofR = -0.5 + fSin2tW;
198   }                                               205   }
199   else if( pName == "nu_tau") // vmg: nu_tau a    206   else if( pName == "nu_tau") // vmg: nu_tau as nu_mu ???
200   {                                               207   {
201     cofL = -0.5 + fSin2tW;                        208     cofL = -0.5 + fSin2tW;
202     cofR = fSin2tW;                               209     cofR = fSin2tW;
203   }                                               210   }
204   else if( pName == "anti_nu_tau")                211   else if( pName == "anti_nu_tau")
205   {                                               212   {
206     cofL = fSin2tW;                               213     cofL = fSin2tW;
207     cofR = -0.5 + fSin2tW;                        214     cofR = -0.5 + fSin2tW;
208   }                                               215   }
209   else                                            216   else
210   {                                               217   {
211     return result;                                218     return result;
212   }                                               219   }
213   xi = 0.5*electron_mass_c2/energy;               220   xi = 0.5*electron_mass_c2/energy;
214                                                   221 
215   cofL2 = cofL*cofL;                              222   cofL2 = cofL*cofL;
216   cofR2 = cofR*cofR;                              223   cofR2 = cofR*cofR;
217   cofLR = cofL*cofR;                              224   cofLR = cofL*cofR;
218                                                   225 
219   // cofs of Tkin/Enu 3rd equation                226   // cofs of Tkin/Enu 3rd equation
220                                                   227 
221   G4double a = cofR2/3.;                          228   G4double a = cofR2/3.;
222   G4double b = -(cofR2+cofLR*xi);                 229   G4double b = -(cofR2+cofLR*xi);
223   G4double c = cofL2+cofR2;                       230   G4double c = cofL2+cofR2;
224                                                   231 
225   G4double xMax  = 1./(1. + xi);                  232   G4double xMax  = 1./(1. + xi);
226   G4double xMax2 = xMax*xMax;                     233   G4double xMax2 = xMax*xMax;
227   G4double xMax3 = xMax*xMax2;                    234   G4double xMax3 = xMax*xMax2;
228                                                   235 
229   G4double d  = -( a*xMax3 + b*xMax2 + c*xMax     236   G4double d  = -( a*xMax3 + b*xMax2 + c*xMax );
230            d *= G4UniformRand();                  237            d *= G4UniformRand();
231                                                   238 
232   // G4cout<<a<<"   "<<b<<"   "<<c<<"   "<<d<<    239   // G4cout<<a<<"   "<<b<<"   "<<c<<"   "<<d<<G4endl<<G4endl;
233                                                   240 
234   // cofs of the incomplete 3rd equation          241   // cofs of the incomplete 3rd equation
235                                                   242 
236   G4double p  = c/a;                              243   G4double p  = c/a;
237            p -= b*b/a/a/3.;                       244            p -= b*b/a/a/3.;
238   G4double q  = d/a;                              245   G4double q  = d/a;
239            q -= b*c/a/a/3.;                       246            q -= b*c/a/a/3.;
240            q += 2*b*b*b/a/a/a/27.;                247            q += 2*b*b*b/a/a/a/27.;
241                                                   248 
242                                                   249 
243   // cofs for the incomplete colutions            250   // cofs for the incomplete colutions
244                                                   251 
245   G4double D  = p*p*p/3./3./3.;                   252   G4double D  = p*p*p/3./3./3.;
246            D += q*q/2./2.;                        253            D += q*q/2./2.;
247                                                   254 
248      // G4cout<<"D = "<<D<<G4endl;                255      // G4cout<<"D = "<<D<<G4endl;
249      // D = -D;                                   256      // D = -D;
250      // G4complex A1 = G4complex(- q/2., std::    257      // G4complex A1 = G4complex(- q/2., std::sqrt(-D) );
251      // G4complex A  = std::pow(A1,1./3.);        258      // G4complex A  = std::pow(A1,1./3.);
252                                                   259 
253      // G4complex B1 = G4complex(- q/2., -std:    260      // G4complex B1 = G4complex(- q/2., -std::sqrt(-D) );
254      // G4complex B  = std::pow(B1,1./3.);        261      // G4complex B  = std::pow(B1,1./3.);
255                                                   262 
256   G4double A1 = - q/2. + std::sqrt(D);            263   G4double A1 = - q/2. + std::sqrt(D);
257   G4double A = std::pow(A1,1./3.);                264   G4double A = std::pow(A1,1./3.);
258                                                   265 
259   G4double B1 = - q/2. - std::sqrt(D);            266   G4double B1 = - q/2. - std::sqrt(D);
260   G4double B = std::pow(-B1,1./3.);               267   G4double B = std::pow(-B1,1./3.);
261            B = -B;                                268            B = -B;
262                                                   269 
263   // roots of the incomplete 3rd equation         270   // roots of the incomplete 3rd equation
264                                                   271 
265   G4complex y1 =  A + B;                          272   G4complex y1 =  A + B;
266   // G4complex y2 = -0.5*(A + B) + 0.5*std::sq    273   // G4complex y2 = -0.5*(A + B) + 0.5*std::sqrt(3.)*(A - B)*G4complex(0.,1.);
267   // G4complex y3 = -0.5*(A + B) - 0.5*std::sq    274   // G4complex y3 = -0.5*(A + B) - 0.5*std::sqrt(3.)*(A - B)*G4complex(0.,1.);
268                                                   275  
269   G4complex x1 = y1 - b/a/3.;                     276   G4complex x1 = y1 - b/a/3.;
270   // G4complex x2 = y2 - b/a/3.;                  277   // G4complex x2 = y2 - b/a/3.;
271   // G4complex x3 = y3 - b/a/3.;                  278   // G4complex x3 = y3 - b/a/3.;
272                                                   279 
273   // G4cout<<"re_x1 = "<<real(x1)<<"; re_x2 =     280   // G4cout<<"re_x1 = "<<real(x1)<<"; re_x2 = "<<real(x2)<<"; re_x3 = "<<real(x3)<<G4endl;
274   // G4cout<<"im_x1 = "<<imag(x1)<<"; im_x2 =     281   // G4cout<<"im_x1 = "<<imag(x1)<<"; im_x2 = "<<imag(x2)<<"; im_x3 = "<<imag(x3)<<G4endl<<G4endl;
275                                                   282 
276   result = real(x1)*energy;                       283   result = real(x1)*energy;
277                                                   284 
278   return result;                                  285   return result;
279 }                                                 286 }
280                                                   287 
281 //                                                288 //
282 //                                                289 //
283 ///////////////////////////                       290 ///////////////////////////
284                                                   291