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

(variance)^2 template More...

#include <VarianceTemplate.hpp>

+ Inheritance diagram for MBMesquite::VarianceTemplate:
+ Collaboration diagram for MBMesquite::VarianceTemplate:

Public Member Functions

MESQUITE_EXPORT VarianceTemplate (QualityMetric *qm=0)
MESQUITE_EXPORT VarianceTemplate (const VarianceTemplate &copy)
 copy constructor
virtual MESQUITE_EXPORT ~VarianceTemplate ()
virtual MESQUITE_EXPORT bool evaluate (EvalType type, PatchData &pd, double &value_out, bool free, MsqError &err)
 Evaluate objective function for specified patch.
virtual MESQUITE_EXPORT bool evaluate_with_gradient (EvalType type, PatchData &pd, double &value_out, std::vector< Vector3D > &grad_out, MsqError &err)
 Evaluate objective function and gradient for specified patch.
virtual MESQUITE_EXPORT bool evaluate_with_Hessian_diagonal (EvalType type, PatchData &pd, double &value_out, std::vector< Vector3D > &grad_out, std::vector< SymMatrix3D > &hess_diag_out, MsqError &err)
 Evaluate objective function and diagonal blocks of Hessian for specified patch.
virtual MESQUITE_EXPORT
ObjectiveFunction
clone () const
 Create copy with same state.
virtual MESQUITE_EXPORT void clear ()

Private Member Functions

void accumulate (double sum, double sqr_sum, size_t count, EvalType type, double &result_sum, double &result_sqr, size_t &global_count)
 Handle EvalType for all eval functions, return OF value.

Private Attributes

size_t mCount
double mSum
double mSqrSum
size_t saveCount
double saveSum
double saveSqrSum
std::vector< size_t > qmHandles
std::vector< size_t > mIndices
std::vector< Vector3DmGradient
std::vector< Vector3DgradSum
std::vector< SymMatrix3DmHessDiag
std::vector< SymMatrix3DhessSum

Detailed Description

(variance)^2 template

This class implements an objective function that is the variance of the quality metric evalutations.

Definition at line 47 of file VarianceTemplate.hpp.


Constructor & Destructor Documentation

Definition at line 51 of file VarianceTemplate.hpp.

References clear().

Referenced by clone().

copy constructor

Define a copy constructor because the compiler-provided default one would also copy the temporary arrays, which would be a waste of time.

Definition at line 63 of file VarianceTemplate.hpp.

        : ObjectiveFunctionTemplate( copy ), mCount( copy.mCount ), mSum( copy.mSum ), mSqrSum( copy.mSqrSum ),
          saveCount( copy.saveCount ), saveSum( copy.saveSum ), saveSqrSum( copy.saveSqrSum )
    {
    }

Definition at line 70 of file VarianceTemplate.hpp.

{}

Member Function Documentation

void MBMesquite::VarianceTemplate::accumulate ( double  sum,
double  sqr_sum,
size_t  count,
EvalType  type,
double &  result_sum,
double &  result_sqr,
size_t &  global_count 
) [private]

Handle EvalType for all eval functions, return OF value.

This function implements the common handling of the EvalType argument for all forms of the 'evaluate' method.

NOTE: This function modifies accumulated values depenending on the value of EvalType.

Parameters:
sumThe sum over the current patch
sqr_sumThe sum of squares over the current patch
countThe number of qm evaluations for the current patch
typeThe evaluation type passed to 'evaluate'
global_countThe total, accumulated number of QM evaluations
result_sumThe sum term of the variance
result_sqrThe sum of squares term of the variance

Definition at line 55 of file VarianceTemplate.cpp.

References MBMesquite::ObjectiveFunction::ACCUMULATE, MBMesquite::ObjectiveFunction::CALCULATE, mCount, mSqrSum, mSum, MBMesquite::ObjectiveFunction::SAVE, saveCount, saveSqrSum, saveSum, moab::sum(), MBMesquite::ObjectiveFunction::TEMPORARY, and MBMesquite::ObjectiveFunction::UPDATE.

Referenced by evaluate(), evaluate_with_gradient(), and evaluate_with_Hessian_diagonal().

