MOAB: Mesh Oriented datABase  (version 5.4.1)
MBMesquite::NonGradient Class Reference

#include <NonGradient.hpp>

+ Inheritance diagram for MBMesquite::NonGradient:
+ Collaboration diagram for MBMesquite::NonGradient:

Public Member Functions

MESQUITE_EXPORT NonGradient (ObjectiveFunction *of)
MESQUITE_EXPORT NonGradient (ObjectiveFunction *of, MsqError &err)
virtual MESQUITE_EXPORT ~NonGradient ()
virtual MESQUITE_EXPORT std::string get_name () const
 Get string name for use in diagnostic and status output.
virtual MESQUITE_EXPORT PatchSetget_patch_set ()
MESQUITE_EXPORT bool project_gradient () const
MESQUITE_EXPORT void project_gradient (bool yesno)
int getDimension ()
double getThreshold ()
int getMaxNumEval ()
double getSimplexDiameterScale ()
void setDimension (int dimension)
void setThreshold (double threshold)
void setMaxNumEval (int maxNumEval)
void setExactPenaltyFunction (bool exactPenalty)
void setSimplexDiameterScale (double newScaleDiameter)
void getRowSum (int numRow, int numCol, std::vector< double > &matrix, std::vector< double > &rowSum)
bool testRowSum (int numRow, int numCol, double *matrix, double *rowSum)
double evaluate (double localArray[], PatchData &pd, MsqError &err)
int initSimplex (double edgeLength)
 edgeLenght is a length scale for the initial polytope.
void printPatch (const PatchData &pd, MsqError &err)
 Generic patch helper function only used by NonGradient.
int getPatchDimension (const PatchData &pd, MsqError &err)
 Generic patch helper function only used by NonGradient.
MESQUITE_EXPORT void set_debugging_level (int level)

Public Attributes

std::vector< double > simplex
 matrix stored by column as a std::vector
std::vector< double > height

Protected Member Functions

virtual MESQUITE_EXPORT void initialize (PatchData &pd, MsqError &err)
virtual MESQUITE_EXPORT void optimize_vertex_positions (PatchData &pd, MsqError &err)
virtual MESQUITE_EXPORT void initialize_mesh_iteration (PatchData &pd, MsqError &err)
virtual MESQUITE_EXPORT void terminate_mesh_iteration (PatchData &pd, MsqError &err)
virtual MESQUITE_EXPORT void cleanup ()

Private Member Functions

double amotry (std::vector< double > &, std::vector< double > &, double *, int, double, PatchData &, MsqError &err)
void printSimplex (std::vector< double > &, std::vector< double > &)
 NonGradient (const NonGradient &pd)
NonGradientoperator= (const NonGradient &pd)

Private Attributes

bool projectGradient
int mDimension
double mThreshold
double mTolerance
int mMaxNumEval
int mNonGradDebug
bool mUseExactPenaltyFunction
double mScaleDiameter

Detailed Description

This is an implementation of a derivative-free optimization algorithm Commonly referred to as the 'amoeba'. This implementation only works on patches containing one free vertex.

Definition at line 52 of file NonGradient.hpp.


Constructor & Destructor Documentation

Definition at line 72 of file NonGradient.cpp.

References MBMesquite::TerminationCriterion::add_iteration_limit(), MBMesquite::QualityImprover::get_inner_termination_criterion(), and set_debugging_level().

    : VertexMover( of ), PatchSetUser( true ), projectGradient( false ), mDimension( 0 ), mThreshold( 0.0 ),
      mTolerance( 0.0 ), mMaxNumEval( 0 ), mNonGradDebug( 0 ), mUseExactPenaltyFunction( true ), mScaleDiameter( 0.1 )
{
    set_debugging_level( 2 );
    // set the default inner termination criterion
    TerminationCriterion* default_crit = get_inner_termination_criterion();
    if( default_crit == NULL )
    {
        return;
    }
    else
    {
        default_crit->add_iteration_limit( 5 );
    }
}

Definition at line 89 of file NonGradient.cpp.

References MBMesquite::TerminationCriterion::add_iteration_limit(), MBMesquite::QualityImprover::get_inner_termination_criterion(), MBMesquite::MsqError::INVALID_STATE, MSQ_ERRRTN, MSQ_SETERR, and set_debugging_level().

    : VertexMover( of ), PatchSetUser( true ), projectGradient( false ), mDimension( 0 ), mThreshold( 0.0 ),
      mTolerance( 0.0 ), mMaxNumEval( 0 ), mNonGradDebug( 0 ), mUseExactPenaltyFunction( true ), mScaleDiameter( 0.1 )
{
    set_debugging_level( 2 );
    // set the default inner termination criterion
    TerminationCriterion* default_crit = get_inner_termination_criterion();
    if( default_crit == NULL )
    {
        MSQ_SETERR( err )
        ( "QualityImprover did not create a default inner "
          "termination criterion.",
          MsqError::INVALID_STATE );
        return;
    }
    else
    {
        default_crit->add_iteration_limit( 5 );MSQ_ERRRTN( err );
    }
}

Definition at line 61 of file NonGradient.hpp.

{}

Member Function Documentation

double MBMesquite::NonGradient::amotry ( std::vector< double > &  ,
std::vector< double > &  ,
double *  ,
int  ,
double  ,
PatchData ,
MsqError err 
) [private]

Definition at line 200 of file NonGradient.cpp.

References evaluate(), getDimension(), mNonGradDebug, and MSQ_PRINT.

Referenced by optimize_vertex_positions().

{
    int numRow = getDimension();
    // int numCol = numRow + 1;
    std::vector< double > ptry( numRow );  // does this make sense?
    double fac1 = ( 1.0 - fac ) / static_cast< double >( numRow );
    double fac2 = fac1 - fac;
    for( int row = 0; row < numRow; row++ )
    {
        ptry[row] = psum[row] * fac1 - p_simplex[row + ihi * numRow] * fac2;
    }
    if( mNonGradDebug >= 3 )
    {
        std::cout << "Try ";
    }
    MSQ_PRINT( 3 )( "Try" );
    double ytry = evaluate( &ptry[0], pd, err );  // value at trial point
    if( mNonGradDebug >= 3 )
    {
        std::cout << ytry << std::endl;
    }
    MSQ_PRINT( 3 )( "yTry" );
    if( ytry < p_height[ihi] )  // better than highest (worst)
    {
        p_height[ihi] = ytry;  // swap ihi and ytry
        for( int row = 0; row < numRow; row++ )
        {
            psum[row] += ( ptry[row] - p_simplex[row + ihi * numRow] );
            p_simplex[row + ihi * numRow] = ptry[row];
        }
    }
    return ytry;
}
void MBMesquite::NonGradient::cleanup ( ) [protected, virtual]

Implements MBMesquite::VertexMover.

Definition at line 562 of file NonGradient.cpp.

References mNonGradDebug.

{
    if( mNonGradDebug >= 2 )
    {
        std::cout << " - Executing NonGradient::iteration_end()" << std::endl;
    }
    // MSQ_PRINT(2)("\n - Executing NonGradient::iteration_end() \n");
}
double MBMesquite::NonGradient::evaluate ( double  localArray[],
PatchData pd,
MsqError err 
)

Definition at line 162 of file NonGradient.cpp.

References MBMesquite::TerminationCriterion::accumulate_inner(), MBMesquite::MsqError::BARRIER_VIOLATED, MBMesquite::MsqError::clear(), MBMesquite::MsqError::error_code(), MBMesquite::OFEvaluator::evaluate(), MBMesquite::QualityImprover::get_inner_termination_criterion(), MBMesquite::VertexMover::get_objective_function_evaluator(), MSQ_ERRZERO, MSQ_SETERR, mUseExactPenaltyFunction, NOT_IMPLEMENTED, MBMesquite::PatchData::num_free_vertices(), MBMesquite::PatchData::set_vertex_coordinates(), MBMesquite::PatchData::snap_vertex_to_domain(), value(), and MBMesquite::PatchData::vertex_by_index().

Referenced by amotry(), and optimize_vertex_positions().

