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

Calculates the L_p objective function raised to the pth power. That is, sums the p_th powers of (the absolute value of) the quality metric values. More...

#include <LPtoPTemplate.hpp>

+ Inheritance diagram for MBMesquite::LPtoPTemplate:
+ Collaboration diagram for MBMesquite::LPtoPTemplate:

Public Member Functions

MESQUITE_EXPORT LPtoPTemplate (QualityMetric *, short, MsqError &)
MESQUITE_EXPORT LPtoPTemplate (short, QualityMetric *)
virtual MESQUITE_EXPORT ~LPtoPTemplate ()
virtual MESQUITE_EXPORT void clear ()
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 bool evaluate_with_Hessian (EvalType type, PatchData &pd, double &value_out, std::vector< Vector3D > &grad_out, MsqHessian &Hessian_out, MsqError &err)
 Evaluate objective function and Hessian for specified patch.
virtual MESQUITE_EXPORT
ObjectiveFunction
clone () const
 Create copy with same state.
MESQUITE_EXPORT void set_dividing_by_n (bool d_bool)

Private Member Functions

double get_value (double power_sum, size_t count, EvalType type, size_t &global_count, MsqError &err)

Private Attributes

short pVal
 The metric value entries are raised to the pVal power.
bool dividingByN
size_t mCount
double mPowSum
size_t saveCount
double savePowSum
std::vector< size_t > qmHandles
std::vector< size_t > mIndices
std::vector< Vector3DmGradient
std::vector< SymMatrix3DmDiag
std::vector< Matrix3DmHessian

Detailed Description

Calculates the L_p objective function raised to the pth power. That is, sums the p_th powers of (the absolute value of) the quality metric values.

Todo:
MB. Suggestions made by Todd Munson: a) There is an inconsistent use of fabs. The hessian evaluation when using the one norm does not take the absolute value, while the gradient does. b) The analytic gradient and hessian evaluations are incorrect when the quality metric changes sign due to taking the absolute value. The negative of the element gradient and hessian also needs to be taken. c) Done. The analytic gradient and hessian evaluations are incorrect when the negate flag is set to -1. The negative of the element gradient and hessian also needs to be taken in this case. d) The malloc in the concrete_eval routine should be removed.

Definition at line 67 of file LPtoPTemplate.hpp.


Constructor & Destructor Documentation

LPtoPTemplate::LPtoPTemplate ( QualityMetric qualitymetric,
short  Pinput,
MsqError err 
)

Definition at line 47 of file LPtoPTemplate.cpp.

References clear(), dividingByN, MBMesquite::MsqError::INVALID_ARG, MSQ_SETERR, and pVal.

Referenced by clone().

    : ObjectiveFunctionTemplate( qualitymetric )
{
    pVal = Pinput;
    if( pVal < 1 )
    {
        MSQ_SETERR( err )( "P_VALUE must be greater than 0.", MsqError::INVALID_ARG );
        return;
    }

    dividingByN = false;

    clear();
}

Definition at line 62 of file LPtoPTemplate.cpp.

References clear().

    : ObjectiveFunctionTemplate( qm ), pVal( P ), dividingByN( false )
{
    clear();
}

Definition at line 77 of file LPtoPTemplate.cpp.

{}

Member Function Documentation

void LPtoPTemplate::clear ( ) [virtual]

Clear any values accumulated for BCD-related eval calls

Implements MBMesquite::ObjectiveFunction.

Definition at line 68 of file LPtoPTemplate.cpp.

References mCount, mPowSum, saveCount, and savePowSum.

Referenced by LPtoPTemplate().

{
    mCount     = 0;
    mPowSum    = 0;
    saveCount  = 0;
    savePowSum = 0;
}
ObjectiveFunction * LPtoPTemplate::clone ( ) const [virtual]

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.

Definition at line 79 of file LPtoPTemplate.cpp.

References LPtoPTemplate().

{
    return new LPtoPTemplate( *this );
}
bool LPtoPTemplate::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.

Definition at line 139 of file LPtoPTemplate.cpp.