{
    switch( type )
    {
        case CALCULATE:
            result_sum   = sum;
            result_sqr   = sqr_sum;
            global_count = count;
            break;

        case ACCUMULATE:
            result_sum   = mSum += sum;
            result_sqr   = mSqrSum += sqr_sum;
            global_count = mCount += count;
            break;

        case SAVE:
            saveSum      = sum;
            saveSqrSum   = sqr_sum;
            saveCount    = count;
            result_sum   = mSum;
            result_sqr   = mSqrSum;
            global_count = mCount;
            break;

        case UPDATE:
            mSum -= saveSum;
            mSqrSum -= saveSqrSum;
            mCount -= saveCount;
            result_sum = mSum += saveSum = sum;
            result_sqr = mSqrSum += saveSqrSum = sqr_sum;
            global_count = mCount += saveCount = count;
            break;

        case TEMPORARY:
            result_sum   = mSum - saveSum + sum;
            result_sqr   = mSqrSum - saveSqrSum + sqr_sum;
            global_count = mCount + count - saveCount;
            break;
    }
}

Clear any values accumulated for BCD-related eval calls

Implements MBMesquite::ObjectiveFunction.

Definition at line 47 of file VarianceTemplate.cpp.

References mCount, mSqrSum, mSum, saveCount, saveSqrSum, and saveSum.

Referenced by VarianceTemplate().

{
    mCount = 0;
    mSum = mSqrSum = 0;
    saveCount      = 0;
    saveSum = saveSqrSum = 0;
}

Create copy with same state.

Create a new instance of the objective function that is a copy of the callee with the same accumulated values, parameters, etc.

Implements MBMesquite::ObjectiveFunction.

Reimplemented in MBMesquite::StdDevTemplate.

Definition at line 42 of file VarianceTemplate.cpp.

References VarianceTemplate().

{
    return new VarianceTemplate( *this );
}
bool MBMesquite::VarianceTemplate::evaluate ( EvalType  type,
PatchData pd,
double &  value_out,
bool  free,
MsqError err 
) [virtual]

Evaluate objective function for specified patch.

Either evaluate the objective function over the passed patch or update the accumulated, global objective function value for changes in the passed patch, depending on the value of the EvalType.

Parameters:
typeEvaluation type.
pdThe patch.
value_outThe passed-back value of the objective fuction.
freeIf true, incorporate the quality metric values only for those metric evaluations that depend on at least one free vertex
Returns:
false if any QualityMetric evaluation returned false, true otherwise.

Implements MBMesquite::ObjectiveFunction.

Reimplemented in MBMesquite::StdDevTemplate.

Definition at line 103 of file VarianceTemplate.cpp.

References MBMesquite::ObjectiveFunction::ACCUMULATE, accumulate(), MBMesquite::QualityMetric::evaluate(), MBMesquite::QualityMetric::get_evaluations(), MBMesquite::QualityMetric::get_negate_flag(), MBMesquite::ObjectiveFunctionTemplate::get_quality_metric(), MBMesquite::QualityMetric::get_single_pass(), MSQ_CHKERR, MSQ_ERRFALSE, n, qmHandles, moab::sum(), and value().

{
    QualityMetric* qm = get_quality_metric();
    if( type == ObjectiveFunction::ACCUMULATE )
        qm->get_single_pass( pd, qmHandles, free, err );
    else
        qm->get_evaluations( pd, qmHandles, free, err );
    MSQ_ERRFALSE( err );

    // calculate OF value for just the patch
    std::vector< size_t >::const_iterator i;
    double value, sum = 0.0, sqr = 0.0;
    for( i = qmHandles.begin(); i != qmHandles.end(); ++i )
    {
        bool result = qm->evaluate( pd, *i, value, err );
        if( MSQ_CHKERR( err ) || !result ) return false;

        sum += value;
        sqr += value * value;
    }

    // get overall OF value, update member data, etc.
    size_t n;
    accumulate( sum, sqr, qmHandles.size(), type, sum, sqr, n );

    // don't divide by zero
    if( n < 1 )
    {
        value_out = 0.0;
        return true;
    }

    const double rms = sqr / n;
    const double avg = sum / n;
    value_out        = qm->get_negate_flag() * ( rms - avg * avg );
    return true;
}
bool MBMesquite::VarianceTemplate::evaluate_with_gradient ( EvalType  eval_type,
PatchData pd,
double &  OF_val,
std::vector< Vector3D > &  grad,
MsqError err 
) [virtual]

Evaluate objective function and gradient for specified patch.

Either evaluate the objective function over the passed patch or update the accumulated, global objective function value for changes in the passed patch, depending on the value of the EvalType.

The default implementation of this function will use the value-only variation of the evaluate method and numerical approximation to calculate gradients. Whenever possible, objective function implementations should provide more efficient analyical gradient calculations.

Parameters:
typeEvaluation type.
pdThe patch.
value_outThe passed-back value of the objective fuction.
grad_outThe gradient of the OF wrt the coordinates of each *free* vertex in the patch.
Returns:
false if any QualityMetric evaluation returned false, true otherwise.

Numerically Calculates the gradient of the ObjectiveFunction for the free vertices in the patch. Returns 'false' if the patch is outside of a required feasible region, returns 'ture' otherwise. The behavior of the function depends on the value of the boolean useLocalGradient. If useLocalGradient is set to 'true', compute_numerical_gradient creates a sub-patch around a free vertex, and then perturbs that vertex in one of the coordinate directions. Only the ObjectiveFunction value on the local sub-patch is used in the computation of the gradient. Therefore, useLocalGradient should only be set to 'true' for ObjectiveFunctions which can use this method. Unless the concrete ObjectiveFunction sets useLocalGradient to 'true' in its constructor, the value will be 'false'. In this case, the objective function value for the entire patch is used in the calculation of the gradient. This is computationally expensive, but it is numerically correct for all (C_1) functions.

Parameters:
pdPatchData on which the gradient is taken.
gradArray of Vector3D of length the number of vertices used to store gradient.
OF_valwill be set to the objective function value.

Reimplemented from MBMesquite::ObjectiveFunction.

Reimplemented in MBMesquite::StdDevTemplate.

Definition at line 141 of file VarianceTemplate.cpp.

References accumulate(), MBMesquite::QualityMetric::evaluate_with_gradient(), MBMesquite::QualityMetric::get_evaluations(), MBMesquite::QualityMetric::get_negate_flag(), MBMesquite::ObjectiveFunctionTemplate::get_quality_metric(), gradSum, mGradient, mIndices, MSQ_CHKERR, MSQ_ERRFALSE, n, MBMesquite::PatchData::num_free_vertices(), MBMesquite::OF_FREE_EVALS_ONLY, qmHandles, moab::sum(), and value().

{
    QualityMetric* qm = get_quality_metric();
    qm->get_evaluations( pd, qmHandles, OF_FREE_EVALS_ONLY, err );
    MSQ_ERRFALSE( err );

    // zero gradient data
    grad_out.clear();  // store sum of metric * gradient of metric, and later OF gradient
    grad_out.resize( pd.num_free_vertices(), Vector3D( 0.0, 0.0, 0.0 ) );
    gradSum.clear();  // store sum of gradients of metrics
    gradSum.resize( pd.num_free_vertices(), Vector3D( 0.0, 0.0, 0.0 ) );

    // calculate OF value and gradient for just the patch
    Matrix3D op;
    std::vector< size_t >::const_iterator i;
    double value, sum = 0.0, sqr = 0.0;
    for( i = qmHandles.begin(); i != qmHandles.end(); ++i )
    {
        bool result = qm->evaluate_with_gradient( pd, *i, value, mIndices, mGradient, err );
        if( MSQ_CHKERR( err ) || !result ) return false;
        if( fabs( value ) < DBL_EPSILON ) continue;

        sum += value;
        sqr += value * value;

        for( size_t j = 0; j < mIndices.size(); ++j )
        {
            const size_t r = mIndices[j];
            gradSum[r] += mGradient[j];
            mGradient[j] *= value;
            grad_out[r] += mGradient[j];
        }
    }

    // update accumulated values (if doing BCD)
    size_t n;
    accumulate( sum, sqr, qmHandles.size(), type, sum, sqr, n );

    // don't divide by zero
    if( n < 1 )
    {
        value_out = 0.0;
        grad_out.clear();
        grad_out.resize( pd.num_free_vertices(), Vector3D( 0.0, 0.0, 0.0 ) );
        return true;
    }

    // calculate OF value
    const double rms = sqr / n;
    const double avg = sum / n;
    value_out        = qm->get_negate_flag() * ( rms - avg * avg );

    // Finish calculation of gradient and Hessian
    const double f = qm->get_negate_flag() * 2.0 / n;
    for( size_t k = 0; k < pd.num_free_vertices(); ++k )
    {
        gradSum[k] *= avg;
        grad_out[k] -= gradSum[k];
        grad_out[k] *= f;
    }

    return true;
}
bool MBMesquite::VarianceTemplate::evaluate_with_Hessian_diagonal ( EvalType  type,
PatchData pd,
double &  value_out,
std::vector< Vector3D > &  grad_out,
std::vector< SymMatrix3D > &  hess_diag_out,
MsqError err 
) [virtual]

