MOAB: Mesh Oriented datABase  (version 5.4.1)
ConjugateGradient.cpp
Go to the documentation of this file.
00001 /* *****************************************************************
00002     MESQUITE -- The Mesh Quality Improvement Toolkit
00003 
00004     Copyright 2004 Sandia Corporation and Argonne National
00005     Laboratory.  Under the terms of Contract DE-AC04-94AL85000
00006     with Sandia Corporation, the U.S. Government retains certain
00007     rights in this software.
00008 
00009     This library is free software; you can redistribute it and/or
00010     modify it under the terms of the GNU Lesser General Public
00011     License as published by the Free Software Foundation; either
00012     version 2.1 of the License, or (at your option) any later version.
00013 
00014     This library is distributed in the hope that it will be useful,
00015     but WITHOUT ANY WARRANTY; without even the implied warranty of
00016     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00017     Lesser General Public License for more details.
00018 
00019     You should have received a copy of the GNU Lesser General Public License
00020     (lgpl.txt) along with this library; if not, write to the Free Software
00021     Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00022 
00023     [email protected], [email protected], [email protected],
00024     [email protected], [email protected], [email protected]
00025 
00026   ***************************************************************** */
00027 /*!
00028   \file   ConjugateGradient.cpp
00029   \brief
00030 
00031   The Conjugate Gradient class is a concrete vertex mover
00032   which performs conjugate gradient minimizaiton.
00033 
00034   \author Michael Brewer
00035   \date   2002-06-19
00036 */
00037 
00038 #include "ConjugateGradient.hpp"
00039 #include <cmath>
00040 #include "MsqDebug.hpp"
00041 #include "MsqTimer.hpp"
00042 //#include "MsqFreeVertexIndexIterator.hpp"
00043 
00044 namespace MBMesquite
00045 {
00046 
00047 extern int get_parallel_rank();
00048 extern int get_parallel_size();
00049 
00050 std::string ConjugateGradient::get_name() const
00051 {
00052     return "ConjugateGradient";
00053 }
00054 
00055 PatchSet* ConjugateGradient::get_patch_set()
00056 {
00057     return PatchSetUser::get_patch_set();
00058 }
00059 
00060 ConjugateGradient::ConjugateGradient( ObjectiveFunction* objective )
00061     : VertexMover( objective ), PatchSetUser( true ), pMemento( NULL ), conjGradDebug( 0 )
00062 {
00063 }
00064 
00065 ConjugateGradient::ConjugateGradient( ObjectiveFunction* objective, MsqError& err )
00066     : VertexMover( objective ), PatchSetUser( true ), pMemento( NULL ), conjGradDebug( 0 )
00067 {
00068     // Michael:: default to global?
00069     set_debugging_level( 0 );
00070     // set the default inner termination criterion
00071     TerminationCriterion* default_crit = get_inner_termination_criterion();
00072     if( default_crit == NULL )
00073     {
00074         MSQ_SETERR( err )
00075         ( "QualityImprover did not create a default inner "
00076           "termination criterion.",
00077           MsqError::INVALID_STATE );
00078         return;
00079     }
00080     else
00081     {
00082         default_crit->add_iteration_limit( 5 );MSQ_ERRRTN( err );
00083     }
00084 }
00085 
00086 ConjugateGradient::~ConjugateGradient()
00087 {
00088     // Checks that cleanup() has been called.
00089     assert( pMemento == NULL );
00090 }
00091 
00092 void ConjugateGradient::initialize( PatchData& pd, MsqError& err )
00093 {
00094     if( get_parallel_size() )
00095     {
00096         MSQ_DBGOUT( 2 ) << "\nP[" << get_parallel_rank() << "] "
00097                         << "o   Performing Conjugate Gradient optimization.\n";
00098     }
00099     else
00100     {
00101         MSQ_DBGOUT( 2 ) << "\no   Performing Conjugate Gradient optimization.\n";
00102     }
00103     pMemento = pd.create_vertices_memento( err );
00104 }
00105 
00106 void ConjugateGradient::initialize_mesh_iteration( PatchData& /*pd*/, MsqError& /*err*/ ) {}
00107 
00108 /*!Performs Conjugate gradient minimization on the PatchData, pd.*/
00109 void ConjugateGradient::optimize_vertex_positions( PatchData& pd, MsqError& err )
00110 {
00111     // pd.reorder();
00112 
00113     MSQ_FUNCTION_TIMER( "ConjugateGradient::optimize_vertex_positions" );
00114 
00115     Timer c_timer;
00116     size_t num_vert = pd.num_free_vertices();
00117     if( num_vert < 1 )
00118     {
00119         MSQ_DBGOUT( 1 ) << "\nEmpty free vertex list in ConjugateGradient\n";
00120         return;
00121     }
00122     /*
00123         //zero out arrays
00124       int zero_loop=0;
00125       while(zero_loop<arraySize){
00126         fGrad[zero_loop].set(0,0,0);
00127         pGrad[zero_loop].set(0,0,0);
00128         fNewGrad[zero_loop].set(0,0,0);
00129         ++zero_loop;
00130       }
00131     */
00132 
00133     // get OF evaluator
00134     OFEvaluator& objFunc = get_objective_function_evaluator();
00135 
00136     size_t ind;
00137     // Michael cull list:  possibly set soft_fixed flags here
00138 
00139     // MsqFreeVertexIndexIterator free_iter(pd, err);  MSQ_ERRRTN(err);
00140 
00141     double f = 0;
00142     // Michael, this isn't equivalent to CUBIT because we only want to check
00143     // the objective function value of the 'bad' elements
00144     // if invalid initial patch set an error.
00145     bool temp_bool = objFunc.update( pd, f, fGrad, err );
00146     assert( fGrad.size() == num_vert );
00147     if( MSQ_CHKERR( err ) ) return;
00148     if( !temp_bool )
00149     {
00150         MSQ_SETERR( err )
00151         ( "Conjugate Gradient not able to get valid gradient "
00152           "and function values on intial patch.",
00153           MsqError::INVALID_MESH );
00154         return;
00155     }
00156     double grad_norm = MSQ_MAX_CAP;
00157 
00158     if( conjGradDebug > 0 )
00159     {
00160         MSQ_PRINT( 2 )( "\nCG's DEGUB LEVEL = %i \n", conjGradDebug );
00161         grad_norm = Linf( arrptr( fGrad ), fGrad.size() );
00162         MSQ_PRINT( 2 )( "\nCG's FIRST VALUE = %f,grad_norm = %f", f, grad_norm );
00163         MSQ_PRINT( 2 )( "\n   TIME %f", c_timer.since_birth() );
00164         grad_norm = MSQ_MAX_CAP;
00165     }
00166 
00167     // Initializing pGrad (search direction).
00168     pGrad.resize( fGrad.size() );
00169     for( ind = 0; ind < num_vert; ++ind )
00170         pGrad[ind] = ( -fGrad[ind] );
00171 
00172     int j      = 0;            // total nb of step size changes ... not used much
00173     int i      = 0;            // iteration counter
00174     unsigned m = 0;            //
00175     double alp = MSQ_MAX_CAP;  // alp: scale factor of search direction
00176                                // we know inner_criterion is false because it was checked in
00177                                // loop_over_mesh before being sent here.
00178     TerminationCriterion* term_crit = get_inner_termination_criterion();
00179 
00180     // while ((i<maxIteration && alp>stepBound && grad_norm>normGradientBound)
00181     //     && !inner_criterion){
00182     while( !term_crit->terminate() )
00183     {
00184         ++i;
00185         // std::cout<<"\Michael delete i = "<<i;
00186         int k = 0;
00187         alp   = get_step( pd, f, k, err );
00188         j += k;
00189         if( conjGradDebug > 2 )
00190         {
00191             MSQ_PRINT( 2 )( "\n  Alp initial, alp = %20.18f", alp );
00192         }
00193 
00194         // if alp == 0, revert to steepest descent search direction
00195         if( alp == 0 )
00196         {
00197             for( m = 0; m < num_vert; ++m )
00198             {
00199                 pGrad[m] = ( -fGrad[m] );
00200             }
00201             alp = get_step( pd, f, k, err );
00202             j += k;
00203             if( conjGradDebug > 1 )
00204             {
00205                 MSQ_PRINT( 2 )( "\n CG's search direction reset." );
00206                 if( conjGradDebug > 2 ) MSQ_PRINT( 2 )( "\n  Alp was zero, alp = %20.18f", alp );
00207             }
00208         }
00209         if( alp != 0 )
00210         {
00211             pd.move_free_vertices_constrained( arrptr( pGrad ), num_vert, alp, err );MSQ_ERRRTN( err );
00212 
00213             if( !objFunc.update( pd, f, fNewGrad, err ) )
00214             {
00215                 MSQ_SETERR( err )
00216                 ( "Error inside Conjugate Gradient, vertices moved "
00217                   "making function value invalid.",
00218                   MsqError::INVALID_MESH );
00219                 return;
00220             }
00221             assert( fNewGrad.size() == (unsigned)num_vert );
00222 
00223             if( conjGradDebug > 0 )
00224             {
00225                 grad_norm = Linf( arrptr( fNewGrad ), num_vert );
00226                 MSQ_PRINT( 2 )
00227                 ( "\nCG's VALUE = %f,  iter. = %i,  grad_norm = %f,  alp = %f", f, i, grad_norm, alp );
00228                 MSQ_PRINT( 2 )( "\n   TIME %f", c_timer.since_birth() );
00229             }
00230             double s11 = 0;
00231             double s12 = 0;
00232             double s22 = 0;
00233             // free_iter.reset();
00234             // while (free_iter.next()) {
00235             //  m=free_iter.value();
00236             for( m = 0; m < num_vert; ++m )
00237             {
00238                 s11 += fGrad[m] % fGrad[m];
00239                 s12 += fGrad[m] % fNewGrad[m];
00240                 s22 += fNewGrad[m] % fNewGrad[m];
00241             }
00242 
00243             // Steepest Descent (takes 2-3 times as long as P-R)
00244             // double bet=0;
00245 
00246             // Fletcher-Reeves (takes twice as long as P-R)
00247             // double bet = s22/s11;
00248 
00249             // Polack-Ribiere
00250             double bet;
00251             if( !divide( s22 - s12, s11, bet ) ) return;  // gradient is zero
00252             // free_iter.reset();
00253             // while (free_iter.next()) {
00254             //  m=free_iter.value();
00255             for( m = 0; m < num_vert; ++m )
00256             {
00257                 pGrad[m] = ( -fNewGrad[m] + ( bet * pGrad[m] ) );
00258                 fGrad[m] = fNewGrad[m];
00259             }
00260             if( conjGradDebug > 2 )
00261             {
00262                 MSQ_PRINT( 2 )
00263                 ( " \nSEARCH DIRECTION INFINITY NORM = %e", Linf( arrptr( fNewGrad ), num_vert ) );
00264             }
00265 
00266         }  // end if on alp == 0
00267 
00268         term_crit->accumulate_patch( pd, err );MSQ_ERRRTN( err );
00269         term_crit->accumulate_inner( pd, f, arrptr( fGrad ), err );MSQ_ERRRTN( err );
00270     }  // end while
00271     if( conjGradDebug > 0 )
00272     {
00273         MSQ_PRINT( 2 )( "\nConjugate Gradient complete i=%i ", i );
00274         MSQ_PRINT( 2 )( "\n-  FINAL value = %f, alp=%4.2e grad_norm=%4.2e", f, alp, grad_norm );
00275         MSQ_PRINT( 2 )( "\n   FINAL TIME %f", c_timer.since_birth() );
00276     }
00277 }
00278 
00279 void ConjugateGradient::terminate_mesh_iteration( PatchData& /*pd*/, MsqError& /*err*/ )
00280 {
00281     //  cout << "- Executing ConjugateGradient::iteration_complete()\n";
00282 }
00283 
00284 void ConjugateGradient::cleanup()
00285 {
00286     //  cout << "- Executing ConjugateGradient::iteration_end()\n";
00287     fGrad.clear();
00288     pGrad.clear();
00289     fNewGrad.clear();
00290     // pMemento->~PatchDataVerticesMemento();
00291     delete pMemento;
00292     pMemento = NULL;
00293 }
00294 
00295 //! Computes a distance to move vertices given an initial position and search direction (stored in
00296 //! data member pGrad).
00297 /*!Returns alp, the double which scales the search direction vector
00298   which when added to the old nodal positions yields the new nodal
00299   positions.*/
00300 /*!\todo Michael NOTE:  ConjugateGradient::get_step's int &j is only
00301   to remain consisitent with CUBIT for an initial test.  It can be
00302   removed.*/
00303 
00304 double ConjugateGradient::get_step( PatchData& pd, double f0, int& j, MsqError& err )
00305 {
00306     // get OF evaluator
00307     OFEvaluator& objFunc = get_objective_function_evaluator();
00308 
00309     size_t num_vertices = pd.num_free_vertices();
00310     // initial guess for alp
00311     double alp = 1.0;
00312     int jmax   = 100;
00313     double rho = 0.5;
00314     // feasible=false implies the mesh is not in the feasible region
00315     bool feasible = false;
00316     int found     = 0;
00317     // f and fnew hold the objective function value
00318     double f    = 0;
00319     double fnew = 0;
00320     // Counter to avoid infinitly scaling alp
00321     j = 0;
00322     // save memento
00323     pd.recreate_vertices_memento( pMemento, err );
00324     // if we must check feasiblility
00325     // while step takes mesh into infeasible region and ...
00326     while( j < jmax && !feasible && alp > MSQ_MIN )
00327     {
00328         ++j;
00329         pd.set_free_vertices_constrained( pMemento, arrptr( pGrad ), num_vertices, alp, err );
00330         feasible = objFunc.evaluate( pd, f, err );
00331         if( err.error_code() == err.BARRIER_VIOLATED )
00332             err.clear();  // barrier violation does not represent an actual error here
00333         MSQ_ERRZERO( err );
00334         // if not feasible, try a smaller alp (take smaller step)
00335         if( !feasible )
00336         {
00337             alp *= rho;
00338         }
00339     }  // end while ...
00340 
00341     // if above while ended due to j>=jmax, no valid step was found.
00342     if( j >= jmax )
00343     {
00344         MSQ_PRINT( 2 )( "\nFeasible Point Not Found" );
00345         return 0.0;
00346     }
00347     // Message::print_info("\nOriginal f %f, first new f = %f, alp = %f",f0,f,alp);
00348     // if new f is larger than original, our step was too large
00349     if( f >= f0 )
00350     {
00351         j = 0;
00352         while( j < jmax && found == 0 )
00353         {
00354             ++j;
00355             alp *= rho;
00356             pd.set_free_vertices_constrained( pMemento, arrptr( pGrad ), num_vertices, alp, err );
00357             // Get new obj value
00358             // if patch is now invalid, then the feasible region is  convex or
00359             // we have an error.  For now, we assume an error.
00360             if( !objFunc.evaluate( pd, f, err ) )
00361             {
00362                 MSQ_SETERR( err )
00363                 ( "Non-convex feasiblility region found.", MsqError::INVALID_MESH );
00364             }
00365             pd.set_to_vertices_memento( pMemento, err );
00366             MSQ_ERRZERO( err );
00367             // if our step has now improved the objective function value
00368             if( f < f0 )
00369             {
00370                 found = 1;
00371             }
00372         }  //   end while j less than jmax
00373            // Message::print_info("\nj = %d found = %d f = %20.18f f0 = %20.18f\n",j,found,f,f0);
00374            // if above ended because of j>=jmax, take no step
00375         if( found == 0 )
00376         {
00377             // Message::print_info("alp = %10.8f, but returning zero\n",alp);
00378             alp = 0.0;
00379             return alp;
00380         }
00381 
00382         j = 0;
00383         // while shrinking the step improves the objFunc value further,
00384         // scale alp down.  Return alp, when scaling once more would
00385         // no longer improve the objFunc value.
00386         while( j < jmax )
00387         {
00388             ++j;
00389             alp *= rho;
00390             // step alp in search direction from original positions
00391             pd.set_free_vertices_constrained( pMemento, arrptr( pGrad ), num_vertices, alp, err );
00392             MSQ_ERRZERO( err );
00393 
00394             // get new objective function value
00395             if( !objFunc.evaluate( pd, fnew, err ) ) MSQ_SETERR( err )
00396             ( "Non-convex feasiblility region found while "
00397               "computing new f.",
00398               MsqError::INVALID_MESH );
00399             if( fnew < f )
00400             {
00401                 f = fnew;
00402             }
00403             else
00404             {
00405                 // Reset the vertices to original position
00406                 pd.set_to_vertices_memento( pMemento, err );
00407                 MSQ_ERRZERO( err );
00408                 alp /= rho;
00409                 return alp;
00410             }
00411         }
00412         // Reset the vertices to original position and return alp
00413         pd.set_to_vertices_memento( pMemento, err );
00414         MSQ_ERRZERO( err );
00415         return alp;
00416     }
00417     // else our new f was already smaller than our original
00418     else
00419     {
00420         j = 0;
00421         // check to see how large of step we can take
00422         while( j < jmax && found == 0 )
00423         {
00424             ++j;
00425             // scale alp up (rho must be less than 1)
00426             alp /= rho;
00427             // step alp in search direction from original positions
00428             pd.set_free_vertices_constrained( pMemento, arrptr( pGrad ), num_vertices, alp, err );
00429             MSQ_ERRZERO( err );
00430 
00431             feasible = objFunc.evaluate( pd, fnew, err );
00432             if( err.error_code() == err.BARRIER_VIOLATED )
00433                 err.clear();  // evaluate() error does not represent an actual problem here
00434             MSQ_ERRZERO( err );
00435             if( !feasible )
00436             {
00437                 alp *= rho;
00438                 // Reset the vertices to original position and return alp
00439                 pd.set_to_vertices_memento( pMemento, err );
00440                 MSQ_ERRZERO( err );
00441                 return alp;
00442             }
00443             if( fnew < f )
00444             {
00445                 f = fnew;
00446             }
00447             else
00448             {
00449                 found = 1;
00450                 alp *= rho;
00451             }
00452         }
00453 
00454         // Reset the vertices to original position and return alp
00455         pd.set_to_vertices_memento( pMemento, err );
00456         MSQ_ERRZERO( err );
00457         return alp;
00458     }
00459 }
00460 
00461 /*!Quadratic one-dimensional line search.*/
00462 /*
00463 double ConjugateGradient::get_step(PatchData &pd,double f0,int &j,
00464                                    MsqError &err){
00465   const double CGOLD = 0.3819660;
00466   const double ZEPS = 1.0e-10;
00467   int n=pd.num_free_vertices();
00468   MsqVertex* vertices=pd.get_vertex_array(err);
00469   double a,b,d,etemp,fb,fu,fv,fw,fx,p,q,r,tol,tol1,tol2,u,v,w,x,xm;
00470   double e=0.0;
00471   d=0.0;
00472   tol=.001;
00473   int iter, maxiter;
00474   maxiter=100;
00475   a=0;
00476   b=.125;
00477   int m=0;
00478   fb=f0-1.0;
00479   iter=0;
00480   //find b such that a b 'should' bracket the min
00481   while (fb<=f0 && iter<maxiter){
00482     ++iter;
00483     b*=2.0;
00484     for(m=0;m<n;++m){
00485       mCoord[m]=mCoord[m] + (b*pGrad[m]);
00486       vertices[m]=(mCoord[m]);
00487     }
00488     fb=objFunc->evaluate(pd,err);
00489   }
00490   iter=0;
00491   x=w=v=(b/2.0);
00492   for(m=0;m<n;++m){
00493     mCoord[m]=mCoord[m] + (x*pGrad[m]);
00494     vertices[m]=(mCoord[m]);
00495   }
00496   fw=fv=fx=objFunc->evaluate(pd,err);
00497   for(iter=0;iter<maxiter;++iter){
00498       //Message::print_info("a=%f,b=%f,x=%f,iter=%i\n",a,b,x,iter);
00499     xm=(a+b)*.5;
00500     tol2=2.0*(tol1=tol*fabs(x)+ZEPS);
00501     if(fabs(x-xm)<= (tol2-0.5*(b-a))){
00502       return x;
00503     }
00504     if(fabs(e)>tol1){
00505       r=(x-w)*(fx-fv);
00506       q=(x-v)*(fx-fw);
00507       p=(x-v)*q-(x-w)*r;
00508       q=2.0*(q-r);
00509       if(q>0.0)
00510         p=-p;
00511       q=fabs(q);
00512       etemp=e;
00513       e=d;
00514       if(fabs(p)>=fabs(0.5*q*etemp)||(p<=q*(a-x))||(p>=q*(b-x))){
00515         d=CGOLD*(e=(x>=xm?a-x:b-x));
00516       }
00517       else{
00518         d=p/q;
00519         u=x+d;
00520         if(u-a<tol2||b-u<tol2)
00521         {
00522           if(tol1<0.0)
00523             d=x-xm;
00524           else
00525             d=xm-x;
00526         }
00527       }
00528     }
00529 
00530     else{
00531       d=CGOLD*(e=(x>=xm?a-x:b-x));
00532     }
00533     if(tol<0.0)
00534       u=(fabs(d)>=tol1?x+d:x-d);
00535     else
00536       u=(fabs(d)>=tol1?x+d:x+d);
00537     for(m=0;m<n;++m){
00538       mCoord[m]=mCoord[m] + (u*pGrad[m]);
00539       vertices[m]=(mCoord[m]);
00540     }
00541     fu=objFunc->evaluate(pd,err);
00542     if(fu<fx){
00543       if(u>=x)
00544         a=x;
00545       else
00546         b=x;
00547       v=w;
00548       w=x;
00549       x=u;
00550       fv=fw;
00551       fw=fx;
00552       fx=fu;
00553     }
00554     else{
00555       if(u<x)
00556         a=u;
00557       else
00558         b=u;
00559       if(fu<=fw||w==x){
00560         v=w;
00561         w=u;
00562         fv=fw;
00563         fw=fu;
00564       }
00565       else if (fu<=fv||v==x||v==w){
00566         v=u;
00567         fv=fu;
00568       }
00569     }
00570   }
00571   for(m=0;m<n;++m){
00572     vertices[m]=(mCoord[m]);
00573   }
00574     //PRINT_WARNING("TOO MANY ITERATIONS IN QUADRATIC LINE SEARCH");
00575   return x;
00576 }
00577 */
00578 
00579 }  // namespace MBMesquite
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines