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

kalbach.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: 1873 $
00029  * $Date: 2006-05-17 08:45:49 -0700 (Wed, 17 May 2006) $
00030  * $Author: dbrown $
00031  * $Id: kalbach.cpp 1873 2006-05-17 15:45:49Z dbrown $
00032  * 
00033  * ******** fete: From ENDF To ENDL *********
00034  */
00035 
00036 // implementation of the class three_d_kalbach
00037 
00038 #include <cmath>
00039 
00040 #include "kalbach.hpp"
00041 #include "convert.hpp"
00042 #include "messaging.hpp"
00043 
00044 extern ENDLClass ENDL;
00045 
00046 // ********* for class Kalbach_data *************
00047 // ------------ Kalbach_data::read_data ----------
00048 void Kalbach_data::read_data( mf6_file& inFile, int num_E_out )
00049 // read in the data for one incident energy
00050 {
00051   string sbuff;
00052   Kalbach_triple triple;
00053  
00054   for ( int count = 0; count < num_E_out; ++count )
00055   {
00056     inFile >> sbuff;
00057     // for the particle energy in MeV
00058     triple.E_out = ENDL.eV2MeV*stod( sbuff );
00059     inFile >> sbuff;
00060     triple.Prob_E_out = stod( sbuff );
00061     inFile >> sbuff;
00062     triple.compound_frac = stod( sbuff );
00063     insert(end(), triple);
00064   }
00065   inFile.ignore(66,'\n');  // Must clear rest of line after using >>
00066 }
00067 // ----------- Kalbach_data::find_data -----------------
00068 Kalbach_data::iterator Kalbach_data::find_data(double E_out)
00069 // Find E_out in the data list and interpolate the data.
00070 // *** This routine is writen for histograms; for linear ***
00071 // *** interpolation, write new code. ***
00072 {
00073   const double eps = 1.0e-6;  // a tolerance on real numbers
00074   Kalbach_data::iterator data_ptr = end();
00075   for(--data_ptr; ; --data_ptr)
00076   {
00077     if(( E_out >= data_ptr->E_out ) ||
00078       ( (data_ptr == begin()) && (E_out >= (1 - eps)*data_ptr->E_out)) )
00079     {
00080       return data_ptr;
00081     }
00082 
00083     // Is E_out below the first one?
00084     if(data_ptr == begin())
00085     {
00086       Warning("Kalbach_data::find_data",pastenum("E_out too small in Kalbach: ",E_out));
00087       return data_ptr;
00088     }
00089   }
00090 }
00091 
00092 // ********* for class one_d_Kalbach ************
00093 // ----------- one_d_Kalbach::one_E_row ----------------
00094 void one_d_Kalbach::one_E_row( double min_E, double max_E )
00095 // for a given mu make new links for secondary energies
00096 // min_E <= E <= max_E
00097 {
00098   // a data check
00099   if( min_E > max_E )
00100   {
00101     SevereError( "one_d_Kalbach::one_E_row", "energies out of order" );
00102   }
00103 
00104   // check for forward in lab frame, backward in c-m frame
00105   double E_lab_trans = mu( ) * mu( ) * map_->E_transl;  // min c-m energy
00106 
00107   if( ( mu( ) > 0.0 ) && ( max_E < E_lab_trans ) )
00108   {
00109     // E_cm is only decreasing
00110     cm_E_decr( min_E, max_E );
00111   }
00112   else if( ( mu( ) > 0.0 ) && ( min_E < E_lab_trans ) &&
00113        ( max_E > E_lab_trans ) )
00114   {
00115     // E_cm switches direction
00116     cm_E_decr( min_E, E_lab_trans );
00117     cm_E_incr( E_lab_trans, max_E );
00118   }
00119   else
00120   {
00121     // E_cm is increasing
00122     cm_E_incr( min_E, max_E );
00123   }
00124 }
00125 // ----------- one_d_Kalbach::cm_E_incr -----------------
00126 // For this mu and this range of lab energies the cm energies increase.
00127 // Spefically, we are in the region E_lab > map_->E_transl - E_cm
00128 void one_d_Kalbach::cm_E_incr( double min_E, double max_E )
00129 {
00130   Pair lab_Emu; // (energy, cosine) in lab coordinates
00131   Pair cm_Emu;  // (energy, cosine) in center-of-mass coordinates
00132   // for new links
00133   double E_out;
00134   double Prob;
00135   dd_link e_out_link;
00136   double dE;    // to widen the histogram
00137 
00138   lab_Emu.mu = mu( );
00139 
00140   // For given E_cm and mu_lab the algebra gives 2 values of E_lab.
00141   // We want the larger one if mu_lab > 0:
00142   int sigma = ( mu( ) > 0.0 ) ? 1 : -1;
00143 
00144   // insert the first point
00145   dE = ENDL_EPSILON( min_E );
00146   lab_Emu.E =  min_E + dE;
00147   // in center-of-mass coordinates
00148   cm_Emu = map_->lab_to_cm( lab_Emu );
00149   // E_f0_r_ptr points to the Kalbach data
00150   E_f0_r_ptr = data_ptr__->find_data( max(cm_Emu.E,0.0) );
00151   // insert this data
00152   Prob = f( lab_Emu.E );
00153   e_out_link = dd_link( lab_Emu.E, Prob );
00154   insert( end( ), e_out_link );
00155 
00156   // identify the largest cm outgoing energy
00157   dE = ENDL_EPSILON( max_E );
00158   lab_Emu.E = max_E - dE;
00159   cm_Emu = map_->lab_to_cm( lab_Emu );
00160   double E_cm_max = cm_Emu.E;
00161 
00162   // The Kalbach data is given as a histogram in cm coordinates.
00163   // Loop until we get to E_cm_max
00164   Kalbach_data::iterator next_data;
00165   for( next_data = E_f0_r_ptr; ; )
00166   {
00167     ++next_data;  // the next Kalbach triplet
00168     one_d_Kalbach::iterator this_link;
00169     one_d_Kalbach::iterator next_link;
00170 
00171     // Are we finished?
00172     if( ( next_data == data_ptr__->end( ) ) ||
00173         ( next_data->E_out > E_cm_max ) )
00174     {
00175       // insert a last link to finish the list
00176       this_link = end( );
00177       --this_link;
00178 
00179       // insert the new last link
00180       Prob = f( max_E );
00181       e_out_link = dd_link( max_E, Prob );
00182       insert( end( ), e_out_link );
00183 
00184       next_link = end( );
00185       --next_link;
00186       thicken( this_link, next_link );
00187       break;
00188     }
00189 
00190     // point to the current last link
00191     this_link = end( );
00192     --this_link;
00193 
00194     // get the end of the histogram
00195     E_out = map_->get_lab_E( next_data->E_out, mu(), sigma );
00196     dE = ENDL_EPSILON( E_out );
00197     // insert the new last link
00198     Prob = f( E_out - dE );
00199     e_out_link = dd_link( E_out - dE, Prob );
00200     insert( end( ), e_out_link );
00201 
00202     next_link = end( );
00203     --next_link;
00204     thicken( this_link, next_link );
00205 
00206     // start the next step of the histogram
00207     E_f0_r_ptr = next_data;
00208     Prob = f( E_out + dE );
00209     e_out_link = dd_link( E_out + dE, Prob );
00210     insert( end( ), e_out_link );
00211   }
00212 
00213   // We currently do not do thickening with respect to mu,
00214   // but this is the place to set up the equiprabable bins
00215   // if and when this thickening gets implemented.
00216   //    if( bin_interp )
00217   //    {
00218   //      get_bins( );
00219   //    }
00220 
00221 }
00222 // ----------- one_d_Kalbach::cm_E_decr -----------------
00223 // For this mu and this range of lab energies the cm energies decrease.
00224 // Specifically, we are in the region E_lab < map_->E_transl - E_cm
00225 void one_d_Kalbach::cm_E_decr( double min_E, double max_E )
00226 {
00227   Pair lab_Emu; // (energy, cosine) in lab coordinates
00228   Pair cm_Emu;  // (energy, cosine) in center-of-mass coordinates
00229   // for new links
00230   double E_out;
00231   double Prob;
00232   dd_link e_out_link;
00233   double dE;    // to widen the histogram
00234 
00235   lab_Emu.mu = mu( );
00236 
00237   // For given E_cm and mu_lab the algebra gives 2 values of E_lab.
00238   // We want the smaller one
00239   int sigma = -1;
00240 
00241   // A first cut at a shift away from the singularities
00242   dE = ENDL_EPSILON( min_E );
00243   // We may need something smaller
00244   double E_diff = max_E - min_E;
00245   if( 4*dE > E_diff ){
00246     dE = E_diff/4.0;
00247   }
00248 
00249   // insert the first point
00250   lab_Emu.E =  min_E + dE;
00251   // in center-of-mass coordinates
00252   cm_Emu = map_->lab_to_cm( lab_Emu );
00253   // E_f0_r_ptr points to the Kalbach data
00254   E_f0_r_ptr = data_ptr__->find_data( max(0.0,cm_Emu.E) );
00255   // insert this data
00256   Prob = f( lab_Emu.E );
00257   e_out_link = dd_link( lab_Emu.E, Prob );
00258   insert( end( ), e_out_link );
00259 
00260   // identify the smallest cm outgoing energy
00261   dE = ENDL_EPSILON( max_E );
00262   lab_Emu.E = max_E - dE;
00263   cm_Emu = map_->lab_to_cm( lab_Emu );
00264   double E_cm_min = cm_Emu.E;
00265 
00266   // The Kalbach data is given as a histogram in cm coordinates.
00267   // Loop until we get to E_cm_min
00268   Kalbach_data::iterator this_data;
00269   for( this_data = E_f0_r_ptr; ; )
00270   {
00271     one_d_Kalbach::iterator this_link;
00272     one_d_Kalbach::iterator next_link;
00273 
00274     // Are we finished?
00275     if( ( this_data == data_ptr__->begin( ) ) ||
00276         ( this_data->E_out < E_cm_min ) )
00277     {
00278       // insert a last link to finish the list
00279       this_link = end( );
00280       --this_link;
00281 
00282       // insert the new last link
00283       Prob = f( max_E );
00284       e_out_link = dd_link( max_E, Prob );
00285       insert( end( ), e_out_link );
00286 
00287       next_link = end( );
00288       --next_link;
00289       thicken( this_link, next_link );
00290       break;
00291     }
00292 
00293     // point to the current last link
00294     this_link = end( );
00295     --this_link;
00296 
00297     // get the start of the histogram
00298     E_out = map_->get_lab_E( this_data->E_out, mu(), sigma );
00299     dE = ENDL_EPSILON( E_out );
00300     // insert the new last link
00301     Prob = f( E_out - dE );
00302     e_out_link = dd_link( E_out - dE, Prob );
00303     insert( end( ), e_out_link );
00304 
00305     next_link = end( );
00306     --next_link;
00307     thicken( this_link, next_link );
00308 
00309     // start the next step of the histogram (the previous cm energy)
00310     --this_data;  // the next Kalbach triplet
00311     E_f0_r_ptr = this_data;
00312     Prob = f( E_out + dE );
00313     e_out_link = dd_link( E_out + dE, Prob );
00314     insert( end( ), e_out_link );
00315   }
00316 }
00317 // ----------- one_d_Kalbach::copy_params -----------------
00318 // Copy the parameters
00319 void one_d_Kalbach::copy_params( Param *params )
00320 {
00321   params__ = params;
00322 }
00323 
00324 // ----------- one_d_Kalbach::copy_data -----------------
00325 // Copy the parameters
00326 void one_d_Kalbach::copy_data( list< Kalbach_data >::iterator data_ptr )
00327 {
00328   data_ptr__ = data_ptr;
00329 }
00330 
00331 // ----------- one_d_Kalbach::get_a ---------------
00332 // Evaluate the Kalbach a function for the slope value
00333 double one_d_Kalbach::get_a( double E_out )
00334 {
00335   return Get_a( E_out, *params__ );
00336 }
00337 
00338 // ----------- one_d_Kalbach::f -----------------
00339 double one_d_Kalbach::f( double E )
00340 // Evaluate the Kalbach probability function for
00341 // cosine mu and energy E_out (laboratory coordinates).
00342 {
00343   Pair lab_Emu; // (energy, cosine) in lab coordinates
00344   Pair cm_Emu;  // (energy, cosine) in center-of-mass coordinates
00345   double Prob;  // return value
00346 
00347   lab_Emu.mu = mu( );
00348   lab_Emu.E = E;
00349 
00350   // in center-of-mass coordinates
00351   cm_Emu = map_->lab_to_cm( lab_Emu );
00352 
00353   if( E_f0_r_ptr == NULL )
00354   {
00355     E_f0_r_ptr = data_ptr__->find_data( max(0.0,cm_Emu.E) );
00356   }
00357 
00358   double slope_a = get_a( cm_Emu.E );
00359   // DAB 2/3/05 -- factor of 0.5 missing from older version of ENDF manual
00360   Prob = 0.5*(E_f0_r_ptr->Prob_E_out*( slope_a/sinh( slope_a ) )*
00361        ( cosh( slope_a*cm_Emu.mu ) +
00362      E_f0_r_ptr->compound_frac*sinh( slope_a*cm_Emu.mu ) ));  
00363 
00364   // in laboratory coordinates
00365   Prob *= map_->J_cm_to_lab( cm_Emu.E, lab_Emu.E );
00366   return Prob;
00367 }
00368 
00369 // ********* for class two_d_Kalbach *************
00370 // ----------- two_d_Kalbach::new_one_d ------
00371 two_d_Kalbach::iterator
00372   two_d_Kalbach::new_one_d( two_d_Kalbach::iterator where )
00373 {
00374   one_d_Kalbach new_link;
00375   insert( where, new_link );
00376 
00377   // point to the new link
00378   two_d_Kalbach::iterator link_ptr = where;
00379   --link_ptr;
00380 
00381   // set up the parameters
00382   link_ptr->copy_params( params_ );
00383   link_ptr->copy_map( &map );
00384   link_ptr->copy_data( data_ptr_ );
00385 
00386   return link_ptr;
00387 }
00388 // ----------- two_d_Kalbach::set_params ------
00389 void two_d_Kalbach::set_params(
00390   list< Kalbach_data >::iterator data_ptr, Param *params )
00391 {
00392   data_ptr_ = data_ptr;
00393   params_ = params;
00394 }
00395 
00396 // ********* for class three_d_Kalbach *************
00397 // ----------- three_d_Kalbach::three_d_Kalbach ------
00398 three_d_Kalbach::three_d_Kalbach( int targ_ZA, int proj_ZA,
00399   int eject_ZA, int this_mult )
00400 // constructor with given target, projectile, and ejectum.
00401 {
00402   // this model makes no sense for natural targets
00403   if( targ_ZA % 1000 == 0 )
00404   {
00405     SevereError( "three_d_Kalbach::three_d_Kalbach",
00406          " the Kalbach model makes no sense for natural targets" );
00407   }
00408 
00409   this_mult_ = this_mult;
00410 
00411   A.set(targ_ZA);
00412   a.set(proj_ZA);
00413   b.set(eject_ZA);
00414   C.set(targ_ZA + proj_ZA);
00415   B.set(C.ZA - eject_ZA);
00416 
00417   // get the separation energies
00418   S_a = get_S(A, a);
00419 
00420   // We need 5 parameters for the Kalbach function
00421   Params.get_space( 5 );
00422   Params[0] = 0.0;  // e_a depends on the incident energy
00423   Params[1] = ( B.AWR + b.AWR ) / B.AWR;  // needed to calculate epsilon_b
00424   Params[2] = get_S( B, b );  //  S_b;
00425 
00426   // Params[3] is the Kalbach M_a factor
00427   switch(a.ZA)
00428   {
00429   case 1:  // neutron
00430     Params[3] = 1.0;
00431     break;
00432   case 1001:  // proton
00433     Params[3] = 1.0;
00434     break;
00435   case 1002:  // deuteron
00436     Params[3] = 1.0;
00437     break;
00438   case 1003:  // tritron
00439     Params[3] = 0.0;
00440     break;
00441   case 2003:  // helium 3
00442     Params[3] = 0.0;
00443     break;
00444   case 2004:  // alpha
00445     Params[3] = 0.0;
00446     break;
00447   default:
00448     Unimplemented("three_d_Kalbach::three_d_Kalbach",
00449         pastenum("M for particle ",a.ZA)+" not implemented");
00450   }
00451 
00452   // Params[4] is the Kalbach m_b factor
00453   switch(b.ZA)
00454   {
00455   case 1:  // neutron
00456     Params[4] = 0.5;
00457     break;
00458   case 1001:  // proton
00459     Params[4] = 1.0;
00460     break;
00461   case 1002:  // deuteron
00462     Params[4] = 1.0;
00463     break;
00464   case 1003:  // tritron
00465     Params[4] = 1.0;
00466     break;
00467   case 2003:  // helium 3
00468     Params[4] = 1.0;
00469     break;
00470   case 2004:  // alpha
00471     Params[4] = 2.0;
00472     break;
00473   default:
00474     Unimplemented("three_d_Kalbach::three_d_Kalbach",
00475         pastenum("M for particle ",b.ZA)+" not implemented");
00476   }
00477 
00478 }
00479 // ----------- three_d_Kalbach::expand_data -----------------
00480 void three_d_Kalbach::expand_data(int num_E_in, mf6_file& inFile)
00481 // read in all of the ENDF/B-VI data and expand it to a list
00482 // num_E_in is the number of data sets
00483 {
00484   // read in the data
00485   read_data( inFile, num_E_in );
00486   if(!ENDL.write_file )
00487   {
00488     return;
00489   }
00490 
00491   if( Global.Value( "Kalbach_i10" ) > 0 )
00492   {
00493     // Check the endep calculation of <E'>
00494     check_i10( );
00495   }
00496 
00497   // make a link for each incident energy
00498   for( Data_ptr = triples.begin( ); Data_ptr != triples.end( );
00499     ++Data_ptr )
00500   {
00501     one_e_in();
00502   }
00503 }
00504 // ----------- three_d_Kalbach::one_e_in -----------------
00505 void three_d_Kalbach::one_e_in( )
00506 // handle the data for one incident energy
00507 {
00508   // this E_1 corresponds to a new link in the 3d list
00509   two_d_Kalbach e_in_link;
00510   e_in_link.E_in( ) = Data_ptr->E_in;
00511   insert( end( ), e_in_link );
00512 
00513   // make a pointer to the new 3-d link
00514   three_d_Kalbach::iterator this_link = end();
00515   --this_link;
00516 
00517   // copy the parameters
00518   this_link->set_params( Data_ptr, &Params );
00519 
00520   // Expand this data into a 3-d link
00521   expand_E_in( this_link );
00522 
00523   // add a link to the cosine list
00524   one_d_table new_cos;
00525   cosines.insert(cosines.end(), new_cos);
00526 
00527   // expand this cosine link
00528   two_d_list<one_d_table>::iterator cos_ptr = cosines.end();
00529   --cos_ptr;
00530   this_link->make_cos_link( Data_ptr->E_in, cos_ptr );
00531 
00532 }
00533 
00534 // ----------- three_d_Kalbach::expand_E_in -----------------
00535 void three_d_Kalbach::expand_E_in( three_d_Kalbach::iterator e_in_link )
00536 // expand the double differential data for 1 incident energy
00537 {
00538   // center-of-mass incident energy plus separation energy of particle.
00539   // This e_a is used later by the get_a routine
00540   Params[0] = e_in_link->E_in()*A.AWR/(A.AWR + a.AWR) + S_a;
00541 
00542   // For ENDL we need to go through cosines first and then E_out.
00543   // But the ENDF/B6 data has E_out first.
00544 
00545   // for the lowest particle energy in MeV
00546   Kalbach_data::iterator this_data = Data_ptr->begin( );
00547   double first_E_out = this_data->E_out;
00548 
00549   // for the highest particle energy in MeV
00550   this_data = Data_ptr->end( );
00551   --this_data;
00552   double last_E_out = this_data->E_out;
00553 
00554   // to map between lab and center-of-mass coordinates
00555   e_in_link->map.set_map( &A, &a, &b, e_in_link->E_in() );
00556 
00557   // identify the geometry of the (E', mu) region in lab coordinates.
00558   e_in_link->get_geom( first_E_out, last_E_out );
00559   if( e_in_link->geom_.fastest == Both_ways )
00560   {
00561     // The fastest outgoing particles can have mu < 0
00562     if( e_in_link->geom_.slowest == Straight )
00563     {
00564       // there are outgoing particles with CM speed 0
00565       e_in_link->full_range( last_E_out );
00566     }
00567     else if( e_in_link->geom_.slowest == Forward )
00568     {
00569       // some outgoing particles are only forward in the lab frame
00570       e_in_link->both_ways_hole( first_E_out, last_E_out );
00571     }
00572     else if( e_in_link->geom_.slowest == Transition )
00573     {
00574       // the slowest outgoing particles move at the CM translational speed
00575       e_in_link->both_ways_trans( first_E_out, last_E_out );
00576     }
00577     else
00578     {
00579       // all outgoing particles move faster than the CM translational speed
00580       e_in_link->all_both_ways( first_E_out, last_E_out );
00581     }
00582   }
00583   else if( e_in_link->geom_.fastest == Transition )
00584   {
00585     // the fastest outgoing particles move at the CM translational speed
00586     if( e_in_link->geom_.slowest == Straight )
00587     {
00588       // there are outgoing particles with CM speed 0
00589       e_in_link->full_trans( last_E_out );
00590     }
00591     else
00592     {
00593       // the slowest outgoing particles have postive CM speed
00594       e_in_link->trans_hole( first_E_out, last_E_out );
00595     }
00596   }
00597   else
00598   {
00599     // the fastest outgoing particles all move forward in the lab frame
00600     if( e_in_link->geom_.slowest == Straight )
00601     {
00602       // there are outgoing particles with CM speed 0
00603       e_in_link->full_forward( 0.0, last_E_out );
00604     }
00605     else
00606     {
00607       // the slowest outgoing particles have postive CM speed
00608       e_in_link->forward_hole( first_E_out, last_E_out );
00609     }
00610   }
00611 }
00612 
00613 // ----------- three_d_Kalbach::read_data -----------------
00614 void three_d_Kalbach::read_data( mf6_file& inFile, int num_E_in )
00615 // read in all of the ENDF/B-VI data
00616 // num_E_in is the number of data sets
00617 {
00618   Kalbach_data one_set;  // Kalbach data for 1 incident energy
00619   list<Kalbach_data>::iterator data_ptr;  // point to the data
00620 
00621   double ZERO;  // zero
00622   double E_1;   // incident neutron energy
00623   int ND;   // number of discrete energies
00624   int NA;   // number of triples            //FIXME!!!
00625   int NW;   // number of words in the record (3*NA)     //FIXME!!!
00626   int NEP;  // number of secondary energies
00627 
00628   string linebuff;
00629   for(int count = 0; count < num_E_in; ++count)
00630   {
00631     getline( inFile, linebuff ); //read next line
00632     read_ddiiii(&linebuff, 
00633           &ZERO, &E_1,   // incident energy
00634           &ND,           // number of discrete energies
00635           &NA,           // number of angular parameters (1)
00636           &NW,           // number of table entries
00637           &NEP);         // number of E_out values
00638     if(ND > 0)
00639     {
00640       SevereError("three_d_Kalbach::read_data","Implement discrete energies for Kalbach");
00641     }
00642 
00643     // build the list for this incident energy (in MeV)
00644     one_set.E_in = E_1*ENDL.eV2MeV;
00645     triples.insert(triples.end(), one_set);
00646 
00647     // point to the new link
00648     data_ptr = triples.end();
00649     --data_ptr;
00650     // fill this list
00651     data_ptr->read_data( inFile, NEP );
00652   }
00653 }
00654 /*
00655 // --------- three_d_Kalbach::fill_in_data -----------------
00656 void three_d_Kalbach::fill_in_data();
00657 // If the emission jumps from all forward to all directions,
00658 // insert a smoother transition.
00659 {
00660   // The translational energy from the initial collisiton is lambda*E_in
00661   double lambda = (A.AWR)*(b.AWR)/((A.AWR + a.AWR)*(A.AWR + a.AWR));
00662 
00663   // look for a jump from all forward
00664   for(list<Kalbach_data>::iterator data_ptr = triples.begin();
00665     data_ptr != triples.end(); ++data_ptr)
00666   {
00667     Kalbach_data::iterator last_set = data_ptr.end();
00668     --last_set;
00669 
00670   }
00671 
00672 }
00673 */
00674 // ----------- three_d_Kalbach::get_S -----------------
00675 double three_d_Kalbach::get_S(Nuclei& targ_A, Nuclei& proj_a)
00676 // evaluate the Kalbach S function for the separation energy
00677 {
00678   double S;
00679 
00680   double I_a;
00681   switch(proj_a.ZA)
00682   {
00683   case 1:  // neutron
00684     I_a = 0.0;
00685     break;
00686   case 1001:  // proton
00687     I_a = 0.0;
00688     break;
00689   case 1002:  // deuteron
00690     I_a = 2.22;
00691     break;
00692   case 1003:  // tritron
00693     I_a = 8.48;
00694     break;
00695   case 2003:  // helium 3
00696     I_a = 7.72;
00697     break;
00698   case 2004:  // alpha
00699     I_a = 28.3;
00700     break;
00701   default:
00702     Unimplemented("three_d_Kalbach::get_S",
00703         pastenum("particle ",proj_a.ZA)+" not implemented in get_S");
00704   }
00705 
00706   const double coef_1 = 15.68;
00707   const double coef_2 = 28.07;
00708   const double coef_3 = 18.56;
00709   const double coef_4 = 33.22;
00710   const double coef_5 = 0.717;
00711   const double coef_6 = 1.211;
00712 
00713   S = coef_1*(C.A - targ_A.A) -
00714       coef_2*((C.N - C.Z)*(C.N - C.Z)/C.A -
00715               (targ_A.N - targ_A.Z)*(targ_A.N - targ_A.Z)/
00716                  targ_A.A) -
00717       coef_3*(exp(2.0/3.0*log(1.0*C.A)) -
00718               exp(2.0/3.0*log(1.0*targ_A.A))) +
00719       coef_4*((C.N - C.Z)*(C.N - C.Z)/
00720                  exp(4.0/3.0*log(1.0*C.A)) -
00721               (targ_A.N - targ_A.Z)*(targ_A.N - targ_A.Z)/
00722                  exp(4.0/3.0*log(1.0*targ_A.A))) -
00723       coef_5*(C.Z*C.Z/exp(1.0/3.0*log(1.0*C.A)) -
00724               targ_A.Z*targ_A.Z/exp(1.0/3.0*log(1.0*targ_A.A))) +
00725       coef_6*(C.Z*C.Z/C.A -
00726               targ_A.Z*A.Z/targ_A.A) - I_a;
00727 
00728   return S;
00729 }
00730 // ----------- three_d_Kalbach::renorm -----------------
00731 void three_d_Kalbach::renorm( ){
00732   for (three_d_list< two_d_Kalbach >::iterator it=begin();it!=end();++it) it->renorm();
00733 }
00734 
00735 // ----------- three_d_Kalbach::check_i10 -----------------
00736 void three_d_Kalbach::check_i10( )
00737 // Write an i=10 file of average energy for outgoing particle.
00738   // What we want is
00739   // (1/2) * \int_\Omega E_b P( E_b, \mu ) \, dE_b \, d\mu
00740   // where everything is in laboratory coordinates.  But we are
00741   // given P is center-of-mass coordinates.  It is better to do
00742   // the integration in center-of-mass coordinates, for 2 reasons.
00743   // First, this avoids having to compute the Jacobian of the mapping
00744   // to laboratory coordinates.  Second, in center-of-mass coordinates
00745   // the domain of integration is a rectangle, -1 <= \mu <= 1, and
00746   // 0 <= E_b <= E_{b \max}.  From equations (6.4) and (6.5) of the
00747   // ENDF/B-VI manual, the resulting integral is
00748   // (1/2) * \int_0^{E_b\max} f_0 \int_{-1}^1 [E_b + \gamma^2 E_a +
00749   //      2 \gamma \mu \sqrt{E_a E_b}] [a / sinh(a) ]
00750   //      [cosh(a \mu) + r sinh(a \mu)] \, d\mu \, dE_b.
00751   // Here, \gamma = sqrt{a.AWR * b.AWR}/(A.AWR + a.AWR), the values of
00752   // f_0 and r are given in the mf6 file, and a is from the formula on
00753   // page 6.5 of the  ENDF/B-VI manual.
00754   //
00755   // We can do the integration over \mu by hand, and the result is
00756   // \int_0^{E_b\max} f_0 { E_b + \gamma^2 E_a +
00757   //   [ 2 r \gamma \sqrt{E_a E_b} / sinh(a) ][ cosh(a) - (1/a)sinh(a) ] dE_b.
00758   // Our routine performs this last integration.
00759   // (Note that the integrand looks as if it has a 0/0 singularity
00760   // at a = 0, but this is removable since cosh(a) - (1/a)sinh(a)
00761   // has Taylor series a^2/3 + a^4/30 + ....)
00762 {
00763   one_d_table i10_data;
00764   double gamma = sqrt( a.AWR * b.AWR )/( A.AWR + a.AWR );
00765 
00766   // for adaptive quadrature
00767   quad_list adapt_quad( av_energy );
00768 
00769   // make a link for each incident energy
00770   for( Data_ptr = triples.begin( ); Data_ptr != triples.end();
00771     ++Data_ptr )
00772   {
00773     dd_link new_link;
00774     new_link.x = Data_ptr->E_in;
00775     double E_a_coef = gamma*gamma*new_link.x; // \gamma^2 E_a
00776     double Er_sum = 0;    // used to integrate sqrt(E)*(coth(a) - 1/a)
00777     double E0_sum = 0;    // used to integrate f_0 * { E_b + \gamma^2 E_a }
00778     double norm_sum = 0;  // used to integrate f_0
00779     Kalbach_data::iterator this_data = Data_ptr->begin();
00780     double prev_E_out = this_data->E_out;
00781     double prev_f0 = this_data->Prob_E_out;
00782     double prev_r = this_data->compound_frac;
00783     double integral_coef = 2 * prev_f0 * gamma * sqrt( new_link.x );
00784 
00785     // loop over the (E_out, f0, r) triples
00786     for( ++this_data ; this_data != Data_ptr->end( ); ++this_data )
00787     {
00788       double E_out = this_data->E_out;
00789       // the data is histograms
00790       norm_sum += prev_f0 * ( E_out - prev_E_out );
00791 
00792       E0_sum += prev_f0 * ( E_out - prev_E_out ) *
00793                 ( 0.5*( E_out + prev_E_out ) + E_a_coef );
00794 
00795       Er_sum += integral_coef * prev_r *
00796           adapt_quad.Simp_quad( prev_E_out, E_out, Params, 1.0e-8 );
00797 
00798       prev_f0 = this_data->Prob_E_out;
00799       prev_r = this_data->compound_frac;
00800       prev_E_out = E_out;
00801     }
00802     new_link.y = this_mult_ * ( E0_sum + Er_sum ) / norm_sum;
00803     i10_data.insert( i10_data.end( ), new_link );
00804   }
00805   i10_data.write_endl( 10 );
00806 }
00807 
00808 // *************** Get_a routine *********************
00809 double Get_a( double E_out, Param& params )
00810 // evaluate the Kalbach a function for the slope value
00811 // Note that:
00812   // params[0] = e_a
00813   // params[1] = ( B.AWR + b.AWR ) / B.AWR
00814   // params[2] = S_b
00815   // params[3] = M_a
00816   // params[4] = m_b
00817 {
00818   double slope_a;
00819   double e_a = params[0];
00820   if( e_a == 0.0 )
00821   {
00822     slope_a = 0.0;
00823   }
00824   else
00825   {
00826     // center-of-mass energy plus separation energy for exit particle
00827     double e_b = E_out*params[1] + params[2];
00828 
00829     const double E_t1 = 130.0;  // MeV
00830     double R_1 = ( e_a < E_t1 ) ? e_a : E_t1;
00831     double X_1 = R_1*e_b/e_a;
00832 
00833     const double E_t3 = 41.0;  // MeV
00834     double R_3 = ( e_a < E_t3 ) ? e_a : E_t3;
00835     double X_3 = R_3*e_b/e_a;
00836 
00837     const double C_1 = 0.04;   // per MeV
00838     const double C_2 = 1.8e-6;  // per MeV^3
00839     const double C_3 = 6.7e-7;  // per MeV^4
00840 
00841     slope_a = X_1*(C_1 + C_2*X_1*X_1) +
00842        C_3*params[3]*params[4]*X_3*X_3*X_3;
00843   }
00844   return slope_a;
00845 }
00846 
00847 // **************** average energy integrand *****************
00848 // integrand for the average energy integral
00849 double av_energy( double E_b, Param& params )
00850 {
00851   double slope_a = Get_a( E_b, params );
00852 
00853   double ans;
00854   ans = ( slope_a == 0.0 ) ? 0.0 :
00855     sqrt( slope_a ) * ( cosh( slope_a ) / sinh( slope_a ) - 
00856                        ( 1.0 / slope_a ) );
00857 
00858   return ans;
00859 }

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