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

list_2d.hpp

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: 1882 $
00029  * $Date: 2006-05-24 15:20:29 -0700 (Wed, 24 May 2006) $
00030  * $Author: dbrown $
00031  * $Id: list_2d.hpp 1882 2006-05-24 22:20:29Z dbrown $
00032  * 
00033  * ******** fete: From ENDF To ENDL *********
00034  */
00035 
00036 // header for the class two_d_list
00037 
00038 #ifndef TWO_D_LIST_CLASS
00039 #define TWO_D_LIST_CLASS
00040 
00041 #include "distrib_1d.hpp"
00042 #include "endl_formats.hpp"
00043 #include "record_types.hpp"
00044 #include "ENDF_file.hpp"
00045 #include "global_params.hpp"
00046 #include "endl_precision.hpp"
00047 
00048 using namespace std;
00049 
00050 extern ENDLClass ENDL;
00051 extern GlobalParameterClass Global;
00052 
00053 // ----------- class two_d_list -----------------
00054 //! 2d list template composed of 2d links
00055 template <class one_d> class two_d_list : public list<one_d>
00056 {
00057 
00058 private:
00059   bool debug_on;        //!< for checking the 2-d interpolation
00060   double tol_2d;        //!< tolerance for 2-d interpolation
00061   double cut_off_2d;    //!< below this value ignore tolerance
00062   double dE_tol;        //!< minimum relative difference in incident energy
00063   double _tag;          //!< used with double differential data
00064                         //! switch for type of interpolation to be used
00065   bool _unit_base;      // if true, thicken using unit-based interpolation;
00066                         // if false, thicken based on equiprobable bins
00067 
00068 protected:
00069   dd_list _weight;      //!< weight for this model
00070 
00071 public:
00072   typedef typename list<one_d>::iterator two_d_iterator;
00073 
00074   //! Default constructor
00075   two_d_list()
00076   {
00077     // debug_on = true;  // for checking the 2-d interpolation
00078     debug_on = false;
00079 
00080     // set the parameters
00081     tol_2d = Global.Value("tol_2d");
00082     cut_off_2d = Global.Value("cut_off_2d");
00083     dE_tol = Global.Value("dE_tol");
00084     _tag=0.0;
00085     _unit_base=false;
00086   }
00087 
00088   //! The usual tag for this list
00089   inline double& E_in( )
00090   {
00091     return _tag;
00092   }
00093 
00094   //! The tag for mf14 gamma data
00095   inline double& E_gamma( )
00096   {
00097     return _tag;
00098   }
00099 
00100   //! Contains the interpolation regions
00101   vector<int> NBT;
00102 
00103   //! Contains the interpolation type for each region
00104   vector<int> INT;
00105 
00106   //! Simple printing utility
00107   void print()
00108   {
00109   // loop through the list
00110     for(two_d_iterator list_ptr = two_d_list::begin();
00111       list_ptr != two_d_list::end(); ++list_ptr)
00112     {
00113       cout << "Energy " << ENDL.data( list_ptr->E_in( ) ) << endl;
00114       list_ptr->print();
00115     }
00116     cout << endl;
00117   }
00118 
00119   // ----------- rescale -----------------
00120   //! Scale the energies of the emitted particle from [0, 1]
00121   void rescale()
00122   {
00123   // loop through the list
00124     for(two_d_iterator list_ptr = two_d_list::begin();
00125       list_ptr != two_d_list::end(); ++list_ptr)
00126     {
00127       list_ptr->rescale();
00128     }
00129   }
00130 
00131   // ----------- renorm -----------------
00132   //! Normalize probabilities such that integrated probability is unity
00133   void renorm()
00134   {
00135   // loop through the list
00136     for(two_d_iterator list_ptr = two_d_list::begin();
00137       list_ptr != two_d_list::end(); ++list_ptr)
00138     {
00139       list_ptr->renorm();
00140     }
00141   }
00142 
00143   // ----------- set_up_endl -----------------
00144   //! set up and open an ENDL file
00145   void set_up_endl( int i_number, fstream& endl_file )
00146   {
00147     ENDL.set_I_number(i_number);
00148     string file_name = ENDL.file_name;
00149 
00150     if(two_d_list::empty())
00151     {
00152       Warning("two_d_list::set_up_endl","The file "+file_name+" is empty.");
00153       return;
00154     }
00155 
00156     if ( ENDL.new_file() )
00157     {
00158       endl_file.open(file_name.c_str(),ios::out);
00159       Info("two_d_list::write_endl","Opening ENDL file "+file_name);
00160     }
00161     else
00162     {
00163       endl_file.open(file_name.c_str(),ios::out|ios::app);
00164       Info("two_d_list::write_endl","Appending ENDL file "+file_name);
00165     }
00166 
00167     //Shove the header info in
00168     endl_file << ENDL.header_line_1 << endl;
00169     endl_file << ENDL.header_line_2 << endl;
00170     endl_file.setf(ios::scientific,ios::floatfield);
00171   }
00172 
00173   // ----------- close_file -----------------
00174   //! close an ENDL file
00175   void close_file( fstream& endl_file )
00176   {
00177     //Put on end of file/section line into output file
00178     endl_file << ENDL.eof_line << endl;
00179 
00180     //Finally, we close the output file
00181     Info("two_d_list::close_file","Closing ENDL file "+ENDL.file_name);
00182     endl_file.close();
00183   }
00184 
00185   // ----------- write_endl -----------------
00186   //! Open an ENDL file and write the data.
00187   void write_endl( int i_number )
00188   {
00189     // set up the file
00190     fstream endl_file;
00191     set_up_endl( i_number, endl_file );
00192 
00193     // print the data
00194     out_data( i_number, endl_file, 0.0 );
00195 
00196     // finish writing and close the file
00197     close_file( endl_file );
00198   }
00199 
00200   // ----------- out_data -----------------
00201   // The list holds angular distributions if i_number = 1
00202   // and isotropic energy distributions if i_number = 4.
00203   // L_order is the Legendre order for i4 data
00204   //! Writes out the ENDL data
00205   void out_data( int i_number, fstream& endl_file, double L_order )
00206   {
00207   
00208     // loop over the data and write it out
00209     for(two_d_iterator link = two_d_list::begin(); 
00210         link != two_d_list::end(); ++link)
00211     {
00212       // write out
00213       for(dd_list::iterator prob_link = link->begin();
00214           prob_link != link->end(); ++prob_link)
00215       {
00216         switch(i_number)
00217         {
00218           case 1:  // angular distributions
00219             endl_file << ENDL.data( link->E_in(), 
00220               prob_link->mu(), prob_link->y )
00221               << endl;
00222             break;
00223           case 4:  // energy distributions
00224             endl_file << ENDL.data( link->E_in(), 
00225               prob_link->E_out(), L_order, prob_link->y )
00226              <<endl;
00227             break;
00228           default:
00229             SevereError("two_d_list::out_data",
00230               pastenum("Improper i_number: ",i_number));
00231         }
00232       }
00233     }
00234   }
00235 
00236   // ----------- mirror -----------------
00237   //! Used to reflect all of the lists mu -> -mu
00238   void mirror()
00239   {
00240   // loop through the list
00241     for(two_d_iterator list = two_d_list::begin();
00242       list != two_d_list::end(); ++list)
00243     {
00244       list->mirror();
00245     }   
00246   }
00247 
00248   // ----------- widen_jumps -----------------
00249   //! Widens the jumps for each member 1-d list and for 
00250   //! the main 2d list
00251   void widen_jumps( double cluster_min=0.0 )
00252   {
00253     // loop through the 1d lists that make up the two_d_list first
00254     for( two_d_iterator list1d = two_d_list::begin( );
00255       list1d != two_d_list::end( ); ++list1d )
00256     {
00257       list1d->widen_jumps( 0.0 );
00258     }   
00259 
00260     // Now for the main 2d list
00261 
00262     // scale dE by the smoother used for histograms
00263     double local_eps, jump_width;
00264 
00265     // container to hold links to elements we need to widen
00266     two_d_iterator cluster_start = two_d_list::begin();
00267     two_d_iterator cluster_end = two_d_list::begin();
00268     int cluster_size=1;
00269 
00270     for( two_d_iterator this_link = two_d_list::begin(); 
00271       this_link != two_d_list::end( ); ++this_link )
00272     {
00273       // get the local energy epsilon (smallest dE can see in output file)
00274       local_eps = ENDL_EPSILON( this_link->tag );
00275       jump_width = ENDL_JUMP_WIDTH( this_link->tag );
00276 
00277       // this_link is too close to the cluster_end, so add it to the cluster
00278       if( abs( this_link->tag - cluster_end->tag ) <= local_eps )
00279       {
00280         cluster_end = this_link;
00281         ++cluster_size;
00282       }
00283       // this_link is far from the cluster_end, so let's process the cluster
00284       else 
00285       {
00286         if ( cluster_size > 1 ) 
00287           widen_cluster( cluster_start, cluster_end, cluster_size, jump_width, cluster_min );
00288         // clear the clusters
00289         cluster_start = ++cluster_end;
00290         cluster_size  = 1;
00291       }
00292     }
00293     // process the last cluster since doesn't get processed otherwise
00294     if ( cluster_size > 1 ) 
00295       widen_cluster( cluster_start, cluster_end, cluster_size, jump_width, cluster_min );
00296   }
00297 
00298 // ----------- widen_cluster -----------------
00299 //! Utility code used by widen_jumps
00300 void widen_cluster( two_d_iterator cluster_start, two_d_iterator cluster_end, 
00301   int cluster_size, double jump_width, double cluster_min )
00302 {
00303 
00304   if( ( cluster_size == 1 ) || ( cluster_start==cluster_end ) )  return;
00305   if ( cluster_size>4 ) 
00306   {
00307     Warning("two_d_list::widen_cluster",pastenum("Found ",cluster_size)+
00308       pastenum(" links in a cluster, check ENDL_EPSILON and ENDL_JUMP_WIDTH values or the data is broken at tag=",
00309       cluster_start->tag)); 
00310   }
00311 
00312   if ( cluster_start->tag <= cluster_min+jump_width )
00313   {
00314     two_d_iterator cluster_middle = cluster_start; 
00315     for (int ijump = 1; ijump<=cluster_size; ++ijump)
00316     {
00317       ++cluster_middle; 
00318       cluster_middle->tag += static_cast<double>(ijump)*jump_width;
00319     }
00320   }
00321   else
00322   {
00323     switch (cluster_size)
00324     {
00325       case 2:
00326         cluster_start->tag -= 0.5*jump_width;
00327         cluster_end->tag   += 0.5*jump_width;
00328         break;
00329       case 3:
00330         cluster_start->tag -= jump_width;
00331         cluster_end->tag   += jump_width;
00332         break;
00333       case 4:
00334         {
00335         two_d_iterator cluster_left = cluster_start;
00336         two_d_iterator cluster_right = cluster_end;
00337         ++cluster_left;--cluster_right;
00338         cluster_start->tag -= 1.5*jump_width;
00339         cluster_left->tag  -= 0.5*jump_width;
00340         cluster_right->tag += 0.5*jump_width;
00341         cluster_end->tag   += 1.5*jump_width;
00342         }
00343         break;
00344       default:
00345         Warning("two_d_list::widen_cluster",pastenum("Found ",cluster_size)+
00346         pastenum(" links in a cluster, check ENDL_EPSILON and ENDL_JUMP_WIDTH values or the data is broken at tag=",
00347         cluster_start->tag)); 
00348         return;
00349     }
00350   }
00351 }
00352 
00353   // ----------- weight_jumps -----------------
00354   //! Insert 1-d lists where the weight has jumps
00355   void weight_jumps()
00356   {
00357     one_d new_dist;  // for inserting new lists
00358 
00359     // widen the jumps in the weights
00360     _weight.widen_jumps( );
00361 
00362     // look for jumps in _weight
00363     two_d_iterator next_list = two_d_list::begin( );
00364     dd_list::iterator wt_ptr = _weight.begin( );
00365     for( ++wt_ptr ;
00366      ( wt_ptr != _weight.end( ) ) && ( next_list != two_d_list::end( ) );
00367         ++wt_ptr )
00368     {
00369       double wt_E = wt_ptr->x;
00370       // find the list with next larger incident energy
00371       for( ; next_list->E_in( ) < wt_E; ++next_list )
00372       {
00373         if( next_list == two_d_list::end( ) )
00374         {
00375           // We went past the maximum energy
00376           break;
00377         }
00378       }
00379 
00380       // we might already have a list at essentinally this energy
00381       double local_eps = ENDL_EPSILON( wt_E );
00382       if( next_list->E_in( ) < wt_E + 1.1*local_eps )
00383       {
00384         next_list->E_in( ) = wt_E;
00385         next_list->weight = wt_ptr->y;
00386       }
00387       else
00388       {
00389         // insert a new list
00390         two_d_iterator prev_list = next_list;
00391         --prev_list;
00392         insert( next_list, new_dist );
00393         two_d_iterator mid_list = prev_list;
00394         ++mid_list;
00395         mid_list->list_interp( wt_E, *prev_list, *next_list );
00396         mid_list->weight = wt_ptr->y;
00397       }
00398     }   
00399   }
00400 
00401   // ----------- thinit -----------------
00402   //! Thins the member lists to within some tolerance
00403   void thinit()
00404   {
00405     for( two_d_iterator list_ptr = two_d_list::begin(); 
00406       list_ptr != two_d_list::end(); ++list_ptr )
00407     {
00408       list_ptr->thinit( );
00409     }
00410   }
00411 
00412   // ----------- thicken -----------------
00413   //! Expands a list to some tolerance
00414   void thicken()
00415   {
00416     bool interp_OK;
00417     two_d_iterator list = two_d_list::begin();
00418     two_d_iterator next_list = list;
00419     ++next_list;
00420 
00421     // loop through the 1-d lists
00422     for( ; next_list != two_d_list::end(); ++list, next_list = list, ++next_list)
00423     {
00424       // keep checking until the interval passes the test
00425       for(;;)
00426       {
00427         interp_OK = check_interp(list);
00428 
00429         if(interp_OK)
00430         {
00431       break;
00432         }
00433         // *** a safely check ***
00434         const int Max_Size = 2000;
00435         if(two_d_list::size() >= Max_Size)
00436         {
00437           SevereError("two_d_list::thicken",
00438             pastenum("got ", Max_Size)+" links in 2-d thicken routine");
00439         }
00440       }
00441     }
00442   }
00443 
00444   // ----------- set_interp -----------------
00445   //! Sets the interpolation type for each list
00446   void set_interp()
00447   {
00448     int count = 0;
00449     int INT_index = 0;
00450     // loop through the list
00451     for(two_d_iterator this_link = two_d_list::begin();
00452       this_link != two_d_list::end(); ++this_link, ++count )
00453     {
00454       // do we change interpolation types?
00455       if ( count >= NBT[ INT_index ] )
00456       {
00457         ++INT_index;
00458         if ( INT_index >= INT.size() )
00459         {
00460           SevereError("two_d_list::set_interp"," INT_index overflow");
00461         }
00462       }
00463       this_link->interp_type = INT[ INT_index ];
00464     }   
00465   }
00466 
00467 
00468   // ----------- chop_highE -----------------
00469   //! Truncates list to the maximum incident energy, E_max
00470   void chop_highE( double Max_energy )
00471   {
00472     one_d new_list;
00473     two_d_iterator last_ptr = two_d_list::end();
00474     --last_ptr;
00475      two_d_iterator prev_ptr = last_ptr;
00476     --prev_ptr;
00477     insert( last_ptr, new_list );
00478     two_d_iterator new_ptr = last_ptr;
00479     --new_ptr;  // point to the new list
00480     new_ptr->list_interp( Max_energy, *prev_ptr, *last_ptr );
00481     erase( last_ptr, two_d_list::end() );  // delete the list above Max_energy
00482   }
00483 
00484   // ----------- check_interp -----------------
00485   //! Checks the accuracy of linear interpolation
00486   // between left_link and ++left_link
00487   // return: true if (worst relative error) <= tol_2d
00488   // if not, insert a new 1-d list good to within tol_1d
00489   bool check_interp(two_d_iterator left_list)
00490   {
00491     // is there a next list?
00492     two_d_iterator next_list = left_list;
00493     ++next_list;
00494     if( next_list == two_d_list::end( ) )
00495     {
00496       SevereError("two_d_list::check_interp",
00497         "attempt to check past the end of the list");
00498     }
00499 
00500     // lin-lin interpolation is always good
00501     if(( next_list->interp_type == 2 )||( next_list->interp_type == 22 ))
00502     {
00503       return( true );
00504     }
00505 
00506     // set width of jump as smallest dE we can see in output file
00507     double min_tol = ENDL_EPSILON( next_list->E_in( ) );
00508 
00509     // accept the link if dE is very small, but not too small
00510     if((next_list->E_in() - left_list->E_in() < dE_tol*next_list->E_in()) &&
00511         (next_list->E_in() - left_list->E_in() < min_tol) )
00512     {
00513       return(true);
00514     }
00515 
00516     // set the energy for an intermediate list
00517     double E_in;
00518     double dE;  // for histogram data
00519     if( next_list->interp_type == 1 ) // histogram data
00520     {
00521       two_d_iterator next_next_list = next_list;
00522       ++next_next_list;
00523       double dE_L = next_list->E_in( ) - left_list->E_in( );
00524       double dE_R = ( next_next_list == two_d_list::end( ) ) ?
00525         next_list->E_in( ) :
00526         next_next_list->E_in( ) - next_list->E_in( );
00527       dE = ( dE_L < dE_R ) ? dE_L : dE_R;
00528 
00529       if(dE <= 0.0)
00530       {
00531         left_list->print( );
00532         next_list->print( );
00533         next_next_list->print( );
00534         SevereError("two_d_list::check_interp","energies out of order");
00535       }
00536 
00537       E_in = next_list->E_in( );
00538       // if next_list->E_in is too close to next_next_list->E_in,
00539       // erase next_list
00540       if( dE_R <= min_tol )
00541       {
00542         erase( next_list );
00543         next_list = next_next_list;
00544       }
00545     }
00546     else
00547     {
00548       E_in = 0.5*( left_list->E_in( ) + next_list->E_in( ) );
00549     }
00550     // insert an intermediate list with incident energy E_in
00551     one_d mid_dist;
00552     insert(next_list, mid_dist);
00553 
00554     // expand the intermediate list and check it
00555     two_d_iterator mid_list = next_list;
00556     --mid_list;
00557     mid_list->E_in( ) = E_in;
00558 
00559     // get the weight for the model if this makes sense
00560     if ( ENDL.F ==5 )
00561     {
00562       mid_list->weight = _weight.evaluate( E_in );
00563     }
00564 
00565     if( next_list->interp_type == 1 ) // histogram data
00566     {
00567       // copy the next list
00568       mid_list->copy( *left_list );
00569       if( E_in > ENDL.Max_E_in - min_tol )
00570       {
00571         mid_list->E_in( ) = ENDL.Max_E_in;  // reset the incident energy tag
00572         erase( next_list );
00573       }
00574       else
00575       {
00576         mid_list->E_in( ) = E_in - min_tol;  // reset the incident energy tag
00577         next_list->interp_type = 2;  // now use lin-lin
00578         next_list->E_in( ) += min_tol;
00579       }
00580       mid_list->interp_type = 2;
00581       return( true );
00582     }
00583 
00584     mid_list->list_interp(E_in, *left_list, *next_list);
00585 
00586     bool is_OK;
00587     if(mid_list->bin_interp)
00588     {
00589       is_OK = mid_list->check_equiprob(*left_list, *next_list, tol_2d);
00590     }
00591     else
00592     {
00593       is_OK = mid_list->check_interp_2d(tol_2d, cut_off_2d);
00594     }
00595 
00596     if(is_OK)
00597     {
00598       // we don't need this list
00599       erase(mid_list);
00600     }
00601     else
00602     {
00603       if(debug_on)
00604       {
00605         cout << endl << " insert list " << endl;
00606         mid_list->print();
00607       }
00608     }
00609     return is_OK;
00610   }
00611 
00612   // ----------- collect_Ein -----------------
00613   //! Ensures that E_incident contains every incident energy
00614   void collect_Ein( list< double >& E_incident )
00615   {
00616     // DAB 2/16/05 -- previous algorithm did not always pick up 
00617     // all energies in the list. 
00618     for( two_d_iterator list_ptr = two_d_list::begin( ); 
00619          list_ptr != two_d_list::end( ); ++list_ptr ){
00620             E_incident.push_back(list_ptr->E_in());
00621       }
00622     E_incident.sort();   //STL gizmo to sort the list
00623     E_incident.unique(); //STL gizmo to remove duplicates
00624   }
00625 
00626   // ----------- copy_data -----------------
00627   //! copies the weights and the tabular data
00628   template <class One_d> void copy_data( two_d_list<One_d>& copy_from )
00629   {
00630     one_d new_dist;
00631 
00632     for( typename two_d_list<One_d>::iterator list_ptr = copy_from.begin( );
00633        list_ptr != copy_from.end( ); ++list_ptr )
00634     {
00635       insert( two_d_list::end( ), new_dist );
00636 
00637       // insert the data
00638       two_d_iterator target_ptr = two_d_list::end( );
00639       --target_ptr;
00640       target_ptr->E_in( ) = list_ptr->E_in( );
00641       target_ptr->weight = list_ptr->weight;
00642       target_ptr->interp_type = target_ptr->interp_type;
00643       for( dd_list::iterator link_ptr = list_ptr->begin();
00644         link_ptr != list_ptr->end(); ++link_ptr )
00645       {
00646         dd_link new_link = dd_link( link_ptr->x, link_ptr->y );
00647         target_ptr->insert( target_ptr->end(), new_link );
00648       }
00649     }
00650   }
00651 
00652   // ----------- fill_in_list -----------------
00653   //! Insures that 2d list has a link for each "E_in" in fill_from
00654   void fill_in_list( list< double >& fill_from )
00655   {
00656     const double fill_tol = Global.Value( "mf5_tol" );  // tolerance for mf5 data
00657 
00658     one_d mid_dist;  // for adding new distributions
00659     two_d_iterator mid_list;
00660     dd_link XYdata;  // for links in our lists
00661 
00662     two_d_iterator list_ptr = two_d_list::begin( );
00663     list< double >::iterator fill_ptr = fill_from.begin( );
00664     double list_E = list_ptr->E_in( );
00665     double fill_E = *fill_ptr;
00666     if( fill_E < ( 1 - fill_tol )*list_E )
00667     {
00668       insert( list_ptr, mid_dist );
00669       // make an initial list with 2 zero links
00670       mid_list = list_ptr;
00671       --mid_list;
00672       mid_list->weight = 0.0;
00673       mid_list->E_in( ) = fill_E;
00674       XYdata.x = 0.0;
00675       XYdata.y = 0.0;
00676       mid_list->insert( mid_list->end( ), XYdata );
00677 
00678       XYdata.x = 1.0e-8;
00679       mid_list->insert( mid_list->end( ), XYdata );
00680       list_ptr = mid_list; // set list_ptr to the start
00681       list_E = fill_E;
00682       ++fill_ptr;
00683     }
00684     else if( fill_E < ( 1 + fill_tol )*list_E )
00685     {
00686       // the lists start at the same incident energies, go to the next
00687       ++fill_ptr;
00688     }
00689     // Go through the fill_from list
00690     two_d_iterator prev_ptr = list_ptr;
00691     ++list_ptr;
00692     for( ; fill_ptr != fill_from.end();  )
00693     {
00694       // test for failure
00695       if( list_ptr == two_d_list::end() )
00696       {
00697         SevereError("two_d_list::fill_in_list",
00698             "Unexpected end of 2d list in fill_in_list");
00699       }
00700       list_E = list_ptr->E_in();
00701       fill_E = *fill_ptr;
00702       if ( fill_E < list_E )
00703       {
00704                 // insert a link at energy fill_E
00705         insert(list_ptr, mid_dist);
00706         // expand the intermediate list
00707         mid_list = list_ptr;
00708         --mid_list;
00709         mid_list->list_interp(fill_E, *prev_ptr, *list_ptr);
00710         if( ENDL.F == 5 )
00711                 {
00712           mid_list->weight = _weight.evaluate( fill_E );
00713                 }
00714         ++fill_ptr;
00715       }
00716       else if ( fill_E > list_E )
00717       {
00718         // go on to the next link
00719         prev_ptr = list_ptr;
00720         ++list_ptr;
00721       }
00722       else
00723       {
00724         // go on to the next link in both lists
00725         prev_ptr = list_ptr;
00726         ++list_ptr;
00727         ++fill_ptr;
00728       }
00729     }
00730   }
00731   // ----------- use_weight -----------------
00732   //! Scales the 1-d lists by their weights
00733   void use_weight( )
00734   {
00735     for( two_d_iterator list_ptr = two_d_list::begin(); 
00736       list_ptr != two_d_list::end(); ++list_ptr )
00737     {
00738       // do not zero out energy distributions
00739       if( list_ptr->weight > 0.0 )
00740       {
00741         *list_ptr *= list_ptr->weight;
00742       }
00743     }
00744   }
00745 
00746 };
00747 
00748 // ************* sum_lists function *****************
00749 // ----------- sum_lists -----------------
00750 //! Insures that 2d list has a link for each E_in in to_add
00751 template <class one_d, class One_d> void sum_lists(
00752   two_d_list<one_d>& sum, two_d_list<One_d>& to_add )
00753 {
00754   const double fill_tol = Global.Value( "mf5_tol" );  // tolerance for mf5 data
00755   typename two_d_list<one_d>::iterator sum_ptr = sum.begin();
00756   typename two_d_list<One_d>::iterator add_ptr = to_add.begin();
00757   for( ;
00758     sum_ptr != sum.end(), add_ptr != to_add.end();
00759     ++sum_ptr, ++add_ptr )
00760   {
00761     if ( abs(sum_ptr->E_in() - add_ptr->E_in()) >
00762       2*fill_tol*add_ptr->E_in() )
00763     {
00764       SevereError("two_d_list::sum_lists",
00765         pastenum("energy mismatch in sum_lists; E_in: ",sum_ptr->E_in())+
00766         pastenum(" ",add_ptr->E_in()));
00767     }
00768     // check for zero weight
00769     if( add_ptr->weight == 0.0 )
00770     {
00771       // usually do nothing
00772       if( sum_ptr->weight == 0.0 )
00773       {
00774         // take the average
00775         *sum_ptr += *add_ptr;
00776         // remove duplicates but not for gammas
00777         if( ( ENDL.F < 12 ) || ( ENDL.F > 15 ) )
00778         {
00779           sum_ptr->thinit();
00780         }
00781       }
00782     }
00783     else if( sum_ptr->weight == 0.0 )
00784     {
00785       // replace *sum_ptr data by *add_ptr data
00786       // save the next pointer and erase this one
00787       typename two_d_list<one_d>::iterator next_ptr = sum_ptr;
00788       ++next_ptr;
00789       sum.erase( sum_ptr );
00790       // insert a new one_d list
00791       one_d new_one_d;
00792       sum.insert( next_ptr, new_one_d );
00793       // reset sum_ptr
00794       sum_ptr = next_ptr;
00795       --sum_ptr;
00796       // copy the *add_ptr data
00797       sum_ptr->dd_list::copy( *add_ptr );
00798       sum_ptr->E_in( ) = add_ptr->E_in( );
00799       sum_ptr->weight = add_ptr->weight;
00800     }
00801     else
00802     {
00803       // add the lists with weight
00804       *add_ptr *= add_ptr->weight;
00805       *sum_ptr += *add_ptr;
00806       // remove duplicates but not for gammas
00807       if( ( ENDL.F < 12 ) || ( ENDL.F > 15 ) )
00808       {
00809         sum_ptr->thinit();
00810       }
00811       sum_ptr->weight += add_ptr->weight;
00812     }
00813   }
00814 }
00815 // ************* sum_lines function *****************
00816 // ----------- sum_lines -----------------
00817 //! Add 2 lists as discrete lines, no interpolation
00818 template <class one_d, class One_d> void sum_lines(
00819   two_d_list<one_d>& sum, two_d_list<One_d>& to_add )
00820 {
00821   const double fill_tol = Global.Value( "mf5_tol" );  // tolerance for mf5 data
00822   typename two_d_list<one_d>::iterator sum_ptr = sum.begin();
00823   typename two_d_list<One_d>::iterator add_ptr = to_add.begin();
00824   for( ;
00825     sum_ptr != sum.end(), add_ptr != to_add.end();
00826     ++sum_ptr, ++add_ptr )
00827   {
00828     if ( abs(sum_ptr->E_in() - add_ptr->E_in()) >
00829       2*fill_tol*add_ptr->E_in() )
00830     {
00831       SevereError("two_d_list::sum_lines",
00832         pastenum("energy mismatch in sum_lines; E_in: ",sum_ptr->E_in())+
00833         pastenum(" ",add_ptr->E_in()));
00834     }
00835     // add the lists
00836     sum_ptr->line_sum( *add_ptr, fill_tol );
00837   }
00838 }
00839 
00840 // *********** class two_d_table ****************
00841 //! Class that holds 2d tabular data as a function of incident energy
00842 // ----------- class two_d_table -----------------
00843 // derive this class from two_d_list
00844 class two_d_table : public two_d_list< one_d_table >
00845 {
00846 };
00847 
00848 
00849 #endif

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