{
    if( pd.num_free_vertices() > 1 )
    {
        MSQ_SETERR( err )
        ( "Only one free vertex per patch implemented", MsqError::NOT_IMPLEMENTED );
    }
    const size_t vertexIndex = 0;
    Vector3D originalVec     = pd.vertex_by_index( vertexIndex );
    Vector3D pointVec;
    for( int index = 0; index < 3; index++ )
    {
        pointVec[index] = point[index];
    }
    pd.set_vertex_coordinates( pointVec, vertexIndex, err );
    pd.snap_vertex_to_domain( vertexIndex, err );

    OFEvaluator& obj_func           = get_objective_function_evaluator();
    TerminationCriterion* term_crit = get_inner_termination_criterion();

    double value;
    bool feasible = obj_func.evaluate( pd, value, err );         // MSQ_ERRZERO(err);
    term_crit->accumulate_inner( pd, value, 0, err );            // MSQ_CHKERR(err);
    if( err.error_code() == err.BARRIER_VIOLATED ) err.clear();  // barrier violation not an error here
    MSQ_ERRZERO( err );
    pd.set_vertex_coordinates( originalVec, vertexIndex, err );
    if( !feasible && mUseExactPenaltyFunction )
    {  // "value" undefined btw
        double ensureFiniteRtol = .25;
        value                   = DBL_MAX * ensureFiniteRtol;
        // std::cout << " Infeasible patch: " << value << std::endl; printPatch( pd, err );
    }
    return value;
}
std::string MBMesquite::NonGradient::get_name ( ) const [virtual]

Get string name for use in diagnostic and status output.

Implements MBMesquite::Instruction.

Definition at line 62 of file NonGradient.cpp.

{
    return "NonGradient";
}

Reimplemented from MBMesquite::PatchSetUser.

Definition at line 67 of file NonGradient.cpp.

Definition at line 79 of file NonGradient.hpp.

References mDimension.

Referenced by amotry(), optimize_vertex_positions(), and printSimplex().

    {
        return ( mDimension );
    }

Definition at line 91 of file NonGradient.hpp.

References mMaxNumEval.

    {
        return ( mMaxNumEval );
    }

Generic patch helper function only used by NonGradient.

Definition at line 291 of file NonGradient.cpp.

References MBMesquite::PatchData::element_by_index(), MBMesquite::MsqMeshEntity::get_element_type(), MBMesquite::MsqError::INVALID_MESH, MSQ_SETERR, and MBMesquite::PatchData::num_elements().

Referenced by initialize_mesh_iteration().

{
    const size_t numElt = pd.num_elements();
    unsigned edimMax    = 0;  // in case numElt == 0
    for( size_t elementId = 0; elementId < numElt; elementId++ )
    {
        const MsqMeshEntity& element = pd.element_by_index( elementId );
        EntityTopology type          = element.get_element_type();
        unsigned edim                = TopologyInfo::dimension( type );
        if( elementId == 0 )
        {
            edimMax = edim;
        }
        else
        {
            if( edimMax != edim )
            {
                MSQ_SETERR( err )
                ( "A patch has elements of different dimensions", MsqError::INVALID_MESH );
                std::cout << "A patch has elements of different dimensions" << edimMax << "  " << edim << std::endl;
                if( edimMax < edim )
                {
                    edimMax = edim;
                }
            }
        }
    }
    return ( edimMax );
}
void MBMesquite::NonGradient::getRowSum ( int  numRow,
int  numCol,
std::vector< double > &  matrix,
std::vector< double > &  rowSum 
)

Definition at line 147 of file NonGradient.cpp.

Referenced by optimize_vertex_positions().

{
    for( int row = 0; row < numRow; row++ )
    {
        rowSum[row] = 0.;
    }
    for( int col = 0; col < numCol; col++ )
    {
        for( int row = 0; row < numRow; row++ )
        {
            rowSum[row] += matrix[row + col * numRow];
        }
    }
}

Definition at line 95 of file NonGradient.hpp.

References mScaleDiameter.

    {
        return ( mScaleDiameter );
    }

Definition at line 83 of file NonGradient.hpp.

References mThreshold.

    {
        return ( mThreshold );
    }
void MBMesquite::NonGradient::initialize ( PatchData pd,
MsqError err 
) [protected, virtual]

Implements MBMesquite::VertexMover.

Definition at line 321 of file NonGradient.cpp.

{}
void MBMesquite::NonGradient::initialize_mesh_iteration ( PatchData pd,
MsqError err 
) [protected, virtual]

Implements MBMesquite::VertexMover.

Definition at line 323 of file NonGradient.cpp.

References MBMesquite::PatchData::get_minmax_edge_length(), MBMesquite::PatchData::get_vertex_array(), getPatchDimension(), mNonGradDebug, mScaleDiameter, MSQ_PRINT, MSQ_SETERR, NOT_IMPLEMENTED, MBMesquite::PatchData::num_free_vertices(), MBMesquite::MsqError::OUT_OF_MEMORY, scale(), setDimension(), setMaxNumEval(), setThreshold(), and simplex.

{
    unsigned elementDimension = getPatchDimension( pd, err );  // to do: react to error
    unsigned dimension        = elementDimension * pd.num_free_vertices();
    // printPatch( pd, err );
    setDimension( dimension );
    int maxNumEval = 100 * dimension;  // 1. Make this a user parameter
    setMaxNumEval( maxNumEval );
    double threshold = 1.e-10;  // avoid division by zero
    setThreshold( threshold );
    double minEdgeLen = 0.0;
    double maxEdgeLen = 0.0;
    //  double ftol = 0.;
    if( dimension > 0 )
    {
        pd.get_minmax_edge_length( minEdgeLen, maxEdgeLen );
        // ftol = minEdgeLen * 1.e-4; // Turn off Amoeba convergence criterion
        if( mNonGradDebug >= 1 )
        {
            std::cout << "minimum edge length " << minEdgeLen << " maximum edge length " << maxEdgeLen << std::endl;
        }
        MSQ_PRINT( 3 )
        ( "minimum edge length %e    maximum edge length %e\n", minEdgeLen, maxEdgeLen );
    }
    //  setTolerance(ftol);
    unsigned numRow = dimension;
    unsigned numCol = numRow + 1;
    if( numRow * numCol <= simplex.max_size() )
    {
        simplex.assign( numRow * numCol, 0. );  // guard against previous simplex value
        double scale = minEdgeLen * mScaleDiameter;
        ;
        const MsqVertex* coord = pd.get_vertex_array( err );
        if( pd.num_free_vertices() > 1 )
        {
            MSQ_SETERR( err )
            ( "Only one free vertex per patch implemented", MsqError::NOT_IMPLEMENTED );
        }
        size_t index = 0;
        for( unsigned col = 0; col < numCol; col++ )
        {
            for( unsigned row = 0; row < numRow; row++ )
            {
                simplex[row + col * numRow] = coord[index][row];
                if( row == col - 1 )
                {
                    simplex[row + col * numRow] += scale / static_cast< double >( numCol );
                }
            }
        }
    }
    else
    {
        MSQ_SETERR( err )( "Use patch with fewer free vertices", MsqError::OUT_OF_MEMORY );
        if( mNonGradDebug >= 1 )
        {
            std::cout << "ERROR: Too many free vertices in patch" << std::endl;
        }
        // MSQ_PRINT(1)("ERROR: Too many free vertices in patch\n");
    }
}
int MBMesquite::NonGradient::initSimplex ( double  edgeLength)

edgeLenght is a length scale for the initial polytope.

NonGradient& MBMesquite::NonGradient::operator= ( const NonGradient pd) [private]
void MBMesquite::NonGradient::optimize_vertex_positions ( PatchData pd,
MsqError err 
) [protected, virtual]

Implements MBMesquite::VertexMover.

Definition at line 385 of file NonGradient.cpp.

References amotry(), evaluate(), MBMesquite::QualityImprover::get_inner_termination_criterion(), getDimension(), getRowSum(), mNonGradDebug, MSQ_FUNCTION_TIMER, MSQ_SETERR, NOT_IMPLEMENTED, MBMesquite::PatchData::num_free_vertices(), printSimplex(), MBMesquite::PatchData::set_vertex_coordinates(), simplex, MBMesquite::PatchData::snap_vertex_to_domain(), and MBMesquite::TerminationCriterion::terminate().

{
    MSQ_FUNCTION_TIMER( "NonGradient::optimize_vertex_positions" );
    int numRow = getDimension();
    int numCol = numRow + 1;
    std::vector< double > p_height( numCol );

    for( int col = 0; col < numCol; col++ )
    {
        p_height[col] = evaluate( &simplex[col * numRow], pd, err );  //  eval patch stuff
    }
    if( mNonGradDebug > 0 )
    {
        printSimplex( simplex, p_height );
    }

    // standardization
    TerminationCriterion* term_crit = get_inner_termination_criterion();
    // int maxNumEval = getMaxNumEval();
    // double threshold = getThreshold();

    //  double ftol = getTolerance();
    int ilo  = 0;  // height[ilo]<=...
    int inhi = 0;  //...<=height[inhi]<=
    int ihi  = 0;  //<=height[ihi]
    //  double rtol = 2.*ftol;
    double ysave;
    double ytry;
    bool afterEvaluation = false;
    std::vector< double > rowSum( numRow );
    getRowSum( numRow, numCol, simplex, rowSum );
    while( !( term_crit->terminate() ) )
    {

        if( mNonGradDebug > 0 )
        {
            printSimplex( simplex, p_height );
        }
        // std::cout << "rtol " << rtol << " ftol " << ftol << " MesquiteIter " <<
        // term_crit->get_iteration_count() << " Done " << term_crit->terminate() << std::endl;

        if( afterEvaluation )
        {
            // reflect highPt through opposite face
            // p_height[0] may vanish
            /*
                  if( !testRowSum( numRow, numCol, &simplex[0], &rowSum[0]) )
                  {
                    // Before uncommenting here and ...
                    //MSQ_SETERR(err)("Internal check sum test A failed", MsqError::INTERNAL_ERROR);
                    //MSQ_ERRRTN(err);
                  }
            */
            ytry = amotry( simplex, p_height, &rowSum[0], ihi, -1.0, pd, err );
            /*
                  if( !testRowSum( numRow, numCol, &simplex[0], &rowSum[0]) )
                  {
                     // ... here, determine a maxVal majorizing the previous as well as the current
               simplex.
                     //MSQ_SETERR(err)("Internal check sum test B failed",
               MsqError::INTERNAL_ERROR);
                     //MSQ_ERRRTN(err);
                  }
            */

            /*
                  if( p_height[0] == 0.)
                  {
                     MSQ_SETERR(err)("(B) Zero objective function value", MsqError::INTERNAL_ERROR);
                     exit(-1);
                  }
            */
            if( ytry <= p_height[ilo] )
            {
                ytry = amotry( simplex, p_height, &rowSum[0], ihi, -2.0, pd, err );
                if( mNonGradDebug >= 3 )
                {
                    std::cout << "Reflect and Expand from highPt " << ytry << std::endl;
                }
                // MSQ_PRINT(3)("Reflect and Expand from highPt : %e\n",ytry);
            }
            else
            {
                if( ytry >= p_height[inhi] )
                {
                    ysave = p_height[ihi];  // Contract along highPt
                    ytry  = amotry( simplex, p_height, &rowSum[0], ihi, 0.5, pd, err );
                    if( ytry >= ysave )
                    {  // contract all directions toward lowPt
                        for( int col = 0; col < numCol; col++ )
                        {
                            if( col != ilo )
                            {
                                for( int row = 0; row < numRow; row++ )
                                {
                                    rowSum[row] = 0.5 * ( simplex[row + col * numRow] + simplex[row + ilo * numRow] );
                                    simplex[row + col * numRow] = rowSum[row];
                                }
                                p_height[col] = evaluate( &rowSum[0], pd, err );
                                if( mNonGradDebug >= 3 )
                                {
                                    std::cout << "Contract all directions toward lowPt value( " << col
                                              << " ) = " << p_height[col] << " ilo = " << ilo << std::endl;
                                }
                                // MSQ_PRINT(3)("Contract all directions toward lowPt value( %d ) =
                                // %e    ilo = %d\n", col, p_height[col], ilo);
                            }
                        }
                    }
                }
            }     // ytri > h(ilo)
        }         // after evaluation
        ilo = 1;  // conditional operator or inline if
        ihi = p_height[0] > p_height[1] ? ( inhi = 1, 0 ) : ( inhi = 0, 1 );
        for( int col = 0; col < numCol; col++ )
        {
            if( p_height[col] <= p_height[ilo] )
            {
                ilo = col;  // ilo := argmin height
            }
            if( p_height[col] > p_height[ihi] )
            {
                inhi = ihi;
                ihi  = col;
            }
            else  // p_height[ihi] >= p_height[col]
                if( col != ihi && p_height[col] > p_height[inhi] ) inhi = col;
        }
        //    rtol=2.0*fabs( p_height[ihi]-p_height[ilo] )/
        //         ( fabs(p_height[ihi])+fabs(p_height[ilo])+threshold );
        afterEvaluation = true;
    }  //  while not converged

    // Always set to current best mesh.
    {
        if( ilo != 0 )
        {
            double yTemp  = p_height[0];
            p_height[0]   = p_height[ilo];  // height dimension numCol
            p_height[ilo] = yTemp;
            for( int row = 1; row < numRow; row++ )
            {
                yTemp                       = simplex[row];
                simplex[row]                = simplex[row + ilo * numRow];
                simplex[row + ilo * numRow] = yTemp;
            }
        }
    }
    if( pd.num_free_vertices() > 1 )
    {
        MSQ_SETERR( err )
        ( "Only one free vertex per patch implemented", MsqError::NOT_IMPLEMENTED );
    }

    Vector3D newPoint( &simplex[0] );
    size_t vertexIndex = 0;  // fix c.f. freeVertexIndex
    pd.set_vertex_coordinates( newPoint, vertexIndex, err );
    pd.snap_vertex_to_domain( vertexIndex, err );
    if( term_crit->terminate() )
    {
        if( mNonGradDebug >= 1 )
        {
            std::cout << "Optimization Termination OptStatus: Max Iter Exceeded" << std::endl;
        }
        // MSQ_PRINT(1)("Optimization Termination OptStatus: Max Iter Exceeded\n");
    }
}
void MBMesquite::NonGradient::printPatch ( const PatchData pd,
MsqError err 
)