References MBMesquite::ObjectiveFunction::ACCUMULATE, MBMesquite::QualityMetric::evaluate(), MBMesquite::QualityMetric::get_evaluations(), MBMesquite::QualityMetric::get_negate_flag(), MBMesquite::ObjectiveFunctionTemplate::get_quality_metric(), MBMesquite::QualityMetric::get_single_pass(), get_value(), MSQ_CHKERR, MSQ_ERRFALSE, pVal, qmHandles, 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, working_sum = 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;

        double tmp_val = value;
        for( short j = 1; j < pVal; ++j )
            tmp_val *= value;
        working_sum += fabs( tmp_val );
    }

    // get overall OF value, update member data, etc.
    size_t global_count;
    value_out = qm->get_negate_flag() * get_value( working_sum, qmHandles.size(), type, global_count, err );
    //  if (!global_count)
    //    return false;  // invalid mesh
    //  else
    return true;
}
bool LPtoPTemplate::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.

Definition at line 171 of file LPtoPTemplate.cpp.

References dividingByN, MBMesquite::QualityMetric::evaluate_with_gradient(), MBMesquite::QualityMetric::get_evaluations(), MBMesquite::QualityMetric::get_negate_flag(), MBMesquite::ObjectiveFunctionTemplate::get_quality_metric(), get_value(), mGradient, mIndices, MSQ_CHKERR, MSQ_ERRFALSE, MBMesquite::PatchData::num_free_vertices(), MBMesquite::OF_FREE_EVALS_ONLY, pVal, and qmHandles.

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

    // zero gradient
    grad_out.clear();
    grad_out.resize( pd.num_free_vertices(), Vector3D( 0.0, 0.0, 0.0 ) );
    bool qm_bool = true;
    double QM_val;
    OF_val = 0.;
    int p1;

    // calculate OF value and gradient for just the patch
    std::vector< size_t >::const_iterator i;
    for( i = qmHandles.begin(); i != qmHandles.end(); ++i )
    {
        qm_bool = qm->evaluate_with_gradient( pd, *i, QM_val, mIndices, mGradient, err );
        if( MSQ_CHKERR( err ) || !qm_bool ) return false;

        QM_val        = fabs( QM_val );
        double QM_pow = 1.0;
        double factor = qm->get_negate_flag();
        if( pVal == 1 )
            QM_pow = 1.0;
        else
        {
            QM_pow = QM_val;
            for( p1 = 2; p1 < pVal; ++p1 )
                QM_pow *= QM_val;
            factor *= QM_pow * pVal;
        }

        OF_val += QM_pow * QM_val;
        for( size_t j = 0; j < mIndices.size(); ++j )
        {
            mGradient[j] *= factor;
            grad_out[mIndices[j]] += mGradient[j];
        }
    }

    // get overall OF value, update member data, etc.
    size_t global_count;
    OF_val = qm->get_negate_flag() * get_value( OF_val, qmHandles.size(), type, global_count, err );
    //  if (!global_count)
    //    return false;  // invalid mesh

    if( dividingByN && global_count )
    {
        const double inv_n = 1.0 / global_count;
        std::vector< Vector3D >::iterator g;
        for( g = grad_out.begin(); g != grad_out.end(); ++g )
            *g *= inv_n;
    }

    return true;
}
bool LPtoPTemplate::evaluate_with_Hessian ( EvalType  type,
PatchData pd,
double &  value_out,
std::vector< Vector3D > &  grad_out,
MsqHessian Hessian_out,
MsqError err 
) [virtual]

Evaluate objective function and 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 will fail.

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.
Hessian_outThe Hessian of the OF wrt the coordinates of each *free* vertex in the patch.
Returns:
false if any QualityMetric evaluation returned false, true otherwise.

Reimplemented from MBMesquite::ObjectiveFunction.

Definition at line 365 of file LPtoPTemplate.cpp.

References MBMesquite::MsqHessian::add(), dividingByN, MBMesquite::QualityMetric::evaluate_with_Hessian(), MBMesquite::QualityMetric::get_evaluations(), MBMesquite::QualityMetric::get_negate_flag(), MBMesquite::ObjectiveFunctionTemplate::get_quality_metric(), get_value(), MBMesquite::MsqError::INVALID_STATE, mGradient, mHessian, mIndices, MSQ_CHKERR, MSQ_ERRFALSE, MSQ_SETERR, n, MBMesquite::PatchData::num_free_vertices(), MBMesquite::OF_FREE_EVALS_ONLY, MBMesquite::Matrix3D::outer_product(), pVal, qmHandles, MBMesquite::MsqHessian::scale(), and MBMesquite::MsqHessian::zero_out().