Evaluate objective function and diagonal blocks of Hessian for specified patch.

Either evaluate the objective function over the passed patch or update the accumulated, global objective function value for changes in the passed patch, depending on the value of the EvalType.

The default implementation of this function evaluate the entire Hessian and discard non-diagonal portions. Concrete objective functions should provide a more efficient implementation that evaluates and accumulates only the required terms.

Parameters:
typeEvaluation type.
pdThe patch.
value_outThe passed-back value of the objective fuction.
grad_outThe gradient of the OF wrt the coordinates of each *free* vertex in the patch.
hess_diag_outThe diagonal blocks of a Hessian. I.e. Decompose the Hessian into 3x3 submatrices and return only the submatrices (blocks) along the diagonal.
Returns:
false if any QualityMetric evaluation returned false, true otherwise.

Reimplemented from MBMesquite::ObjectiveFunction.

Reimplemented in MBMesquite::StdDevTemplate.

Definition at line 209 of file VarianceTemplate.cpp.

References accumulate(), MBMesquite::QualityMetric::evaluate_with_Hessian_diagonal(), MBMesquite::QualityMetric::get_evaluations(), MBMesquite::QualityMetric::get_negate_flag(), MBMesquite::ObjectiveFunctionTemplate::get_quality_metric(), gradSum, hessSum, mGradient, mHessDiag, mIndices, MSQ_CHKERR, MSQ_ERRFALSE, n, MBMesquite::PatchData::num_free_vertices(), MBMesquite::OF_FREE_EVALS_ONLY, MBMesquite::outer(), qmHandles, moab::sum(), and value().