Generic patch helper function only used by NonGradient.

Definition at line 240 of file NonGradient.cpp.

References MBMesquite::PatchData::get_vertex_array(), mNonGradDebug, MSQ_PRINT, MBMesquite::PatchData::num_corners(), MBMesquite::PatchData::num_elements(), MBMesquite::PatchData::num_free_vertices(), MBMesquite::PatchData::num_nodes(), MBMesquite::PatchData::num_slave_vertices(), and numVert.

{
    if( mNonGradDebug == 0 )
    {
        return;
    }
    const size_t numNode = pd.num_nodes();  // 27,  27 what?
    MSQ_PRINT( 3 )( "Number of Vertices: %d\n", (int)pd.num_nodes() );
    const size_t numVert = pd.num_free_vertices();  // 1
    MSQ_PRINT( 3 )( "Num Free = %d\n", (int)pd.num_free_vertices() );
    const size_t numSlaveVert = pd.num_slave_vertices();  // 0
    const size_t numCoin      = pd.num_corners();         // 64
    const MsqVertex* coord    = pd.get_vertex_array( err );
    MSQ_PRINT( 3 )( "Number of Vertices: %d\n", (int)pd.num_nodes() );

    std::cout << "Patch " << numNode << "  " << numVert << "  " << numSlaveVert << "  " << numCoin << std::endl;
    MSQ_PRINT( 3 )( "\n" );
    std::cout << "Coordinate ";
    std::cout << "         " << std::endl;
    for( size_t index = 0; index < numVert; index++ )
    {
        std::cout << coord[index][0] << "  " << coord[index][1] << "  " << coord[index][2] << std::endl;
    }
    // const size_t numElt = pd.num_elements();
    if( mNonGradDebug >= 3 )
    {
        std::cout << "Number of Elements: " << pd.num_elements() << std::endl;
    }
    MSQ_PRINT( 3 )( "Number of Elements: %d\n", (int)pd.num_elements() );
}
void MBMesquite::NonGradient::printSimplex ( std::vector< double > &  p_simplex,
std::vector< double > &  p_height 
) [private]

Definition at line 271 of file NonGradient.cpp.

References getDimension().

Referenced by optimize_vertex_positions().

{
    int numRow = getDimension();
    int numCol = numRow + 1;
    for( int col = 0; col < numCol; col++ )
    {
        // std::cout << "simplex[ " << col << "]= " ;
        for( int row = 0; row < numRow; row++ )
        {
            std::cout << p_simplex[row + col * numRow] << " ";
        }
        // std::cout << "           "  << p_height[col] << std::endl;
    }
    for( int col = 0; col < numCol; col++ )
    {
        std::cout << p_height[col] << " ";
    }
    std::cout << std::endl;
}

Definition at line 68 of file NonGradient.hpp.

References projectGradient.

    {
        return projectGradient;
    }

Definition at line 74 of file NonGradient.hpp.

References projectGradient.

    {
        projectGradient = yesno;
    }

Obtain diagnostic data during optimization off=level 0, ... level 3 = maximal

Definition at line 137 of file NonGradient.hpp.

References mNonGradDebug.

Referenced by main(), and NonGradient().

    {
        mNonGradDebug = level;
    }
