Main Page | Namespace List | Class Hierarchy | Class List | File List | Class Members | File Members | Related Pages

charge.cpp

Go to the documentation of this file.
00001 /*
00002  * ******** fete: From ENDF To ENDL *********
00003  * 
00004  * Copyright (c) 2006, The Regents of the University of California. 
00005  * All rights reserved.
00006  * 
00007  * Produced at the Lawrence Livermore National Laboratory. 
00008  * Written by David A. Brown, Gerry Hedstrom, Tony Hill
00009  * 
00010  * This file is part of fete v1.0  (UCRL-CODE-218718)
00011  * 
00012  * Please read the COPYING file for "Our Notice and GNU General 
00013  * Public License" in the root of this software distribution.  
00014  * 
00015  * This program is free software; you can redistribute it and/or modify 
00016  * it under the terms of the GNU General Public License (as published by 
00017  * the Free Software Foundation) version 2, dated June 1991. 
00018  * 
00019  * This program is distributed in the hope that it will be useful, 
00020  * but WITHOUT ANY WARRANTY; without even the IMPLIED WARRANTY OF 
00021  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the terms 
00022  * and conditions of the GNU General Public License for more details. 
00023  * 
00024  * You should have received a copy of the GNU General Public License along 
00025  * with this program; if not, write to the Free Software Foundation, Inc., 
00026  * 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 
00027  * 
00028  * $Revision: 1747 $
00029  * $Date:$
00030  * $Author:$
00031  * $Id:$
00032  * 
00033  * ******** fete: From ENDF To ENDL *********
00034  */
00035 
00036 // Implementation of the classes used to translate the charged-particle
00037 // scattering data for identical particles.
00038 
00039 #include "charge.hpp"
00040 #include "endl_formats.hpp"
00041 #include "global_params.hpp"
00042 #include "bdfls_tools.hpp"
00043 
00044 extern ENDLClass ENDL;
00045 extern GlobalParameterClass Global;
00046 extern bdflsClass bdfls;  // for the data constants
00047 
00048 // ******************* class one_d_charge ******************
00049 // ------------ one_d_charge::f -----------------------
00050 double one_d_charge::f( double mu )
00051 {
00052   if( ( mu < -1.0 ) || ( mu >= 1.0 ) ||
00053       ( same_ && ( mu < 0.0 ) ) )
00054   {
00055     SevereError( "one_d_charge::f",
00056          pastenum("mu out of range: ", mu ) );
00057   }
00058 
00059   // the different terms
00060   double coulomb;  // Coulomb portion
00061   double interfere;   // the electron interference term
00062   double nuke;     // the nuclear reaction term
00063 
00064   if( LTP_model == 1 )
00065   {
00066     // are the particles identical?
00067     if( same_ )
00068     {
00069       coulomb = coulomb_same( mu );
00070       interfere = interfere_same( mu );
00071       nuke = nuclear_same( mu );
00072     }
00073     else
00074     {
00075       coulomb = coulomb_diff( mu );
00076       interfere = interfere_diff( mu );
00077       nuke = nuclear_diff( mu );
00078     }
00079     return coulomb - interfere + nuke;
00080   }
00081   else
00082   {
00083     // are the particles identical?
00084     if( same_ )
00085     {
00086       coulomb = coulomb_same( mu );
00087     }
00088     else
00089     {
00090       coulomb = coulomb_diff( mu );
00091     }
00092     return coulomb;
00093   }
00094 }
00095 // ------------ one_d_charge::coulomb_same -----------------------
00096 double one_d_charge::coulomb_same( double mu )
00097 // 1/2 the first (Coulomb) term from (6.8) for identical particles
00098 {
00099   double term = ( 2 * eta * eta /( k_sq * (1-mu)*(1-mu) ) ) *
00100     ( ( 1 + mu*mu ) / ( 1 - mu*mu ) +
00101       ( spin_factor / (2*spin_ + 1) ) *
00102       cos( eta * log( ( 1 + mu ) / ( 1 - mu ) ) ) );
00103 
00104   return 0.5*term;
00105 }
00106 // ------------ one_d_charge::interfere_same -----------------------
00107 double one_d_charge::interfere_same( double mu )
00108 // 1/2 the second (interference) term from (6.12) for identical particles
00109 {
00110   complex< double > exponent( 0.0, eta * log( ( 1 - mu ) / 2 ) );
00111   complex< double > term_1 = ( 1 + mu ) * exp( exponent );
00112 
00113   complex< double > expon_2( 0.0, eta * log( ( 1 + mu ) / 2 ) );
00114   complex< double > term_2 = ( 1 - mu ) * exp( expon_2 );
00115 
00116   // the zeroth-order term
00117   complex< double > sum = 0.5* ( term_1 + term_2 ) * a_coef[0];
00118 
00119   // the first-order term
00120   if( order > 0 )
00121   {
00122     sum += 1.5 * ( term_1 - term_2 ) * a_coef[1] * mu;
00123   }
00124 
00125   // the remaining terms
00126   double sign = 1.0;    // int won't work
00127   double P_prev = 1.0;  // P_0(mu)
00128   double P_this = mu;   // P_1(mu)
00129   double P_next;
00130   for( int count = 2; count <= order; ++count )
00131   {
00132     P_next = ((2*count - 1)*mu*P_this - (count-1)*P_prev)/count;
00133     complex< double > this_coef = a_coef[count] *
00134       (count + 0.5)*( term_1 + sign*term_2 )*P_next;
00135     sum += this_coef;
00136 
00137     // shift the values
00138     P_prev = P_this;
00139     P_this = P_next;
00140     sign *= -1;
00141   }
00142     
00143   return eta / ( 1 - mu*mu ) * (sum.real( ));
00144 }
00145 // ------------ one_d_charge::nuclear_same -----------------------
00146 double one_d_charge::nuclear_same( double mu )
00147 // 1/2 the third (nuclear) term from (6.12) for identical particles
00148 {
00149 
00150   // the zeroth-order term
00151   double sum =b_coef[0] / 2;
00152 
00153   // the remaining (even) terms
00154   double P_prev = 1.0;  // P_0(mu)
00155   double P_this = mu;   // P_1(mu)
00156   double P_next;
00157   for( int count = 2; count <= order; ++count )
00158   {
00159     P_next = ((2*count - 1)*mu*P_this - (count-1)*P_prev)/count;
00160     if( count % 2 == 0 )
00161     {
00162       sum += (count + 0.5)*P_next*b_coef[ count/2 ];
00163     }
00164     // shift the values
00165     P_prev = P_this;
00166     P_this = P_next;
00167   }
00168     
00169   return 0.5*sum;
00170 }
00171 // ------------ one_d_charge::coulomb_diff -----------------------
00172 double one_d_charge::coulomb_diff( double mu )
00173 // the first (Coulomb) term from (6.7) for different particles
00174 {
00175   double term = eta * eta /( k_sq * (1-mu)*(1-mu) );
00176 
00177   return term;
00178 }
00179 // ------------ one_d_charge::interfere_diff -----------------------
00180 double one_d_charge::interfere_diff( double mu )
00181 // the second (interference) term from (6.11) for different particles
00182 {
00183   // the zeroth-order term
00184   complex< double > sum = 0.5 * a_coef[0];
00185 
00186   // the first-order term
00187   if( order > 0 )
00188   {
00189     sum += 1.5 * a_coef[1] * mu;
00190   }
00191 
00192   // the remaining terms
00193   double P_prev = 1.0;  // P_0(mu)
00194   double P_this = mu;   // P_1(mu)
00195   double P_next;
00196   for( int count = 2; count <= order; ++count )
00197   {
00198     P_next = ((2*count - 1)*mu*P_this - (count-1)*P_prev)/count;
00199     complex< double > this_coef = a_coef[count] *
00200       (count + 0.5) * P_next;
00201     sum += this_coef;
00202 
00203     // shift the values
00204     P_prev = P_this;
00205     P_this = P_next;
00206   }
00207   
00208   complex< double > exponent( 0.0, eta * log( ( 1 - mu ) / 2 ) );
00209   sum *= exp( exponent );
00210   return 2 * eta / ( 1 - mu ) * sum.real( );
00211 }
00212 // ------------ one_d_charge::nuclear_diff -----------------------
00213 double one_d_charge::nuclear_diff( double mu )
00214 // the third (nuclear) term from (6.11) for different particles
00215 {
00216 
00217   // the zeroth-order term
00218   double sum = b_coef[0] / 2;
00219 
00220   // the first-order term
00221   if( order > 0 )
00222   {
00223     sum += 1.5 * b_coef[1] * mu;
00224   }
00225 
00226   // the remaining (even) terms
00227   double P_prev = 1.0;  // P_0(mu)
00228   double P_this = mu;   // P_1(mu)
00229   double P_next;
00230   for( int count = 2; count <= 2*order; ++count )
00231   {
00232     P_next = ((2*count - 1)*mu*P_this - (count-1)*P_prev)/count;
00233     sum += (count + 0.5)*P_next*b_coef[ count ];
00234 
00235     // shift the values
00236     P_prev = P_this;
00237     P_this = P_next;
00238   }
00239     
00240   return sum;
00241 }
00242 
00243 // ------------ one_d_charge::initiate -----------------------
00244 void one_d_charge::initiate( mf6_file& inFile, double e_in, double spin )
00245 {
00246   // are the particles the same?
00247   int ZA_in = ENDL.yo_to_za( ENDL.incident_particle );
00248   same_ = ( ZA_in == ENDL.ZA );
00249 
00250   // Calculate the parameter eta
00251   E_in( ) = e_in;
00252   // the charge
00253   int Z_in = ENDL.get_Z( ENDL.incident_particle );
00254   int Z_targ = ENDL.ZA / 1000;
00255 
00256   // square of parameter k from (6.9)
00257   // But the correct formula is
00258   //  k = (A/(A+1)) * sqrt(2*m_1*E) / (h_bar*c)
00259   // We use m_ and E in MeV, h_bar in MeV sh, c in cm/sh
00260   // h_bar = bdfls.nuclear_constants[9] in MeV sh
00261   // c = bdfls.nuclear_constants[8] in cm/sh
00262   double denom = bdfls.nuclear_constants[8] *
00263     bdfls.nuclear_constants[9];
00264   double amu_to_MeV = bdfls.nuclear_constants[5];
00265   double mass_ratio = ENDL.target_mass / ENDL.projectile_mass;
00266   // the 10^{24} is to convert cm^2 to barns
00267   k_sq = ( mass_ratio / ( mass_ratio + 1 ) ) * ENDL.projectile_mass * 
00268     amu_to_MeV * e_in /
00269     ( denom * denom * 1.0e24 );
00270 
00271   // eta parameter from (6.10)
00272   // If h_bar = 1 and c = 1, then e^2 = alpha (fine-structure constant)
00273   // value from physics.nist.gov
00274   double alpha = 7.297352568e-3;
00275   eta = Z_in * Z_targ * alpha *
00276     sqrt( 0.5 * amu_to_MeV * ENDL.projectile_mass / e_in );
00277 
00278   // (-1)^{2*spin}
00279   double ds = ENDL_EPSILON( 1.0 );
00280   // is 2*spin even or odd?
00281   spin_factor = ( static_cast<int>( 2 * spin * (1+ds) ) % 2 == 0 ) ?
00282     1 : -1;
00283 
00284   // save the spin
00285   spin_ = spin;
00286 
00287   // the amount of data depends on the model
00288   if( LTP_model == 1 )
00289   {
00290     a_coef.reserve( order + 1 );
00291     a_coef.resize( order + 1 );
00292     if( same_ )
00293     {
00294       b_coef.reserve( order + 1 );
00295       b_coef.resize( order + 1 );
00296     }
00297     else
00298     {
00299       b_coef.reserve( 2*order + 1 );
00300       b_coef.resize( 2*order + 1 );
00301     }
00302   }
00303   else
00304   {
00305   }
00306 }
00307 // ----------- one_d_charge::read_data -----------
00308 void one_d_charge::read_data( mf6_file& inFile )
00309 {
00310   string linebuff;  // the input line
00311   string strbuff;   // a substring
00312   int num_data;
00313   int num_b;
00314   int count;
00315   int sub_count;
00316 
00317   if( LTP_model == 1 )
00318   {
00319     if( same_ )
00320     {
00321       num_b = order + 1;
00322       num_data = 3*order + 3;
00323     }
00324     else
00325     {
00326       num_b = 2*order + 1;
00327       num_data = 4*order + 3;
00328     }
00329     // read the b-coefficients
00330     for( count = 0; count < num_b; ++count )
00331     {
00332       double next_coef;
00333       sub_count = count % 6;  //  there are 6 entries per line
00334       if( sub_count == 0 )
00335       {
00336         getline( inFile, linebuff); //Read in new line
00337       }
00338       strbuff = linebuff.substr( sub_count*11, 11 );
00339       read_d( &strbuff, &next_coef );
00340       b_coef[ count ] = next_coef;
00341     }
00342     sub_count = ( sub_count + 1 ) % 6;  // we just read the last b
00343     // read the a-coefficients
00344     for(  ; count < num_data; count += 2 )
00345     {
00346       double real_part;
00347       double imag_part;
00348 
00349       if( sub_count == 0 )
00350       {
00351         getline( inFile, linebuff); //Read in new line
00352       }
00353       // the real and imaginary parts may be on different lines
00354       if( sub_count < 5 )
00355       {
00356     // same line
00357         strbuff = linebuff.substr( sub_count*11, 22 );
00358         read_dd( &strbuff, &real_part, &imag_part );
00359         a_coef[ count - num_b ] = complex< double >( real_part, imag_part );
00360       }
00361       else
00362       {
00363     // different lines
00364         strbuff = linebuff.substr( sub_count*11, 11 );
00365         read_d( &strbuff, &real_part );
00366         getline( inFile, linebuff); //Read in new line
00367         strbuff = linebuff.substr( 0, 11 );
00368         read_d( &strbuff, &imag_part );
00369         a_coef[ count - num_b ] = complex< double >( real_part, imag_part );
00370       }
00371 
00372       sub_count = ( sub_count + 2 ) % 6;  //  there are 6 entries per line
00373     }
00374   }
00375 }
00376 // ------------ one_d_charge::read_table -----------------------
00377 void one_d_charge::read_table( mf6_file& inFile, double e_in, int num_pairs )
00378 {
00379   string linebuff;  // the input line
00380   string strbuff;   // a substring
00381 
00382   dd_link link;  // for creating new links
00383 
00384   // read the (mu, probability) pairs---the deviation from Rutherford
00385   for( int count = 0; count < num_pairs; ++count )
00386   {
00387     int sub_count = count % 3;
00388     if( sub_count == 0 )
00389     {
00390       getline( inFile, linebuff); //Read in new line
00391     }
00392     strbuff = linebuff.substr( sub_count*22, 22 );
00393     read_dd( &strbuff, &link.x, &link.y );
00394     deviation.insert( deviation.end( ), link );
00395   }
00396 }
00397 // ----------- one_d_charge::expand -----------
00398 void one_d_charge::expand( )
00399 // Expand a charge model into a linked list.
00400 {
00401   dd_link link;  // for creating new links
00402 
00403   // initialize the list with maximum and minimum mu
00404   double mu = Global.Value( "max_coulomb_mu" );
00405   if( same_ )
00406   {
00407     // initialize the list with 0 and maximum mu
00408     link = dd_link( 0.0, f( 0.0 ) );
00409   }
00410   else
00411   {
00412     // initialize the list with maximum and minimum mu
00413     link = dd_link( -1.0, f( -1.0 ) );
00414   }
00415   insert( end( ), link );
00416 
00417   link = dd_link( mu, f( mu ) );
00418   insert( end( ), link) ;
00419 
00420   thicken( );
00421 
00422   // For identical particles ENDF models only the forward emission;
00423   // the other particle gives the backward emission.  In order to
00424   // avoid the jump discontinuity, we use the average for both particles.
00425   one_d_charge::iterator this_link = begin( );
00426   for( ++this_link; this_link != end( ); ++this_link )
00427   {
00428     link.x = -this_link->x;
00429     link.y = this_link->y;
00430     push_front( link );
00431   }
00432 }
00433 // ----------- one_d_charge::expand_12 -----------
00434 void one_d_charge::expand_12( double scale )
00435 // Expand a charge model into a linked list.
00436 {
00437   // the angular distribution is coulomb + xs*deviation
00438   deviation *= scale;
00439   *this += deviation;
00440 
00441   // check for positivity
00442   for( one_d_charge::iterator this_link = begin( );
00443        this_link != end( ); ++this_link )
00444   {
00445     if( this_link->y < 0.0 )
00446     {
00447       this_link->y = 0.0;
00448     }
00449   }
00450 }
00451 
00452 // ********* for class two_d_charge *************
00453 // ----------- two_d_charge::master -----------------
00454 // Reads in all the data and expands into linearly interpolable pointwise data
00455 void two_d_charge::master( mf6_file& inFile )
00456 {
00457   bool done = false;
00458   xs_read = false;  // cross section data not yet read
00459 
00460   double spin;
00461   int NR;  // number of interpolation regions
00462   int NE;  // number of incident energies
00463 
00464   inFile.read_spin( &spin, &NR, &NE );
00465   inFile.get_regions( NR, NBT, INT );
00466 
00467   // do we need to implement interpolation?
00468   bool bad_interp = false;
00469   int j;
00470   for( j = 0; j < INT.size( ); ++j )
00471   {
00472     if( ( INT[ j ] != 2 ) && ( INT[ j ] != 6 ) )
00473     {
00474       bad_interp = true;
00475       break;
00476     }
00477   }
00478   if ( bad_interp )
00479   {
00480     Warning( "two_d_charge::master",
00481         pastenum( "strange interpolation: ", INT[ j ] ) );
00482   }
00483 
00484   // handle each incident neutron energy
00485   for ( int iNE = 0; iNE < NE; iNE++ )
00486   {
00487     if( !done )
00488     {
00489       one_E_in( inFile, spin, &done );
00490     }
00491   }
00492 
00493   two_d_charge::iterator this_link = begin( );
00494   for(  ; this_link != end( ); ++this_link )
00495   {
00496     // convert the formulas to tables
00497     this_link->expand( );
00498     if( this_link->LTP_model > 2 )
00499     {
00500       // add xs*deviation
00501       double scale = xs.evaluate( this_link->E_in( ) );
00502       this_link->expand_12( scale );
00503     }
00504   }
00505   // get the reaction cross section
00506   get_xs( );
00507 
00508   if( ENDL.write_file ) 
00509   {
00510     ENDL.set_yo( ENDL.incident_particle );
00511     write_endl( 1 ); //write the ENDL file
00512     // do we have charged-particle scattering for a light target?
00513     if ( ( ENDL.T == 2 ) && ( ENDL.ZA <= 2004 ) )
00514     {
00515       ENDL.outgoing_particle = 10 + ENDL.za_to_yo( ENDL.ZA );
00516       mirror( );
00517       write_endl( 1 );
00518     }
00519   }
00520 }
00521 
00522 // ----------- two_d_charge::one_E_in -----------------
00523 // Reads in data for one incident energy
00524 // Sets *done = true if the incident energy is at or above the maximum
00525 void two_d_charge::one_E_in( mf6_file& inFile, double spin,
00526   bool *done )
00527 {
00528   // make a new link for this incident energy
00529   one_d_charge new_link;
00530 
00531   // read the model and the incident energy
00532   double e_in;
00533   int LTP;  // the model
00534   int NL;   // depends on the model
00535   inFile.read_model( &e_in, &LTP, &NL );
00536 
00537   if( !(*done) )
00538   {
00539     // energy in MeV
00540     e_in *= ENDL.eV2MeV;
00541 
00542     // insert the 1-d list as a link
00543     insert( end( ), new_link );
00544     two_d_charge::iterator this_link = end( );
00545     --this_link;
00546 
00547     this_link->LTP_model = LTP;
00548     this_link->initiate( inFile, e_in, spin );
00549 
00550     if( LTP == 1 )
00551     {
00552       // for this model NL is the Legendre order
00553       this_link->order = NL;
00554       // read the Legendre coefficients
00555       this_link->read_data( inFile );
00556     }
00557     else if( LTP > 2 )
00558     {
00559       this_link->read_table( inFile, e_in, NL );
00560       if( LTP == 14 )
00561       {
00562         Unimplemented( "two_d_charge::one_E_in",
00563          "implement lin-log interpolation" );
00564       }
00565       // did we already read the cross section file?
00566       if( !xs_read )
00567       {
00568     read_xs( );
00569       }
00570     }
00571     else
00572     {
00573       Unimplemented( "two_d_charge::one_E_in",
00574          pastenum( "implement model: ", LTP ) );
00575     }
00576   }
00577   else
00578   {
00579     // just read the rest of the data
00580   }
00581 
00582   double dE = ENDL_EPSILON( ENDL.Max_E_in );
00583   if( e_in > ENDL.Max_E_in - dE )
00584   {
00585     *done = true;
00586   }
00587 }
00588 // ----------- two_d_charge::read_xs -----------------
00589 // Reads in the cross section data
00590 void two_d_charge::read_xs( )
00591 {
00592   mf3_file xs_File;   //Instantiate the input file
00593   xs_File.open(3, ENDL.T);
00594 
00595   // get the target ZA and mass
00596   int ZA;
00597   double AWR;
00598   xs_File.read_line1(&ZA, &AWR);
00599 
00600   //read second line
00601   double QM, QI;
00602   int LR, NR, NP;
00603   xs_File.read_line2(
00604           &QM,             // mass difference Q value
00605           &QI,             // reaction Q value for the lowest state
00606           &LR,             // complex or "breakup" flag
00607           &NR,             // number of INT regions 
00608           &NP);            // number of data points
00609   xs_File.get_regions( NR, xs.NBT, xs.INT );
00610   xs.read_data( NP, xs_File, ENDL.Max_E_in );
00611   xs_File.close();
00612 
00613   xs_read = true;
00614 }
00615 // ----------- two_d_charge::get_xs -----------------
00616 // Calculates the cross section
00617 void two_d_charge::get_xs( )
00618 {
00619   one_d_table new_xs;  // Instantiate the cross section object
00620   dd_link link;  // for creating new links
00621 
00622   two_d_charge::iterator this_link = begin( );
00623   for(  ; this_link != end( ); ++this_link )
00624   {
00625     link.x = this_link->E_in( );
00626     link.y = this_link->get_norm( );
00627     if( link.y > 0.0 )
00628     {
00629       this_link->renorm( link.y );
00630     }
00631     new_xs.insert( new_xs.end( ), link );
00632   }
00633 
00634   // print the cross sections
00635   ENDL.set_yo( 0 );
00636   new_xs.write_endl( 0 );
00637 }

Generated on Thu Sep 7 10:30:02 2006 for fete -- From ENDFB6 To ENDL by doxygen 1.3.4