Referenced by ObjectiveFunctionTest::test_compute_ana_hessian_tet(), and ObjectiveFunctionTest::test_compute_ana_hessian_tet_scaled().

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

    // zero gradient and hessian
    grad.clear();
    grad.resize( pd.num_free_vertices(), 0.0 );
    hessian.zero_out();

    double QM_val, QM_pow = 1.0;
    double fac1, fac2;
    Matrix3D elem_outer_product;
    bool qm_bool;
    size_t i, j, n;
    short p;

    // Loops over all elements in the patch.
    OF_val = 0.0;
    std::vector< size_t >::const_iterator k;
    for( k = qmHandles.begin(); k != qmHandles.end(); ++k )
    {
        // Computes \nabla^2 Q(e). Only the free vertices will have non-zero entries.
        qm_bool = qm->evaluate_with_Hessian( pd, *k, QM_val, mIndices, mGradient, mHessian, err );
        if( MSQ_CHKERR( err ) || !qm_bool ) return false;
        QM_val = fabs( QM_val );

        // **** Computes Hessian ****
        const size_t nve = mIndices.size();
        if( pVal == 1 )
        {
            QM_pow = 1.0;
            n      = 0;
            for( i = 0; i < nve; ++i )
            {
                for( j = i; j < nve; ++j )
                {
                    // negate if necessary
                    mHessian[n] *= negate_flag;
                    hessian.add( mIndices[i], mIndices[j], mHessian[n], err );
                    MSQ_ERRFALSE( err );
                    ++n;
                }
            }
            fac1 = 1;
        }
        else if( pVal >= 2 )
        {
            // Computes the coefficients:
            QM_pow = 1.0;
            for( p = 0; p < pVal - 2; ++p )
                QM_pow *= QM_val;
            // 1 - computes p(p-1)Q(e)^{p-2}
            fac2 = pVal * ( pVal - 1 ) * QM_pow;
            // 2 - computes  pQ(e)^{p-1}
            QM_pow *= QM_val;
            fac1 = pVal * QM_pow;

            // fac1 *= qm->get_negate_flag();
            // fac2 *= qm->get_negate_flag();

            n = 0;
            for( i = 0; i < nve; ++i )
            {
                for( j = i; j < nve; ++j )
                {
                    elem_outer_product.outer_product( mGradient[i], mGradient[j] );
                    elem_outer_product *= fac2;
                    mHessian[n] *= fac1;
                    mHessian[n] += elem_outer_product;
                    mHessian[n] *= negate_flag;
                    hessian.add( mIndices[i], mIndices[j], mHessian[n], err );
                    MSQ_ERRFALSE( err );
                    ++n;
                }
            }
        }
        else
        {
            MSQ_SETERR( err )( " invalid P value.", MsqError::INVALID_STATE );
            return false;
        }

        // **** Computes Gradient ****

        // For each vertex in the element ...
        for( i = 0; i < nve; ++i )
        {
            // ... computes p*q^{p-1}*grad(q) ...
            mGradient[i] *= fac1 * negate_flag;
            // ... and accumulates it in the objective function gradient.
            // also scale the gradient by the scaling factor
            assert( mIndices[i] < pd.num_free_vertices() );
            grad[mIndices[i]] += mGradient[i];
        }

        // **** computes Objective Function value \sum_{i=1}^{N_e} |q_i|^P ****
        OF_val += QM_pow * QM_val;
    }

    size_t global_count;
    OF_val = negate_flag * get_value( OF_val, qmHandles.size(), type, global_count, err );
    //  if (!global_count)
    //    return false;  // invalid mesh

    if( dividingByN && global_count )
    {
        const double inv_n = 1.0 / global_count;
        std::vector< Vector3D >::iterator g;
        for( g = grad.begin(); g != grad.end(); ++g )
            *g *= inv_n;
        hessian.scale( inv_n );
    }

    return true;
}
bool LPtoPTemplate::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.

Definition at line 234 of file LPtoPTemplate.cpp.

References dividingByN, MBMesquite::QualityMetric::evaluate_with_Hessian_diagonal(), MBMesquite::QualityMetric::get_evaluations(), MBMesquite::QualityMetric::get_negate_flag(), MBMesquite::ObjectiveFunctionTemplate::get_quality_metric(), get_value(), MBMesquite::MsqError::INVALID_STATE, mDiag, mGradient, mIndices, MSQ_CHKERR, MSQ_ERRFALSE, MSQ_SETERR, MBMesquite::PatchData::num_free_vertices(), MBMesquite::OF_FREE_EVALS_ONLY, pVal, and qmHandles.

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

    // zero gradient and hessian
    grad.clear();
    grad.resize( pd.num_free_vertices(), 0.0 );
    hess_diag.clear();
    hess_diag.resize( pd.num_free_vertices(), 0.0 );

    double QM_val, QM_pow = 1.0;
    double fac1, fac2;
    const double negate_flag = qm->get_negate_flag();
    bool qm_bool;
    size_t i;
    short p;

    // Loops over all elements in the patch.
    OF_val = 0.0;
    std::vector< size_t >::const_iterator k;
    for( k = qmHandles.begin(); k != qmHandles.end(); ++k )
    {
        // Computes \nabla^2 Q(e). Only the free vertices will have non-zero entries.
        qm_bool = qm->evaluate_with_Hessian_diagonal( pd, *k, QM_val, mIndices, mGradient, mDiag, err );
        if( MSQ_CHKERR( err ) || !qm_bool ) return false;
        QM_val = fabs( QM_val );

        // **** Computes Hessian ****
        const size_t nve = mIndices.size();
        if( pVal == 1 )
        {
            QM_pow = 1.0;
            for( i = 0; i < nve; ++i )
            {
                mDiag[i] *= negate_flag;
                hess_diag[mIndices[i]] += mDiag[i];
            }
            fac1 = 1;
        }
        else if( pVal >= 2 )
        {
            // Computes the coefficients:
            QM_pow = 1.0;
            for( p = 0; p < pVal - 2; ++p )
                QM_pow *= QM_val;
            // 1 - computes p(p-1)Q(e)^{p-2}
            fac2 = pVal * ( pVal - 1 ) * QM_pow;
            // 2 - computes  pQ(e)^{p-1}
            QM_pow *= QM_val;
            fac1 = pVal * QM_pow;

            // fac1 *= qm->get_negate_flag();
            // fac2 *= qm->get_negate_flag();

            for( i = 0; i < nve; ++i )
            {
                SymMatrix3D op( mGradient[i] );
                op *= fac2;
                mDiag[i] *= fac1;
                op += mDiag[i];
                op *= negate_flag;
                hess_diag[mIndices[i]] += op;
            }
        }
        else
        {
            MSQ_SETERR( err )( " invalid P value.", MsqError::INVALID_STATE );
            return false;
        }

        // **** Computes Gradient ****

        // For each vertex in the element ...
        for( i = 0; i < nve; ++i )
        {
            // ... computes p*q^{p-1}*grad(q) ...
            mGradient[i] *= fac1 * negate_flag;
            // ... and accumulates it in the objective function gradient.
            // also scale the gradient by the scaling factor
            assert( mIndices[i] < pd.num_free_vertices() );
            grad[mIndices[i]] += mGradient[i];
        }

        // **** computes Objective Function value \sum_{i=1}^{N_e} |q_i|^P ****
        OF_val += QM_pow * QM_val;
    }

    size_t global_count;
    OF_val = negate_flag * get_value( OF_val, qmHandles.size(), type, global_count, err );
    //  if (!global_count)
    //    return false;  // invalid mesh

    if( dividingByN && global_count )
    {
        const double inv_n = 1.0 / global_count;
        for( i = 0; i < pd.num_free_vertices(); ++i )
        {
            grad[i] *= inv_n;
            hess_diag[i] *= inv_n;
        }
    }

    return true;
}
double LPtoPTemplate::get_value ( double  power_sum,
size_t  count,
EvalType  type,
size_t &  global_count,
MsqError err 
) [private]

Definition at line 84 of file LPtoPTemplate.cpp.

References MBMesquite::ObjectiveFunction::ACCUMULATE, MBMesquite::ObjectiveFunction::CALCULATE, dividingByN, mCount, mPowSum, MBMesquite::ObjectiveFunction::SAVE, saveCount, savePowSum, MBMesquite::ObjectiveFunction::TEMPORARY, and MBMesquite::ObjectiveFunction::UPDATE.

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

{
    double result = 0;
    switch( type )
    {
        default:
        case CALCULATE:
            result       = power_sum;
            global_count = count;
            break;

        case ACCUMULATE:
            mPowSum += power_sum;
            mCount += count;
            result       = mPowSum;
            global_count = mCount;
            break;

        case SAVE:
            savePowSum   = power_sum;
            saveCount    = count;
            result       = mPowSum;
            global_count = mCount;
            break;

        case UPDATE:
            mPowSum -= savePowSum;
            mCount -= saveCount;
            savePowSum = power_sum;
            saveCount  = count;
            mPowSum += savePowSum;
            mCount += saveCount;
            result       = mPowSum;
            global_count = mCount;
            break;

        case TEMPORARY:
            result       = mPowSum - savePowSum + power_sum;
            global_count = mCount + count - saveCount;
            break;
    }

    //  if (!global_count)
    //    {
    //      MSQ_SETERR(err)(" global_count is zero, possibly due to an invalid mesh.",
    //      MsqError::INVALID_MESH); return -1;  // result is invalid
    //    }
    if( dividingByN && global_count ) result /= global_count;
    return result;
}

Use set_dividing_by_n to control whether this objective function divides it's final value by the number of metric values used to compute the objective function value. That is, if the associated metric is element based, the obejctive function value is divided by the number of elements. If it is vertex based, the objective function is divided by the number of vertices. If this function is passed 'true', the function value will be scale. If it is passed false, the function value will not be scaled.

Definition at line 121 of file LPtoPTemplate.hpp.

References dividingByN.

Referenced by ObjectiveFunctionTest::test_compute_ana_hessian_tet_scaled(), ObjectiveFunctionTest::test_compute_gradient_LPtoPTemplate_L2_scaled(), ObjectiveFunctionTest::test_diagonal_gradient_LPtoPTemplate_L2_scaled(), ObjectiveFunctionTest::test_hessian_diagonal_LPtoPTemplate_L2_scaled(), and ObjectiveFunctionTest::test_hessian_gradient_LPtoPTemplate_L2_scaled().

    {
        dividingByN = d_bool;
    }

Member Data Documentation

dividingByN is true if we are dividing the objective function by the number of metric values.

Definition at line 133 of file LPtoPTemplate.hpp.

Referenced by evaluate_with_gradient(), evaluate_with_Hessian(), evaluate_with_Hessian_diagonal(), get_value(), LPtoPTemplate(), and set_dividing_by_n().

The number of accumulated entires

Definition at line 135 of file LPtoPTemplate.hpp.

Referenced by clear(), and get_value().

std::vector< SymMatrix3D > MBMesquite::LPtoPTemplate::mDiag [mutable, private]

Temporary storage for qm Hessian diagonal blocks

Definition at line 147 of file LPtoPTemplate.hpp.

Referenced by evaluate_with_Hessian_diagonal().

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

Temporary storage for qm gradient

Definition at line 145 of file LPtoPTemplate.hpp.

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

std::vector< Matrix3D > MBMesquite::LPtoPTemplate::mHessian [mutable, private]

Temporary storage for qm Hessian

Definition at line 149 of file LPtoPTemplate.hpp.

Referenced by evaluate_with_Hessian().

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

Temporary storage for qm vertex indices

Definition at line 143 of file LPtoPTemplate.hpp.

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

The accumulated sum of values

Definition at line 136 of file LPtoPTemplate.hpp.

Referenced by clear(), and get_value().

The metric value entries are raised to the pVal power.

Definition at line 130 of file LPtoPTemplate.hpp.

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

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

Temporary storage for qm sample handles

Definition at line 141 of file LPtoPTemplate.hpp.

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

Saved count from previous patch

Definition at line 137 of file LPtoPTemplate.hpp.

Referenced by clear(), and get_value().

Saved sum from previous patch

Definition at line 138 of file LPtoPTemplate.hpp.

Referenced by clear(), and get_value().

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