void MBMesquite::NonGradient::setDimension ( int  dimension) [inline]

Definition at line 99 of file NonGradient.hpp.

References mDimension.

Referenced by initialize_mesh_iteration().

    {
        mDimension = dimension;
    }
void MBMesquite::NonGradient::setExactPenaltyFunction ( bool  exactPenalty) [inline]

Definition at line 115 of file NonGradient.hpp.

References mUseExactPenaltyFunction.

Referenced by main().

    {
        mUseExactPenaltyFunction = exactPenalty;
    }
void MBMesquite::NonGradient::setMaxNumEval ( int  maxNumEval) [inline]

Definition at line 111 of file NonGradient.hpp.

References mMaxNumEval.

Referenced by initialize_mesh_iteration().

    {
        mMaxNumEval = maxNumEval;
    }
void MBMesquite::NonGradient::setSimplexDiameterScale ( double  newScaleDiameter) [inline]

Definition at line 119 of file NonGradient.hpp.

References mScaleDiameter.

Referenced by main().

    {
        mScaleDiameter = newScaleDiameter;
    }
void MBMesquite::NonGradient::setThreshold ( double  threshold) [inline]

Definition at line 103 of file NonGradient.hpp.

References mThreshold.

Referenced by initialize_mesh_iteration().

    {
        mThreshold = threshold;
    }
void MBMesquite::NonGradient::terminate_mesh_iteration ( PatchData pd,
MsqError err 
) [protected, virtual]

Implements MBMesquite::VertexMover.

Definition at line 553 of file NonGradient.cpp.

References mNonGradDebug.

{
    if( mNonGradDebug >= 2 )
    {
        std::cout << "- Executing NonGradient::iteration_complete()" << std::endl;
    }
    // MSQ_PRINT(2)("\n - Executing NonGradient::iteration_complete() \n");
}
bool MBMesquite::NonGradient::testRowSum ( int  numRow,
int  numCol,
double *  matrix,
double *  rowSum 
)

Definition at line 110 of file NonGradient.cpp.

References mNonGradDebug, and MSQ_PRINT.

{
    bool same = true;
    std::vector< double > rowSum( numRow, 0. );
    double maxVal = 0.;
    for( int col = 0; col < numCol; col++ )
    {
        for( int row = 0; row < numRow; row++ )
        {
            rowSum[row] += matrix[row + col * numRow];
            if( fabs( matrix[row + col * numRow] ) > maxVal )
            {
                maxVal = fabs( matrix[row + col * numRow] );
            }
        }
    }
    double machEps    = 1.e-14 * static_cast< double >( numRow );  // better to use system parameters
    double upperBound = machEps * maxVal + 1.e-15;
    for( int row = 0; row < numRow; row++ )
    {
        if( fabs( rowSum[row] - oldRowSum[row] ) > upperBound )
        {
            same = false;
            if( mNonGradDebug >= 2 )
            {
                std::cout << " Warning: NonGradient Row Sum " << row << " Test: value " << rowSum[row]
                          << " Discrepancy " << rowSum[row] - oldRowSum[row] << " maxVal " << maxVal << " numRow "
                          << numRow << " numCol " << numCol << std::endl;
            }
            MSQ_PRINT( 2 )
            ( "NonGradient Row Sum [%d] Test failure: value %22.15e  difference %22.15e \n", row, rowSum[row],
              rowSum[row] - oldRowSum[row] );
        }
    }
    return ( same );
}

Member Data Documentation

std::vector< double > MBMesquite::NonGradient::height

Definition at line 130 of file NonGradient.hpp.

Definition at line 151 of file NonGradient.hpp.

Referenced by getDimension(), and setDimension().

Definition at line 154 of file NonGradient.hpp.

Referenced by getMaxNumEval(), and setMaxNumEval().

Definition at line 152 of file NonGradient.hpp.

Referenced by getThreshold(), and setThreshold().

Definition at line 153 of file NonGradient.hpp.

Definition at line 158 of file NonGradient.hpp.

Referenced by evaluate(), and setExactPenaltyFunction().

Definition at line 150 of file NonGradient.hpp.

Referenced by project_gradient().

std::vector< double > MBMesquite::NonGradient::simplex

matrix stored by column as a std::vector

Definition at line 129 of file NonGradient.hpp.

Referenced by initialize_mesh_iteration(), and optimize_vertex_positions().

List of all members.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines