MOAB: Mesh Oriented datABase
(version 5.4.1)
|
Base class for concrete quality metrics. More...
#include <QualityMetric.hpp>
Public Types | |
enum | MetricType { VERTEX_BASED, ELEMENT_BASED } |
enum | AveragingMethod { LINEAR, RMS, HMS, SUM, SUM_SQUARED, HARMONIC, LAST_WITH_HESSIAN = HARMONIC, MINIMUM, MAXIMUM, GEOMETRIC, LAST_WITH_GRADIENT = GEOMETRIC, STANDARD_DEVIATION, MAX_OVER_MIN, MAX_MINUS_MIN, SUM_OF_RATIOS_SQUARED } |
Public Member Functions | |
virtual MESQUITE_EXPORT | ~QualityMetric () |
virtual MESQUITE_EXPORT MetricType | get_metric_type () const =0 |
virtual MESQUITE_EXPORT std::string | get_name () const =0 |
virtual MESQUITE_EXPORT int | get_negate_flag () const =0 |
1 if metric should be minimized, -1 if metric should be maximized. | |
virtual MESQUITE_EXPORT void | get_evaluations (PatchData &pd, std::vector< size_t > &handles, bool free_vertices_only, MsqError &err)=0 |
Get locations at which metric can be evaluated. | |
virtual MESQUITE_EXPORT void | get_single_pass (PatchData &pd, std::vector< size_t > &handles, bool free_vertices_only, MsqError &err) |
Get locations at which metric can be evaluated for use in BCD intialization and QualityAssessor. | |
virtual MESQUITE_EXPORT bool | evaluate (PatchData &pd, size_t handle, double &value, MsqError &err)=0 |
Get metric value at a logical location in the patch. | |
virtual MESQUITE_EXPORT bool | evaluate_with_indices (PatchData &pd, size_t handle, double &value, std::vector< size_t > &indices, MsqError &err)=0 |
Get metric value at a logical location in the patch. | |
virtual MESQUITE_EXPORT bool | evaluate_with_gradient (PatchData &pd, size_t handle, double &value, std::vector< size_t > &indices, std::vector< Vector3D > &gradient, MsqError &err) |
Get metric value and gradient at a logical location in the patch. | |
virtual MESQUITE_EXPORT bool | evaluate_with_Hessian_diagonal (PatchData &pd, size_t handle, double &value, std::vector< size_t > &indices, std::vector< Vector3D > &gradient, std::vector< SymMatrix3D > &Hessian_diagonal, MsqError &err) |
Get metric value and gradient at a logical location in the patch. | |
virtual MESQUITE_EXPORT bool | evaluate_with_Hessian (PatchData &pd, size_t handle, double &value, std::vector< size_t > &indices, std::vector< Vector3D > &gradient, std::vector< Matrix3D > &Hessian, MsqError &err) |
Get metric value and deravitives at a logical location in the patch. | |
MESQUITE_EXPORT double | weighted_average_metrics (const double coef[], const double metric_values[], const int &num_values, MsqError &err) |
virtual MESQUITE_EXPORT void | initialize_queue (MeshDomainAssoc *mesh_and_domain, const Settings *settings, MsqError &err) |
Called at start of instruction queue processing. | |
Static Public Member Functions | |
static double | vertex_barrier_function (double det, double delta) |
Escobar Barrier Function for Shape and Other Metrics. | |
static MESQUITE_EXPORT void | remove_fixed_gradients (EntityTopology type, uint32_t fixed_vertices, std::vector< Vector3D > &gradients) |
Remove from vector any gradient terms corresponding to a fixed vertex. | |
static MESQUITE_EXPORT void | remove_fixed_diagonals (EntityTopology type, uint32_t fixed_vertices, std::vector< Vector3D > &gradients, std::vector< SymMatrix3D > &hess_diagonal_blocks) |
Remove from vectors any gradient terms and hessian diagonal blcoks corresponding to a fixed vertex. | |
static MESQUITE_EXPORT void | remove_fixed_hessians (EntityTopology type, uint32_t fixed_vertices, std::vector< Matrix3D > &hessians) |
Remove from vector any Hessian blocks corresponding to a fixed vertex. | |
static MESQUITE_EXPORT uint32_t | fixed_vertex_bitmap (PatchData &pd, const MsqMeshEntity *elem, std::vector< size_t > &free_indices) |
Convert fixed vertex format from list to bit flags. | |
Protected Member Functions | |
QualityMetric () | |
Private Attributes | |
int | feasible |
std::vector< Matrix3D > | tmpHess |
bool | keepFiniteDiffEps |
bool | haveFiniteDiffEps |
double | finiteDiffEps |
Base class for concrete quality metrics.
Definition at line 71 of file QualityMetric.hpp.
AveragingMethod allows you to set how the quality metric values attained at each sample point will be averaged together to produce a single metric value for an element.
Definition at line 303 of file QualityMetric.hpp.
{ LINEAR, //!< the linear average RMS, //!< the root-mean-squared average HMS, //!< the harmonic-mean-squared average SUM, //!< the sum of the values SUM_SQUARED, //!< the sum of the squares of the values HARMONIC, //!< the harmonic average LAST_WITH_HESSIAN = HARMONIC, MINIMUM, //!< the minimum value MAXIMUM, //!< the maximum value GEOMETRIC, //!< the geometric average LAST_WITH_GRADIENT = GEOMETRIC, STANDARD_DEVIATION, //!< the standard deviation squared of the values MAX_OVER_MIN, //!< the maximum value minus the minum value MAX_MINUS_MIN, //!< the maximum value divided by the minimum value SUM_OF_RATIOS_SQUARED //!< (1/(N^2))*(SUM (SUM (v_i/v_j)^2)) };
VERTEX_BASED |
Iterate over vertices to evaluate metric. |
ELEMENT_BASED |
Iterate over elements to evaluate metric. |
Reimplemented in MBMesquite::AWQualityMetric, and MBMesquite::TQualityMetric.
Definition at line 77 of file QualityMetric.hpp.
{ VERTEX_BASED, /**< Iterate over vertices to evaluate metric. */ ELEMENT_BASED /**< Iterate over elements to evaluate metric. */ };
MBMesquite::QualityMetric::QualityMetric | ( | ) | [inline, protected] |
Definition at line 74 of file QualityMetric.hpp.
: keepFiniteDiffEps( false ), haveFiniteDiffEps( false ) {}
virtual MESQUITE_EXPORT MBMesquite::QualityMetric::~QualityMetric | ( | ) | [inline, virtual] |
Definition at line 83 of file QualityMetric.hpp.
{}
virtual MESQUITE_EXPORT bool MBMesquite::QualityMetric::evaluate | ( | PatchData & | pd, |
size_t | handle, | ||
double & | value, | ||
MsqError & | err | ||
) | [pure virtual] |
Get metric value at a logical location in the patch.
Evaluate the metric at one location in the PatchData.
pd | The patch. |
handle | The location in the patch (as passed back from get_evaluations). |
value | The output metric value. |
Implemented in OFTestBadQM, DistTestMetric, MetricLogger, TriTauMetric, ParabolicVertexMetric, OFTestQM, CompareMetric, ConstantElementMetric, LinearVertexMetric, NumericQM, FauxMetric< B >, MBMesquite::MultiplyQualityMetric, MBMesquite::CompareQM, MBMesquite::PowerQualityMetric, MBMesquite::IdealWeightInverseMeanRatio, MBMesquite::UntangleBetaQualityMetric, MBMesquite::TMPQualityMetric, MBMesquite::IdealWeightMeanRatio, MBMesquite::AffineMapMetric, MBMesquite::EdgeLengthRangeQualityMetric, MBMesquite::ScalarMultiplyQualityMetric, MBMesquite::ConditionNumberQualityMetric, MBMesquite::AddQualityMetric, MBMesquite::LocalSizeQualityMetric, MBMesquite::EdgeLengthQualityMetric, MBMesquite::NumericalQM, MBMesquite::VertexConditionNumberQualityMetric, MBMesquite::ScalarAddQualityMetric, MBMesquite::AspectRatioGammaQualityMetric, MBMesquite::ElementPMeanP, MBMesquite::VertexPMeanP, MBMesquite::ElementMaxQM, MBMesquite::VertexMaxQM, MBMesquite::ElementAvgQM, MBMesquite::SizeMetric, and MBMesquite::EdgeLengthMetric.
Referenced by MBMesquite::PMeanPMetric::average(), QualityMetricTester::compare_eval_and_eval_with_indices(), MBMesquite::NonSmoothDescent::compute_function(), MBMesquite::LInfTemplate::evaluate(), MBMesquite::ElementAvgQM::evaluate(), MBMesquite::MaxTemplate::evaluate(), MBMesquite::ElementMaxQM::evaluate(), MBMesquite::VertexMaxQM::evaluate(), MBMesquite::ScalarAddQualityMetric::evaluate(), MBMesquite::NumericalQM::evaluate(), MBMesquite::AddQualityMetric::evaluate(), MBMesquite::ScalarMultiplyQualityMetric::evaluate(), MBMesquite::PatchPowerMeanP::evaluate(), MBMesquite::VarianceTemplate::evaluate(), MBMesquite::LPtoPTemplate::evaluate(), MBMesquite::CompareQM::evaluate(), MBMesquite::PowerQualityMetric::evaluate(), MBMesquite::MultiplyQualityMetric::evaluate(), MBMesquite::PMeanPTemplate::evaluate(), evaluate_with_gradient(), MBMesquite::ElementQM::evaluate_with_indices(), MBMesquite::EdgeQM::evaluate_with_indices(), MBMesquite::MetricWeight::get_weight(), MBMesquite::InverseMetricWeight::get_weight(), test_bad_element(), QualityMetricTester::test_domain_deviation_quality(), test_evaluate(), QualityMetricTester::test_measures_quality(), QualityMetricTester::test_measures_transform(), QualityMetricTester::test_measures_vertex_quality(), QualityMetricTester::test_transform_invariant(), QualityMetricTester::test_type_is_not_supported(), and QualityMetricTester::test_type_is_supported().
bool MBMesquite::QualityMetric::evaluate_with_gradient | ( | PatchData & | pd, |
size_t | handle, | ||
double & | value, | ||
std::vector< size_t > & | indices, | ||
std::vector< Vector3D > & | gradient, | ||
MsqError & | err | ||
) | [virtual] |
Get metric value and gradient at a logical location in the patch.
Evaluate the metric at one location in the PatchData.
pd | The patch. |
handle | The location in the patch (as passed back from get_evaluations). |
value | The output metric value. |
indices | The free vertices that the evaluation is a function of, specified as vertex indices in the PatchData. |
gradient | The gradient of the metric as a function of the coordinates of the free vertices passed back in the indices list. |
Reimplemented in CompareMetric, FauxMetric< B >, MBMesquite::MultiplyQualityMetric, MBMesquite::CompareQM, MBMesquite::PowerQualityMetric, MBMesquite::IdealWeightInverseMeanRatio, MBMesquite::IdealWeightMeanRatio, MBMesquite::AWQualityMetric, MBMesquite::TQualityMetric, MBMesquite::ScalarMultiplyQualityMetric, MBMesquite::AddQualityMetric, MBMesquite::NumericalQM, MBMesquite::ScalarAddQualityMetric, MBMesquite::VertexPMeanP, MBMesquite::ElementPMeanP, and MBMesquite::EdgeLengthMetric.
Definition at line 98 of file QualityMetric.cpp.
References evaluate(), evaluate_with_indices(), finiteDiffEps, MBMesquite::get_delta_C(), haveFiniteDiffEps, INTERNAL_ERROR, keepFiniteDiffEps, MSQ_CHKERR, MSQ_ERRZERO, MSQ_SETERR, MBMesquite::PatchData::set_vertex_coordinates(), value(), and MBMesquite::PatchData::vertex_by_index().
Referenced by MBMesquite::PMeanPMetric::average_with_gradient(), QualityMetricTester::compare_analytical_and_numerical_gradients(), TMPQualityMetricTest< QMType >::compare_analytical_and_numerical_gradients(), QualityMetricTester::compare_eval_with_grad_and_eval_with_diagonal(), QualityMetricTester::compare_eval_with_grad_and_eval_with_hessian(), QualityMetricTester::compare_eval_with_indices_and_eval_with_gradient(), QualityMetricTest::compare_gradient(), MBMesquite::NonSmoothDescent::compute_gradient(), MBMesquite::ScalarAddQualityMetric::evaluate_with_gradient(), MBMesquite::NumericalQM::evaluate_with_gradient(), MBMesquite::PatchPowerMeanP::evaluate_with_gradient(), MBMesquite::VarianceTemplate::evaluate_with_gradient(), MBMesquite::AddQualityMetric::evaluate_with_gradient(), MBMesquite::ScalarMultiplyQualityMetric::evaluate_with_gradient(), MBMesquite::LPtoPTemplate::evaluate_with_gradient(), MBMesquite::PowerQualityMetric::evaluate_with_gradient(), MBMesquite::CompareQM::evaluate_with_gradient(), MBMesquite::PMeanPTemplate::evaluate_with_gradient(), MBMesquite::MultiplyQualityMetric::evaluate_with_gradient(), evaluate_with_Hessian(), QualityMetricTester::test_domain_deviation_gradient(), QualityMetricTester::test_grad_transform_invariant(), QualityMetricTest::test_gradient_constant(), QualityMetricTest::test_gradient_linear(), QualityMetricTest::test_gradient_parabolic(), QualityMetricTester::test_gradient_reflects_quality(), QualityMetricTest::test_gradient_tau(), QualityMetricTester::test_gradient_with_fixed_vertex(), QualityMetricTester::test_ideal_element_zero_gradient(), QualityMetricTester::test_ideal_element_zero_vertex_gradient(), QualityMetricTester::test_type_is_not_supported(), QualityMetricTester::test_type_is_supported(), and QualityMetricTester::test_vertex_gradient_reflects_quality().
{ indices.clear(); bool valid = evaluate_with_indices( pd, handle, value, indices, err ); if( MSQ_CHKERR( err ) || !valid ) return false; if( indices.empty() ) return true; // get initial pertubation amount double delta_C = finiteDiffEps; if( !haveFiniteDiffEps ) { delta_C = get_delta_C( pd, indices, err ); MSQ_ERRZERO( err ); if( keepFiniteDiffEps ) { finiteDiffEps = delta_C; haveFiniteDiffEps = true; } } const double delta_inv_C = 1.0 / delta_C; const int reduction_limit = 15; gradient.resize( indices.size() ); for( size_t v = 0; v < indices.size(); ++v ) { const Vector3D pos = pd.vertex_by_index( indices[v] ); /* gradient in the x, y, z direction */ for( int j = 0; j < 3; ++j ) { double delta = delta_C; double delta_inv = delta_inv_C; double metric_value; Vector3D delta_v( 0, 0, 0 ); // perturb the node and calculate gradient. The while loop is a // safety net to make sure the epsilon perturbation does not take // the element out of the feasible region. int counter = 0; for( ;; ) { // perturb the coordinates of the free vertex in the j direction // by delta delta_v[j] = delta; pd.set_vertex_coordinates( pos + delta_v, indices[v], err ); MSQ_ERRZERO( err ); // compute the function at the perturbed point location valid = evaluate( pd, handle, metric_value, err ); if( !MSQ_CHKERR( err ) && valid ) break; if( ++counter >= reduction_limit ) { MSQ_SETERR( err ) ( "Perturbing vertex by delta caused an inverted element.", MsqError::INTERNAL_ERROR ); return false; } delta *= 0.1; delta_inv *= 10.; } // put the coordinates back where they belong pd.set_vertex_coordinates( pos, indices[v], err ); // compute the numerical gradient gradient[v][j] = ( metric_value - value ) * delta_inv; } // for(j) } // for(indices) return true; }
bool MBMesquite::QualityMetric::evaluate_with_Hessian | ( | PatchData & | pd, |
size_t | handle, | ||
double & | value, | ||
std::vector< size_t > & | indices, | ||
std::vector< Vector3D > & | gradient, | ||
std::vector< Matrix3D > & | Hessian, | ||
MsqError & | err | ||
) | [virtual] |
Get metric value and deravitives at a logical location in the patch.
Evaluate the metric at one location in the PatchData.
pd | The patch. |
handle | The location in the patch (as passed back from get_evaluations). |
value | The output metric value. |
indices | The free vertices that the evaluation is a function of, specified as vertex indices in the PatchData. |
gradient | The gradient of the metric as a function of the coordinates of the free vertices passed back in the indices list. |
Hessian | The Hessian of the metric as a function of the coordinates. The Hessian is passed back as the upper-triangular portion of the matrix in row-major order, where each Matrix3D is the portion of the Hessian with respect to the vertices at the corresponding positions in the indices list. |
Reimplemented in CompareMetric, FauxMetric< B >, MBMesquite::CompareQM, MBMesquite::MultiplyQualityMetric, MBMesquite::IdealWeightInverseMeanRatio, MBMesquite::PowerQualityMetric, MBMesquite::IdealWeightMeanRatio, MBMesquite::AWQualityMetric, MBMesquite::TQualityMetric, MBMesquite::NumericalQM, MBMesquite::ScalarMultiplyQualityMetric, MBMesquite::AddQualityMetric, MBMesquite::ScalarAddQualityMetric, MBMesquite::VertexPMeanP, and MBMesquite::ElementPMeanP.
Definition at line 173 of file QualityMetric.cpp.
References evaluate_with_gradient(), finiteDiffEps, MBMesquite::get_delta_C(), haveFiniteDiffEps, INTERNAL_ERROR, keepFiniteDiffEps, MSQ_CHKERR, MSQ_ERRZERO, MSQ_SETERR, MBMesquite::PatchData::set_vertex_coordinates(), and MBMesquite::PatchData::vertex_by_index().
Referenced by MBMesquite::PMeanPMetric::average_with_Hessian(), PMeanPTemplateTest::check_result(), QualityMetricTester::compare_analytical_and_numerical_hessians(), TMPQualityMetricTest< QMType >::compare_analytical_and_numerical_hessians(), QualityMetricTester::compare_eval_with_diag_and_eval_with_hessian(), QualityMetricTester::compare_eval_with_grad_and_eval_with_hessian(), QualityMetricTester::compare_eval_with_indices_and_eval_with_hessian(), MBMesquite::ScalarAddQualityMetric::evaluate_with_Hessian(), MBMesquite::PatchPowerMeanP::evaluate_with_Hessian(), MBMesquite::AddQualityMetric::evaluate_with_Hessian(), MBMesquite::ScalarMultiplyQualityMetric::evaluate_with_Hessian(), MBMesquite::NumericalQM::evaluate_with_Hessian(), MBMesquite::LPtoPTemplate::evaluate_with_Hessian(), MBMesquite::PowerQualityMetric::evaluate_with_Hessian(), MBMesquite::MultiplyQualityMetric::evaluate_with_Hessian(), MBMesquite::CompareQM::evaluate_with_Hessian(), MBMesquite::PMeanPTemplate::evaluate_with_Hessian(), evaluate_with_Hessian_diagonal(), QualityMetricTest::test_Hessian_constant(), QualityMetricTest::test_Hessian_linear(), QualityMetricTest::test_Hessian_parabolic(), QualityMetricTest::test_Hessian_tau(), QualityMetricTester::test_hessian_transform_invariant(), QualityMetricTester::test_hessian_with_fixed_vertex(), QualityMetricTester::test_ideal_element_positive_definite_Hessian(), QualityMetricTester::test_symmetric_Hessian_diagonal_blocks(), QualityMetricTester::test_type_is_not_supported(), and QualityMetricTester::test_type_is_supported().
{ indices.clear(); gradient.clear(); keepFiniteDiffEps = true; bool valid = evaluate_with_gradient( pd, handle, value, indices, gradient, err ); keepFiniteDiffEps = false; if( MSQ_CHKERR( err ) || !valid ) { haveFiniteDiffEps = false; return false; } if( indices.empty() ) { haveFiniteDiffEps = false; return true; } // get initial pertubation amount double delta_C; if( haveFiniteDiffEps ) { delta_C = finiteDiffEps; } else { delta_C = get_delta_C( pd, indices, err ); MSQ_ERRZERO( err ); } assert( delta_C < 1e30 ); const double delta_inv_C = 1.0 / delta_C; const int reduction_limit = 15; std::vector< Vector3D > temp_gradient( indices.size() ); const int num_hess = indices.size() * ( indices.size() + 1 ) / 2; Hessian.resize( num_hess ); for( unsigned v = 0; v < indices.size(); ++v ) { const Vector3D pos = pd.vertex_by_index( indices[v] ); for( int j = 0; j < 3; ++j ) { // x, y, and z double delta = delta_C; double delta_inv = delta_inv_C; double metric_value; Vector3D delta_v( 0, 0, 0 ); // find finite difference for gradient int counter = 0; for( ;; ) { delta_v[j] = delta; pd.set_vertex_coordinates( pos + delta_v, indices[v], err ); MSQ_ERRZERO( err ); valid = evaluate_with_gradient( pd, handle, metric_value, indices, temp_gradient, err ); if( !MSQ_CHKERR( err ) && valid ) break; if( ++counter >= reduction_limit ) { MSQ_SETERR( err ) ( "Algorithm did not successfully compute element's " "Hessian.\n", MsqError::INTERNAL_ERROR ); haveFiniteDiffEps = false; return false; } delta *= 0.1; delta_inv *= 10.0; } pd.set_vertex_coordinates( pos, indices[v], err ); MSQ_ERRZERO( err ); // compute the numerical Hessian for( unsigned w = 0; w <= v; ++w ) { // finite difference to get some entries of the Hessian Vector3D fd( temp_gradient[w] ); fd -= gradient[w]; fd *= delta_inv; // For the block at position w,v in a matrix, we need the corresponding index // (mat_index) in a 1D array containing only upper triangular blocks. unsigned sum_w = w * ( w + 1 ) / 2; // 1+2+3+...+w unsigned mat_index = w * indices.size() + v - sum_w; Hessian[mat_index][0][j] = fd[0]; Hessian[mat_index][1][j] = fd[1]; Hessian[mat_index][2][j] = fd[2]; } } // for(j) } // for(indices) haveFiniteDiffEps = false; return true; }
bool MBMesquite::QualityMetric::evaluate_with_Hessian_diagonal | ( | PatchData & | pd, |
size_t | handle, | ||
double & | value, | ||
std::vector< size_t > & | indices, | ||
std::vector< Vector3D > & | gradient, | ||
std::vector< SymMatrix3D > & | Hessian_diagonal, | ||
MsqError & | err | ||
) | [virtual] |
Get metric value and gradient at a logical location in the patch.
Evaluate the metric at one location in the PatchData.
pd | The patch. |
handle | The location in the patch (as passed back from get_evaluations). |
value | The output metric value. |
indices | The free vertices that the evaluation is a function of, specified as vertex indices in the PatchData. |
gradient | The gradient of the metric as a function of the coordinates of the free vertices passed back in the indices list. |
Hessian_diagonal | The 3x3 blocks along the diagonal of the Hessian matrix. |
Reimplemented in CompareMetric, MBMesquite::CompareQM, MBMesquite::IdealWeightInverseMeanRatio, MBMesquite::IdealWeightMeanRatio, MBMesquite::AWQualityMetric, MBMesquite::TQualityMetric, MBMesquite::VertexPMeanP, MBMesquite::NumericalQM, and MBMesquite::ElementPMeanP.
Definition at line 273 of file QualityMetric.cpp.
References evaluate_with_Hessian(), MSQ_CHKERR, and tmpHess.
Referenced by MBMesquite::PMeanPMetric::average_with_Hessian_diagonal(), QualityMetricTester::compare_analytical_and_numerical_diagonals(), TMPQualityMetricTest< QMType >::compare_analytical_and_numerical_diagonals(), QualityMetricTester::compare_eval_with_diag_and_eval_with_hessian(), QualityMetricTester::compare_eval_with_grad_and_eval_with_diagonal(), QualityMetricTester::compare_eval_with_indices_and_eval_with_diagonal(), MBMesquite::NumericalQM::evaluate_with_Hessian_diagonal(), MBMesquite::VarianceTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::LPtoPTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::CompareQM::evaluate_with_Hessian_diagonal(), MBMesquite::PMeanPTemplate::evaluate_with_Hessian_diagonal(), QualityMetricTest::test_diagonal_constant(), QualityMetricTest::test_diagonal_linear(), QualityMetricTest::test_diagonal_parabolic(), QualityMetricTest::test_diagonal_tau(), and QualityMetricTester::test_diagonal_with_fixed_vertex().
{ bool rval = evaluate_with_Hessian( pd, handle, value, indices, gradient, tmpHess, err ); if( MSQ_CHKERR( err ) || !rval ) return rval; size_t s = indices.size(); Hessian_diagonal.resize( s ); std::vector< Matrix3D >::const_iterator h = tmpHess.begin(); for( size_t i = 0; i < indices.size(); ++i ) { Hessian_diagonal[i] = h->upper(); h += s--; } return rval; }
virtual MESQUITE_EXPORT bool MBMesquite::QualityMetric::evaluate_with_indices | ( | PatchData & | pd, |
size_t | handle, | ||
double & | value, | ||
std::vector< size_t > & | indices, | ||
MsqError & | err | ||
) | [pure virtual] |
Get metric value at a logical location in the patch.
Evaluate the metric at one location in the PatchData.
pd | The patch. |
handle | The location in the patch (as passed back from get_evaluations). |
value | The output metric value. |
indices | The free vertices that the evaluation is a function of, specified as vertex indices in the PatchData. |
Implemented in TriTauMetric, ParabolicVertexMetric, ConstantElementMetric, CompareMetric, LinearVertexMetric, NumericQM, FauxMetric< B >, MBMesquite::MultiplyQualityMetric, MBMesquite::EdgeQM, MBMesquite::CompareQM, MBMesquite::PowerQualityMetric, MBMesquite::TMPQualityMetric, MBMesquite::AffineMapMetric, MBMesquite::EdgeLengthRangeQualityMetric, MBMesquite::ScalarMultiplyQualityMetric, MBMesquite::AddQualityMetric, MBMesquite::LocalSizeQualityMetric, MBMesquite::EdgeLengthQualityMetric, MBMesquite::NumericalQM, MBMesquite::ElementQM, MBMesquite::VertexConditionNumberQualityMetric, MBMesquite::ScalarAddQualityMetric, MBMesquite::VertexPMeanP, and MBMesquite::VertexMaxQM.
Referenced by MBMesquite::PMeanPMetric::average_with_indices(), QualityMetricTester::compare_eval_and_eval_with_indices(), QualityMetricTester::compare_eval_with_indices_and_eval_with_diagonal(), QualityMetricTester::compare_eval_with_indices_and_eval_with_gradient(), QualityMetricTester::compare_eval_with_indices_and_eval_with_hessian(), QualityMetricTest::compare_indices(), evaluate_with_gradient(), MBMesquite::VertexMaxQM::evaluate_with_indices(), MBMesquite::ScalarAddQualityMetric::evaluate_with_indices(), MBMesquite::NumericalQM::evaluate_with_indices(), MBMesquite::AddQualityMetric::evaluate_with_indices(), MBMesquite::ScalarMultiplyQualityMetric::evaluate_with_indices(), MBMesquite::CompareQM::evaluate_with_indices(), MBMesquite::PowerQualityMetric::evaluate_with_indices(), MBMesquite::MultiplyQualityMetric::evaluate_with_indices(), QualityMetricTester::test_domain_deviation_quality(), QualityMetricTester::test_get_element_indices(), QualityMetricTester::test_get_indices_fixed(), QualityMetricTester::test_get_sample_indices(), QualityMetricTester::test_get_vertex_indices(), QualityMetricTester::test_type_is_not_supported(), and QualityMetricTester::test_type_is_supported().
uint32_t MBMesquite::QualityMetric::fixed_vertex_bitmap | ( | PatchData & | pd, |
const MsqMeshEntity * | elem, | ||
std::vector< size_t > & | free_indices | ||
) | [static] |
Convert fixed vertex format from list to bit flags.
Given list of pointers to fixed vertices as passed to evaluation functions, convert to bit flag format used for many utility functions in this class. Bits correspond to vertices in the canonical vertex order, beginning with the least-significant bit. The bit is cleared for free vertices and set (1) for fixed vertices.
Definition at line 294 of file QualityMetric.cpp.
References MBMesquite::MsqMeshEntity::get_vertex_index_array(), MBMesquite::PatchData::num_free_vertices(), and MBMesquite::MsqMeshEntity::vertex_count().
Referenced by MBMesquite::IdealWeightMeanRatio::evaluate_with_gradient(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_gradient(), MBMesquite::IdealWeightMeanRatio::evaluate_with_Hessian(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_Hessian(), MBMesquite::IdealWeightMeanRatio::evaluate_with_Hessian_diagonal(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_Hessian_diagonal(), and QualityMetricTest::test_fixed_vertex_list().
{ indices.clear(); uint32_t result = ~(uint32_t)0; unsigned num_vtx = elem->vertex_count(); const size_t* vertices = elem->get_vertex_index_array(); indices.clear(); for( unsigned i = 0; i < num_vtx; ++i ) { if( vertices[i] < pd.num_free_vertices() ) { indices.push_back( vertices[i] ); result &= ~(uint32_t)( 1 << i ); } } return result; }
virtual MESQUITE_EXPORT void MBMesquite::QualityMetric::get_evaluations | ( | PatchData & | pd, |
std::vector< size_t > & | handles, | ||
bool | free_vertices_only, | ||
MsqError & | err | ||
) | [pure virtual] |
Get locations at which metric can be evaluated.
Different metrics are evaluated for different things within a patch. For example, an element-based metric will be evaluated once for each element in patch, a vertex-based metric once for each veretx, and a target/sample-point based metric will be evaluated once for each samle point in each element. This method returns a list of handles, one for each location in the patch at which the metric can be evaluated. The handle values are used as input to the evaluate methods.
pd | The patch |
handles | Output list of handles |
free_vertices_only | If true, only pass back evaluation points that depend on at least one free vertex. |
Implemented in CompareMetric, NumericQM, FauxMetric< B >, MBMesquite::MultiplyQualityMetric, MBMesquite::CompareQM, MBMesquite::PowerQualityMetric, MBMesquite::AffineMapMetric, MBMesquite::TMPQualityMetric, MBMesquite::ScalarMultiplyQualityMetric, MBMesquite::AddQualityMetric, MBMesquite::NumericalQM, MBMesquite::EdgeQM, MBMesquite::ScalarAddQualityMetric, MBMesquite::ElementQM, and MBMesquite::VertexQM.
Referenced by QualityMetricTester::compare_analytical_and_numerical_diagonals(), TMPQualityMetricTest< QMType >::compare_analytical_and_numerical_diagonals(), QualityMetricTester::compare_analytical_and_numerical_gradients(), TMPQualityMetricTest< QMType >::compare_analytical_and_numerical_gradients(), QualityMetricTester::compare_analytical_and_numerical_hessians(), TMPQualityMetricTest< QMType >::compare_analytical_and_numerical_hessians(), QualityMetricTester::compare_eval_and_eval_with_indices(), QualityMetricTester::compare_eval_with_diag_and_eval_with_hessian(), QualityMetricTester::compare_eval_with_grad_and_eval_with_diagonal(), QualityMetricTester::compare_eval_with_grad_and_eval_with_hessian(), QualityMetricTester::compare_eval_with_indices_and_eval_with_diagonal(), QualityMetricTester::compare_eval_with_indices_and_eval_with_gradient(), QualityMetricTester::compare_eval_with_indices_and_eval_with_hessian(), MBMesquite::LInfTemplate::evaluate(), MBMesquite::MaxTemplate::evaluate(), MBMesquite::PatchPowerMeanP::evaluate(), MBMesquite::VarianceTemplate::evaluate(), MBMesquite::LPtoPTemplate::evaluate(), MBMesquite::PMeanPTemplate::evaluate(), MBMesquite::PatchPowerMeanP::evaluate_with_gradient(), MBMesquite::VarianceTemplate::evaluate_with_gradient(), MBMesquite::LPtoPTemplate::evaluate_with_gradient(), MBMesquite::PMeanPTemplate::evaluate_with_gradient(), MBMesquite::PatchPowerMeanP::evaluate_with_Hessian(), MBMesquite::LPtoPTemplate::evaluate_with_Hessian(), MBMesquite::PMeanPTemplate::evaluate_with_Hessian(), MBMesquite::VarianceTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::LPtoPTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::PMeanPTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::ScalarAddQualityMetric::get_evaluations(), MBMesquite::NumericalQM::get_evaluations(), MBMesquite::AddQualityMetric::get_evaluations(), MBMesquite::ScalarMultiplyQualityMetric::get_evaluations(), MBMesquite::PowerQualityMetric::get_evaluations(), MBMesquite::CompareQM::get_evaluations(), MBMesquite::MultiplyQualityMetric::get_evaluations(), get_single_pass(), MBMesquite::NonSmoothDescent::init_opt(), MBMesquite::VertexMover::loop_over_mesh(), MBMesquite::QualityAssessor::loop_over_mesh_internal(), test_bad_element(), QualityMetricTester::test_diagonal_with_fixed_vertex(), QualityMetricTester::test_domain_deviation_gradient(), QualityMetricTester::test_domain_deviation_quality(), test_evaluate(), QualityMetricTester::test_get_element_evaluations(), QualityMetricTester::test_get_element_indices(), QualityMetricTester::test_get_in_element_evaluations(), QualityMetricTester::test_get_indices_fixed(), QualityMetricTester::test_get_sample_evaluations(), QualityMetricTester::test_get_sample_indices(), QualityMetricTester::test_get_vertex_evaluations(), QualityMetricTester::test_grad_transform_invariant(), QualityMetricTester::test_gradient_reflects_quality(), QualityMetricTester::test_gradient_with_fixed_vertex(), QualityMetricTester::test_hessian_transform_invariant(), QualityMetricTester::test_hessian_with_fixed_vertex(), QualityMetricTester::test_ideal_element_positive_definite_Hessian(), QualityMetricTester::test_ideal_element_zero_gradient(), QualityMetricTester::test_ideal_element_zero_vertex_gradient(), QualityMetricTester::test_measures_quality(), QualityMetricTester::test_measures_transform(), QualityMetricTester::test_symmetric_Hessian_diagonal_blocks(), QualityMetricTester::test_transform_invariant(), QualityMetricTester::test_type_is_not_supported(), and QualityMetricTester::test_type_is_supported().
virtual MESQUITE_EXPORT MetricType MBMesquite::QualityMetric::get_metric_type | ( | ) | const [pure virtual] |
Implemented in OFTestBadQM, MetricLogger, CompareMetric, OFTestQM, NumericQM, FauxMetric< B >, MBMesquite::CompareQM, MBMesquite::MultiplyQualityMetric, MBMesquite::PowerQualityMetric, MBMesquite::NumericalQM, MBMesquite::ScalarMultiplyQualityMetric, MBMesquite::ElemSampleQM, MBMesquite::ScalarAddQualityMetric, MBMesquite::AddQualityMetric, MBMesquite::EdgeQM, MBMesquite::ElementQM, and MBMesquite::VertexQM.
Referenced by MBMesquite::AddQualityMetric::AddQualityMetric(), MBMesquite::QualityAssessor::find_or_add(), MBMesquite::AddQualityMetric::get_metric_type(), MBMesquite::ScalarAddQualityMetric::get_metric_type(), MBMesquite::ScalarMultiplyQualityMetric::get_metric_type(), MBMesquite::NumericalQM::get_metric_type(), MBMesquite::PowerQualityMetric::get_metric_type(), MBMesquite::MultiplyQualityMetric::get_metric_type(), MBMesquite::CompareQM::get_metric_type(), MBMesquite::MultiplyQualityMetric::MultiplyQualityMetric(), QualityMetricTester::test_get_element_indices(), QualityMetricTester::test_get_indices_fixed(), QualityMetricTester::test_get_sample_indices(), and QualityMetricTester::test_get_vertex_indices().
virtual MESQUITE_EXPORT std::string MBMesquite::QualityMetric::get_name | ( | ) | const [pure virtual] |
Implemented in OFTestBadQM, DistTestMetric, MetricLogger, TriTauMetric, ParabolicVertexMetric, CompareMetric, ConstantElementMetric, OFTestQM, LinearVertexMetric, NumericQM, FauxMetric< B >, MBMesquite::IdealWeightInverseMeanRatio, MBMesquite::UntangleBetaQualityMetric, MBMesquite::CompareQM, MBMesquite::AWQualityMetric, MBMesquite::TQualityMetric, MBMesquite::MultiplyQualityMetric, MBMesquite::IdealWeightMeanRatio, MBMesquite::PowerQualityMetric, FauxMetric< B >, MBMesquite::EdgeLengthRangeQualityMetric, MBMesquite::ConditionNumberQualityMetric, MBMesquite::LocalSizeQualityMetric, MBMesquite::EdgeLengthQualityMetric, MBMesquite::AffineMapMetric, MBMesquite::ScalarMultiplyQualityMetric, MBMesquite::VertexConditionNumberQualityMetric, MBMesquite::NumericalQM, MBMesquite::AspectRatioGammaQualityMetric, MBMesquite::ElementPMeanP, MBMesquite::VertexPMeanP, MBMesquite::ElementMaxQM, MBMesquite::VertexMaxQM, MBMesquite::ScalarAddQualityMetric, MBMesquite::ElementAvgQM, MBMesquite::SizeMetric, MBMesquite::AddQualityMetric, and MBMesquite::EdgeLengthMetric.
Referenced by MBMesquite::AddQualityMetric::get_name(), MBMesquite::ElementAvgQM::get_name(), MBMesquite::ScalarAddQualityMetric::get_name(), MBMesquite::ElementMaxQM::get_name(), MBMesquite::VertexMaxQM::get_name(), MBMesquite::ElementPMeanP::get_name(), MBMesquite::VertexPMeanP::get_name(), MBMesquite::NumericalQM::get_name(), MBMesquite::ScalarMultiplyQualityMetric::get_name(), MBMesquite::PowerQualityMetric::get_name(), MBMesquite::MultiplyQualityMetric::get_name(), MBMesquite::CompareQM::get_name(), MBMesquite::CompareQM::print_stats(), and smooth_mesh().
virtual MESQUITE_EXPORT int MBMesquite::QualityMetric::get_negate_flag | ( | ) | const [pure virtual] |
1 if metric should be minimized, -1 if metric should be maximized.
Implemented in OFTestBadQM, DistTestMetric, MetricLogger, TriTauMetric, ParabolicVertexMetric, CompareMetric, ConstantElementMetric, OFTestQM, LinearVertexMetric, NumericQM, FauxMetric< B >, MBMesquite::IdealWeightInverseMeanRatio, MBMesquite::UntangleBetaQualityMetric, MBMesquite::CompareQM, MBMesquite::MultiplyQualityMetric, MBMesquite::IdealWeightMeanRatio, MBMesquite::PowerQualityMetric, MBMesquite::EdgeLengthRangeQualityMetric, MBMesquite::ConditionNumberQualityMetric, MBMesquite::LocalSizeQualityMetric, MBMesquite::EdgeLengthQualityMetric, MBMesquite::ScalarMultiplyQualityMetric, MBMesquite::VertexConditionNumberQualityMetric, MBMesquite::AffineMapMetric, MBMesquite::TMPQualityMetric, MBMesquite::NumericalQM, MBMesquite::AspectRatioGammaQualityMetric, MBMesquite::ElementPMeanP, MBMesquite::VertexPMeanP, MBMesquite::ElementMaxQM, MBMesquite::ScalarAddQualityMetric, MBMesquite::VertexMaxQM, MBMesquite::ElementAvgQM, MBMesquite::SizeMetric, MBMesquite::AddQualityMetric, and MBMesquite::EdgeLengthMetric.
Referenced by MBMesquite::AddQualityMetric::AddQualityMetric(), MBMesquite::StdDevTemplate::evaluate(), MBMesquite::LInfTemplate::evaluate(), MBMesquite::MaxTemplate::evaluate(), MBMesquite::PatchPowerMeanP::evaluate(), MBMesquite::VarianceTemplate::evaluate(), MBMesquite::LPtoPTemplate::evaluate(), MBMesquite::PMeanPTemplate::evaluate(), MBMesquite::StdDevTemplate::evaluate_with_gradient(), MBMesquite::PatchPowerMeanP::evaluate_with_gradient(), MBMesquite::VarianceTemplate::evaluate_with_gradient(), MBMesquite::LPtoPTemplate::evaluate_with_gradient(), MBMesquite::PMeanPTemplate::evaluate_with_gradient(), MBMesquite::PatchPowerMeanP::evaluate_with_Hessian(), MBMesquite::LPtoPTemplate::evaluate_with_Hessian(), MBMesquite::PMeanPTemplate::evaluate_with_Hessian(), MBMesquite::StdDevTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::VarianceTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::LPtoPTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::PMeanPTemplate::evaluate_with_Hessian_diagonal(), MBMesquite::AddQualityMetric::get_negate_flag(), MBMesquite::ElementAvgQM::get_negate_flag(), MBMesquite::ScalarAddQualityMetric::get_negate_flag(), MBMesquite::ElementMaxQM::get_negate_flag(), MBMesquite::VertexMaxQM::get_negate_flag(), MBMesquite::VertexPMeanP::get_negate_flag(), MBMesquite::ElementPMeanP::get_negate_flag(), MBMesquite::NumericalQM::get_negate_flag(), MBMesquite::ScalarMultiplyQualityMetric::get_negate_flag(), MBMesquite::PowerQualityMetric::get_negate_flag(), MBMesquite::MultiplyQualityMetric::get_negate_flag(), MBMesquite::CompareQM::get_negate_flag(), MBMesquite::MultiplyQualityMetric::MultiplyQualityMetric(), QualityMetricTester::test_domain_deviation_gradient(), QualityMetricTester::test_domain_deviation_quality(), QualityMetricTester::test_gradient_reflects_quality(), QualityMetricTester::test_ideal_element_positive_definite_Hessian(), QualityMetricTester::test_measures_quality(), QualityMetricTester::test_measures_transform(), QualityMetricTester::test_measures_vertex_quality(), and QualityMetricTester::test_vertex_gradient_reflects_quality().
void MBMesquite::QualityMetric::get_single_pass | ( | PatchData & | pd, |
std::vector< size_t > & | handles, | ||
bool | free_vertices_only, | ||
MsqError & | err | ||
) | [virtual] |
Get locations at which metric can be evaluated for use in BCD intialization and QualityAssessor.
For element-based, sample-based, and vertex-based metrics, this function is the same as get_evaluations. For edge-based metrics it returns only a subset of the results for get_evaluations such that each edge in the mesh is visited only once even though it would normally be visited twice when iterating over patches of the mesh. This assumes that no vertex occurs in more than one patch without its MSQ_PATCH_FIXED flag set. This assumption is true for both element-on-vertex and global patches.
pd | The patch |
handles | Output list of handles |
free_vertices_only | If true, only pass back evaluation points that depend on at least one free vertex. |
Reimplemented in MBMesquite::EdgeQM.
Definition at line 48 of file QualityMetric.cpp.
References get_evaluations().
Referenced by MBMesquite::PatchPowerMeanP::evaluate(), MBMesquite::VarianceTemplate::evaluate(), MBMesquite::LPtoPTemplate::evaluate(), MBMesquite::PMeanPTemplate::evaluate(), and MBMesquite::QualityAssessor::loop_over_mesh_internal().
{ get_evaluations( pd, handles, free_vertices_only, err ); }
void MBMesquite::QualityMetric::initialize_queue | ( | MeshDomainAssoc * | mesh_and_domain, |
const Settings * | settings, | ||
MsqError & | err | ||
) | [virtual] |
Called at start of instruction queue processing.
Do any preliminary global initialization, consistency checking, etc. Default implementation does nothing.
Reimplemented in MBMesquite::TMPQualityMetric.
Definition at line 46 of file QualityMetric.cpp.
Referenced by MBMesquite::ObjectiveFunctionTemplate::initialize_queue().
{}
void MBMesquite::QualityMetric::remove_fixed_diagonals | ( | EntityTopology | type, |
uint32_t | fixed_vertices, | ||
std::vector< Vector3D > & | gradients, | ||
std::vector< SymMatrix3D > & | hess_diagonal_blocks | ||
) | [static] |
Remove from vectors any gradient terms and hessian diagonal blcoks corresponding to a fixed vertex.
Remove terms from vector that correspond to fixed vertices.
type | Element type |
fixed_vertices | Bit flags, one per vertex, 1 if vertex is fixed. |
gradients | Array of gradients |
hess_diagonal_blocks | Array of diagonal blocks of Hessian matrix. |
Definition at line 329 of file QualityMetric.cpp.
References corners.
Referenced by MBMesquite::IdealWeightMeanRatio::evaluate_with_Hessian_diagonal(), and MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_Hessian_diagonal().
{ const unsigned num_vertex = TopologyInfo::corners( type ); unsigned r, w; for( r = 0; r < num_vertex && !( fixed & ( 1 << r ) ); ++r ) ; for( w = r++; r < num_vertex; ++r ) { if( !( fixed & ( 1 << r ) ) ) { grads[w] = grads[r]; diags[w] = diags[r]; ++w; } } grads.resize( w ); diags.resize( w ); }
void MBMesquite::QualityMetric::remove_fixed_gradients | ( | EntityTopology | type, |
uint32_t | fixed_vertices, | ||
std::vector< Vector3D > & | gradients | ||
) | [static] |
Remove from vector any gradient terms corresponding to a fixed vertex.
Remove terms from vector that correspond to fixed vertices.
type | Element type |
fixed_vertices | Bit flags, one per vertex, 1 if vertex is fixed. |
gradients | Array of gradients |
Definition at line 312 of file QualityMetric.cpp.
References corners.
Referenced by MBMesquite::IdealWeightMeanRatio::evaluate_with_gradient(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_gradient(), MBMesquite::IdealWeightMeanRatio::evaluate_with_Hessian(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_Hessian(), and QualityMetricTest::test_remove_fixed_gradients().
{ const unsigned num_vertex = TopologyInfo::corners( elem_type ); unsigned r, w; for( r = 0; r < num_vertex && !( fixed & ( 1 << r ) ); ++r ) ; for( w = r++; r < num_vertex; ++r ) { if( !( fixed & ( 1 << r ) ) ) { grads[w] = grads[r]; ++w; } } grads.resize( w ); }
void MBMesquite::QualityMetric::remove_fixed_hessians | ( | EntityTopology | type, |
uint32_t | fixed_vertices, | ||
std::vector< Matrix3D > & | hessians | ||
) | [static] |
Remove from vector any Hessian blocks corresponding to a fixed vertex.
Remove blocks from vector that correspond to fixed vertices.
type | Element type |
fixed_vertices | Bit flags, one per vertex, 1 if vertex is fixed. |
hessians | Array of Hessian blocks (upper trianguler, row-major) |
Definition at line 351 of file QualityMetric.cpp.
References corners.
Referenced by MBMesquite::IdealWeightMeanRatio::evaluate_with_Hessian(), MBMesquite::IdealWeightInverseMeanRatio::evaluate_with_Hessian(), and QualityMetricTest::test_remove_fixed_hessians().
{ const unsigned num_vertex = TopologyInfo::corners( elem_type ); unsigned r, c, i = 0, w = 0; for( r = 0; r < num_vertex; ++r ) { if( fixed & ( 1 << r ) ) { i += num_vertex - r; continue; } for( c = r; c < num_vertex; ++c ) { if( !( fixed & ( 1 << c ) ) ) { hessians[w] = hessians[i]; ++w; } ++i; } } hessians.resize( w ); }
static double MBMesquite::QualityMetric::vertex_barrier_function | ( | double | det, |
double | delta | ||
) | [inline, static] |
Escobar Barrier Function for Shape and Other Metrics.
Definition at line 227 of file QualityMetric.hpp.
double MBMesquite::QualityMetric::weighted_average_metrics | ( | const double | coef[], |
const double | metric_values[], | ||
const int & | num_values, | ||
MsqError & | err | ||
) |
takes an array of coefficients and an array of metrics (both of length num_value) and averages the contents using averaging method 'method'.
Definition at line 375 of file QualityMetric.cpp.
{ // MSQ_MAX needs to be made global? // double MSQ_MAX=1e10; double total_value = 0.0; int i = 0; // if no values, return zero if( num_values <= 0 ) { return 0.0; } for( i = 0; i < num_values; ++i ) { total_value += coef[i] * metric_values[i]; } total_value /= (double)num_values; return total_value; }
int MBMesquite::QualityMetric::feasible [private] |
Definition at line 331 of file QualityMetric.hpp.
double MBMesquite::QualityMetric::finiteDiffEps [private] |
Location for finite difference Hessian code to store this value so that it doesn't need to be recalculated if the gradient calculation is also finite difference
Definition at line 338 of file QualityMetric.hpp.
Referenced by evaluate_with_gradient(), and evaluate_with_Hessian().
bool MBMesquite::QualityMetric::haveFiniteDiffEps [private] |
True if finite difference Hessian code has calculated finiteDiffEps
Definition at line 336 of file QualityMetric.hpp.
Referenced by evaluate_with_gradient(), and evaluate_with_Hessian().
bool MBMesquite::QualityMetric::keepFiniteDiffEps [private] |
True if gradient finite difference calculation should set finiteDiffEps
Definition at line 334 of file QualityMetric.hpp.
Referenced by evaluate_with_gradient(), and evaluate_with_Hessian().
std::vector< Matrix3D > MBMesquite::QualityMetric::tmpHess [private] |
Definition at line 333 of file QualityMetric.hpp.
Referenced by evaluate_with_Hessian_diagonal().