{
    QualityMetric* qm = get_quality_metric();
    qm->get_evaluations( pd, qmHandles, OF_FREE_EVALS_ONLY, err );
    MSQ_ERRFALSE( err );

    // zero gradient and Hessian data
    grad_out.clear();  // store sum of metric * gradient of metric, and later OF gradient
    grad_out.resize( pd.num_free_vertices(), Vector3D( 0.0, 0.0, 0.0 ) );
    gradSum.clear();  // store sum of gradients of metrics
    gradSum.resize( pd.num_free_vertices(), Vector3D( 0.0, 0.0, 0.0 ) );
    hessSum.clear();  // store sum of Hessians of metrics
    hessSum.resize( pd.num_free_vertices(), SymMatrix3D( 0.0 ) );
    hess_diag_out.clear();  // store sum of metric * outer_product(metric gradient), and later OF Hessian
    hess_diag_out.resize( pd.num_free_vertices(), SymMatrix3D( 0.0 ) );

    // calculate OF value and gradient for just the patch
    Matrix3D op;
    std::vector< size_t >::const_iterator i;
    double value, sum = 0.0, sqr = 0.0;
    for( i = qmHandles.begin(); i != qmHandles.end(); ++i )
    {
        bool result = qm->evaluate_with_Hessian_diagonal( pd, *i, value, mIndices, mGradient, mHessDiag, err );
        if( MSQ_CHKERR( err ) || !result ) return false;
        if( fabs( value ) < DBL_EPSILON ) continue;

        sum += value;
        sqr += value * value;

        for( size_t j = 0; j < mIndices.size(); ++j )
        {
            const size_t r = mIndices[j];

            hessSum[r] += mHessDiag[j];
            hess_diag_out[r] += outer( mGradient[j] );
            mHessDiag[j] *= value;
            hess_diag_out[r] += mHessDiag[j];

            gradSum[r] += mGradient[j];
            mGradient[j] *= value;
            grad_out[r] += mGradient[j];
        }
    }

    // update accumulated values (if doing BCD)
    size_t n;
    accumulate( sum, sqr, qmHandles.size(), type, sum, sqr, n );

    // don't divide by zero
    if( n < 1 )
    {
        value_out = 0.0;
        grad_out.clear();
        grad_out.resize( pd.num_free_vertices(), Vector3D( 0.0, 0.0, 0.0 ) );
        hess_diag_out.clear();
        hess_diag_out.resize( pd.num_free_vertices(), SymMatrix3D( 0.0 ) );
        return true;
    }

    // calculate OF value
    const double rms = sqr / n;
    const double avg = sum / n;
    value_out        = qm->get_negate_flag() * ( rms - avg * avg );

    // Finish calculation of gradient and Hessian
    const double inv_n = 1.0 / n;
    const double f     = qm->get_negate_flag() * 2.0 / n;
    for( size_t k = 0; k < pd.num_free_vertices(); ++k )
    {
        hessSum[k] *= avg;
        hess_diag_out[k] -= hessSum[k];
        hess_diag_out[k] -= inv_n * outer( gradSum[k] );
        hess_diag_out[k] *= f;

        gradSum[k] *= avg;
        grad_out[k] -= gradSum[k];
        grad_out[k] *= f;
    }

    return true;
}

Member Data Documentation

std::vector< Vector3D > MBMesquite::VarianceTemplate::gradSum [mutable, private]

Definition at line 132 of file VarianceTemplate.hpp.

Referenced by evaluate_with_gradient(), and evaluate_with_Hessian_diagonal().

std::vector< SymMatrix3D > MBMesquite::VarianceTemplate::hessSum [mutable, private]

Definition at line 134 of file VarianceTemplate.hpp.

Referenced by evaluate_with_Hessian_diagonal().

The number of accumulated entires

Definition at line 120 of file VarianceTemplate.hpp.

Referenced by accumulate(), and clear().

std::vector< Vector3D > MBMesquite::VarianceTemplate::mGradient [mutable, private]

Temporary storage for qm gradient

Definition at line 132 of file VarianceTemplate.hpp.

Referenced by evaluate_with_gradient(), and evaluate_with_Hessian_diagonal().

std::vector< SymMatrix3D > MBMesquite::VarianceTemplate::mHessDiag [mutable, private]

Temporary storage for qm Hessian diagonal data

Definition at line 134 of file VarianceTemplate.hpp.

Referenced by evaluate_with_Hessian_diagonal().

std::vector< size_t > MBMesquite::VarianceTemplate::mIndices [mutable, private]

Temporary storage for qm vertex indices

Definition at line 130 of file VarianceTemplate.hpp.

Referenced by evaluate_with_gradient(), and evaluate_with_Hessian_diagonal().

The running sum of the square of QM values

Definition at line 122 of file VarianceTemplate.hpp.

Referenced by accumulate(), and clear().

The runnnig sum of the qualtiy metric valuse

Definition at line 121 of file VarianceTemplate.hpp.

Referenced by accumulate(), and clear().

std::vector< size_t > MBMesquite::VarianceTemplate::qmHandles [mutable, private]

Temporary storage for qm sample handles

Definition at line 128 of file VarianceTemplate.hpp.

Referenced by evaluate(), evaluate_with_gradient(), and evaluate_with_Hessian_diagonal().

Saved count from previous patch

Definition at line 123 of file VarianceTemplate.hpp.

Referenced by accumulate(), and clear().

Saved sum from previous patch

Definition at line 125 of file VarianceTemplate.hpp.

Referenced by accumulate(), and clear().

Saved sum from previous patch

Definition at line 124 of file VarianceTemplate.hpp.

Referenced by accumulate(), and clear().

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