MOAB: Mesh Oriented datABase
(version 5.4.1)
|
Base class for concrete Objective Functions ObjectiveFunction contains a pointer to a QualityMetric. If the ObjectiveFunction is associated with more than one QualityMetric (i.e., the Objective is a composite, and the composed ObjectiveFunctions are associated with different QualityMetrics), then the QualityMetric pointer is set to NULL.. More...
#include <ObjectiveFunction.hpp>
Public Types | |
enum | EvalType { CALCULATE, ACCUMULATE, SAVE, UPDATE, TEMPORARY } |
Public Member Functions | |
virtual | ~ObjectiveFunction () |
virtual void | initialize_queue (MeshDomainAssoc *mesh_and_domain, const Settings *settings, MsqError &err)=0 |
Called at start of instruction queue processing. | |
virtual bool | initialize_block_coordinate_descent (MeshDomainAssoc *mesh_and_domain, const Settings *settings, PatchSet *user_set, MsqError &err)=0 |
Initial accumulated value for block coordinate descent algorithms. | |
virtual bool | evaluate (EvalType type, PatchData &pd, double &value_out, bool free, MsqError &err)=0 |
Evaluate objective function for specified patch. | |
virtual 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 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 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 ObjectiveFunction * | clone () const =0 |
Create copy with same state. | |
virtual void | clear ()=0 |
virtual int | min_patch_layers () const =0 |
Protected Member Functions | |
double | get_eps (PatchData &pd, EvalType eval_type, double &local_val, int k, size_t vertex_index, MsqError &err) |
Returns eps used in the numerical gradient calculation. | |
Private Member Functions | |
bool | compute_subpatch_numerical_gradient (EvalType type, EvalType get_eps_eval_type, PatchData &pd, double &flocal, Vector3D &grad, MsqError &err) |
Compute numerical approx. of gradient for 1 vertex's coordinates. | |
bool | compute_patch_numerical_gradient (EvalType type, EvalType get_eps_eval_type, PatchData &pd, double &flocal, std::vector< Vector3D > &grad, MsqError &err) |
Base class for concrete Objective Functions ObjectiveFunction contains a pointer to a QualityMetric. If the ObjectiveFunction is associated with more than one QualityMetric (i.e., the Objective is a composite, and the composed ObjectiveFunctions are associated with different QualityMetrics), then the QualityMetric pointer is set to NULL..
Definition at line 76 of file ObjectiveFunction.hpp.
CALCULATE |
Do not modify or use any accumulated value in the calculation of the objective function value. Evaluate the objective function only over the passed patch. Used for NASH-type solutions. |
ACCUMULATE |
Incorporate the evaluation of the passed patch into the accumulated global objective function value. This EvalType is used by the default initialize implemenation, and need not be considered by ObjectiveFunction implementations that provide their own implementation of the initialize method. |
SAVE |
Save the evaluation over the passed patch such that it can later be removed from the accumulated global objective function value. Assume that the current accumulated value already includes the evaluation of this patch. Do *not* modify any accumulated, global value. |
UPDATE |
Assume that the patch passed to this call is a modified version of the patch previously passed with either the SAVE or UPDATE EvalType specified. Update the accumulated global value accordingly (remove previously saved local data from accumulated value, and incorporate evaluation over this patch into the accumulated value.) Save the necessary data such that the results of incorporaring the evaluation of this patch into the accumulated global value can be removed upon a subsequent call with this EvalType. |
TEMPORARY |
For this EvalType, the passed back results should be the same as for the UPDATE EvalType, but any accumulated data or data corresponding to the previous local values should *not* be modified. |
Definition at line 79 of file ObjectiveFunction.hpp.
{ /** Do not modify or use any accumulated value in the * calculation of the objective function value. * Evaluate the objective function only over the passed * patch. Used for NASH-type solutions. */ CALCULATE, /** Incorporate the evaluation of the passed patch into * the accumulated global objective function value. This * EvalType is used by the default initialize implemenation, * and need not be considered by ObjectiveFunction * implementations that provide their own implementation of * the initialize method. */ ACCUMULATE, /** Save the evaluation over the passed patch such that it * can later be removed from the accumulated global objective * function value. Assume that the current accumulated * value already includes the evaluation of this patch. Do * *not* modify any accumulated, global value. */ SAVE, /** Assume that the patch passed to this call is a modified * version of the patch previously passed with either the * SAVE or UPDATE EvalType specified. Update the accumulated * global value accordingly (remove previously saved local * data from accumulated value, and incorporate evaluation * over this patch into the accumulated value.) Save the * necessary data such that the results of incorporaring * the evaluation of this patch into the accumulated global * value can be removed upon a subsequent call with this * EvalType. */ UPDATE, /** For this EvalType, the passed back results should be the * same as for the UPDATE EvalType, but any accumulated data * or data corresponding to the previous local values should * *not* be modified. */ TEMPORARY };
MBMesquite::ObjectiveFunction::~ObjectiveFunction | ( | ) | [virtual] |
Definition at line 49 of file ObjectiveFunction.cpp.
{}
virtual void MBMesquite::ObjectiveFunction::clear | ( | ) | [pure virtual] |
Clear any values accumulated for BCD-related eval calls
Implemented in MBMesquite::PMeanPTemplate, MBMesquite::CompositeOFAdd, NumericalTestOF, MBMesquite::CompositeOFMultiply, MBMesquite::CompositeOFScalarMultiply, MBMesquite::CompositeOFScalarAdd, MBMesquite::VarianceTemplate, DummyOF, FauxObjectiveFunction, MBMesquite::LPtoPTemplate, MBMesquite::MaxTemplate, and MBMesquite::LInfTemplate.
Referenced by MBMesquite::CompositeOFScalarAdd::clear(), MBMesquite::CompositeOFScalarMultiply::clear(), MBMesquite::CompositeOFMultiply::clear(), MBMesquite::CompositeOFAdd::clear(), evaluate_with_gradient(), MBMesquite::ObjectiveFunctionTemplate::initialize_block_coordinate_descent(), ObjectiveFunctionTests::test_clone(), and ObjectiveFunctionTests::test_eval_type().
virtual ObjectiveFunction* MBMesquite::ObjectiveFunction::clone | ( | ) | const [pure 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.
Implemented in MBMesquite::PMeanPTemplate, MBMesquite::LPtoPTemplate, MBMesquite::CompositeOFAdd, MBMesquite::CompositeOFMultiply, NumericalTestOF, MBMesquite::CompositeOFScalarMultiply, MBMesquite::CompositeOFScalarAdd, MBMesquite::VarianceTemplate, MBMesquite::PatchPowerMeanP, DummyOF, FauxObjectiveFunction, MBMesquite::StdDevTemplate, MBMesquite::MaxTemplate, and MBMesquite::LInfTemplate.
Referenced by MBMesquite::CompositeOFScalarAdd::clone(), MBMesquite::CompositeOFScalarMultiply::clone(), MBMesquite::CompositeOFMultiply::clone(), MBMesquite::CompositeOFAdd::clone(), evaluate_with_gradient(), ObjectiveFunctionTests::test_clone(), and CompositeOFTest::test_composite_clone().
bool MBMesquite::ObjectiveFunction::compute_patch_numerical_gradient | ( | EvalType | type, |
EvalType | get_eps_eval_type, | ||
PatchData & | pd, | ||
double & | flocal, | ||
std::vector< Vector3D > & | grad, | ||
MsqError & | err | ||
) | [private] |
Definition at line 124 of file ObjectiveFunction.cpp.
References b, eps, evaluate(), get_eps(), MBMesquite::MsqError::INVALID_STATE, MSQ_CHKERR, MSQ_ERRZERO, MSQ_SETERR, MBMesquite::PatchData::num_free_vertices(), and MBMesquite::OF_FREE_EVALS_ONLY.
Referenced by evaluate_with_gradient().
{ double flocald = 0; double eps = 0; bool b = evaluate( type, pd, flocal, OF_FREE_EVALS_ONLY, err ); if( MSQ_CHKERR( err ) || !b ) { return false; } for( size_t i = 0; i < pd.num_free_vertices(); ++i ) { // loop over the three coords x,y,z for( int j = 0; j < 3; ++j ) { eps = get_eps( pd, subtype, flocald, j, i, err ); MSQ_ERRZERO( err ); if( eps == 0 ) { MSQ_SETERR( err ) ( "Dividing by zero in Objective Functions numerical grad", MsqError::INVALID_STATE ); return false; } grad[i][j] = ( flocald - flocal ) / eps; } } return true; }
bool MBMesquite::ObjectiveFunction::compute_subpatch_numerical_gradient | ( | EvalType | type, |
EvalType | get_eps_eval_type, | ||
PatchData & | pd, | ||
double & | flocal, | ||
Vector3D & | grad, | ||
MsqError & | err | ||
) | [private] |
Compute numerical approx. of gradient for 1 vertex's coordinates.
Compute the numerical approximation of the gradient of the objective function for the coordinates of a single vertex given a patch for which that vertex is the only free vertex.
type | Evaluation type. |
pd | A patch containing a single free vertex and the adjacent elements necessary for complete evaluation of the dependence of the objective fuction on the coordinates of the free vertex. |
flocal | The objective function value for the unmodified subpatch. (Output.) |
grad | The gradient of the OF with respect to the coordinates of the free veretx. (Output.) |
Definition at line 90 of file ObjectiveFunction.cpp.
References b, eps, evaluate(), get_eps(), MBMesquite::MsqError::INVALID_STATE, MSQ_CHKERR, MSQ_ERRZERO, MSQ_SETERR, MBMesquite::PatchData::num_free_vertices(), and MBMesquite::OF_FREE_EVALS_ONLY.
Referenced by evaluate_with_gradient().
{ assert( pd.num_free_vertices() == 1 ); double flocald = 0; double eps = 0; bool b = evaluate( type, pd, flocal, OF_FREE_EVALS_ONLY, err ); if( MSQ_CHKERR( err ) || !b ) { return false; } // loop over the three coords x,y,z for( int j = 0; j < 3; ++j ) { eps = get_eps( pd, subtype, flocald, j, 0, err ); MSQ_ERRZERO( err ); if( eps == 0 ) { MSQ_SETERR( err ) ( "Dividing by zero in Objective Functions numerical grad", MsqError::INVALID_STATE ); return false; } grad[j] = ( flocald - flocal ) / eps; } return true; }
virtual bool MBMesquite::ObjectiveFunction::evaluate | ( | EvalType | type, |
PatchData & | pd, | ||
double & | value_out, | ||
bool | free, | ||
MsqError & | err | ||
) | [pure 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.
type | Evaluation type. |
pd | The patch. |
value_out | The passed-back value of the objective fuction. |
free | If true, incorporate the quality metric values only for those metric evaluations that depend on at least one free vertex |
Implemented in MBMesquite::PMeanPTemplate, NumericalTestOF, MBMesquite::LPtoPTemplate, MBMesquite::CompositeOFAdd, DummyOF, MBMesquite::CompositeOFMultiply, MBMesquite::PatchPowerMeanP, MBMesquite::VarianceTemplate, FauxObjectiveFunction, MBMesquite::CompositeOFScalarMultiply, MBMesquite::CompositeOFScalarAdd, MBMesquite::MaxTemplate, MBMesquite::LInfTemplate, and MBMesquite::StdDevTemplate.
Referenced by ObjectiveFunctionTests::compare_numerical_hessian(), compute_patch_numerical_gradient(), compute_subpatch_numerical_gradient(), MBMesquite::CompositeOFScalarAdd::evaluate(), MBMesquite::CompositeOFScalarMultiply::evaluate(), MBMesquite::CompositeOFMultiply::evaluate(), MBMesquite::CompositeOFAdd::evaluate(), MBMesquite::OFEvaluator::evaluate(), ObjectiveFunctionTests::evaluate_internal(), evaluate_with_gradient(), get_eps(), MBMesquite::ObjectiveFunctionTemplate::initialize_block_coordinate_descent(), ObjectiveFunctionTests::test_clone(), CompositeOFTest::test_composite_clone(), CompositeOFTest::test_eval_fails(), CompositeOFTest::test_evaluate(), ObjectiveFunctionTests::test_handles_invalid_qm(), ObjectiveFunctionTests::test_handles_qm_error(), CompositeOFTest::test_invalid_eval(), ObjectiveFunctionTest::test_max_negate_flag(), ObjectiveFunctionTests::test_negate_flag(), and MBMesquite::OFEvaluator::update().
bool MBMesquite::ObjectiveFunction::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.
type | Evaluation type. |
pd | The patch. |
value_out | The passed-back value of the objective fuction. |
grad_out | The gradient of the OF wrt the coordinates of each *free* vertex in the patch. |
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.
pd | PatchData on which the gradient is taken. |
grad | Array of Vector3D of length the number of vertices used to store gradient. |
OF_val | will be set to the objective function value. |
Reimplemented in MBMesquite::PMeanPTemplate, MBMesquite::LPtoPTemplate, MBMesquite::CompositeOFAdd, DummyOF, MBMesquite::CompositeOFMultiply, MBMesquite::VarianceTemplate, MBMesquite::PatchPowerMeanP, MBMesquite::CompositeOFScalarMultiply, MBMesquite::CompositeOFScalarAdd, and MBMesquite::StdDevTemplate.
Definition at line 180 of file ObjectiveFunction.cpp.
References ACCUMULATE, b, CALCULATE, MBMesquite::MsqError::clear(), clear(), clone(), compute_patch_numerical_gradient(), compute_subpatch_numerical_gradient(), evaluate(), MBMesquite::PatchData::get_subpatch(), layers, min_patch_layers(), MSQ_CHKERR, MSQ_ERRZERO, MBMesquite::PatchData::num_free_vertices(), MBMesquite::OF_FREE_EVALS_ONLY, SAVE, and TEMPORARY.
Referenced by ObjectiveFunctionTests::compare_diagonal_gradient(), ObjectiveFunctionTests::compare_hessian_gradient(), ObjectiveFunctionTests::compare_numerical_gradient(), MBMesquite::OFEvaluator::evaluate(), ObjectiveFunctionTests::evaluate_internal(), MBMesquite::CompositeOFScalarAdd::evaluate_with_gradient(), MBMesquite::CompositeOFScalarMultiply::evaluate_with_gradient(), MBMesquite::CompositeOFMultiply::evaluate_with_gradient(), MBMesquite::CompositeOFAdd::evaluate_with_gradient(), NumericalOFTest::test_changed(), NumericalOFTest::test_gradient_values(), NumericalOFTest::test_handles_eval_failure(), NumericalOFTest::test_handles_eval_false(), ObjectiveFunctionTests::test_handles_invalid_qm(), ObjectiveFunctionTests::test_handles_qm_error(), ObjectiveFunctionTests::test_negate_flag(), NumericalOFTest::test_unchanged(), and MBMesquite::OFEvaluator::update().
{ bool b; grad.resize( pd.num_free_vertices() ); // Fast path for single-free-vertex patch if( pd.num_free_vertices() == 1 ) { const EvalType sub_type = ( eval_type == CALCULATE ) ? CALCULATE : TEMPORARY; b = compute_subpatch_numerical_gradient( eval_type, sub_type, pd, OF_val, grad[0], err ); return !MSQ_CHKERR( err ) && b; } ObjectiveFunction* of = this; std::unique_ptr< ObjectiveFunction > deleter; if( eval_type == CALCULATE ) { of->clear(); b = of->evaluate( ACCUMULATE, pd, OF_val, OF_FREE_EVALS_ONLY, err ); if( err ) { // OF doesn't support BCD type evals, try slow method err.clear(); of->clear(); b = compute_patch_numerical_gradient( CALCULATE, CALCULATE, pd, OF_val, grad, err ); return !MSQ_CHKERR( err ) && b; } else if( !b ) return b; } else { b = this->evaluate( eval_type, pd, OF_val, OF_FREE_EVALS_ONLY, err ); if( MSQ_CHKERR( err ) || !b ) return false; of = this->clone(); deleter = std::unique_ptr< ObjectiveFunction >( of ); } // Determine number of layers of adjacent elements based on metric type. unsigned layers = min_patch_layers(); // Create a subpatch for each free vertex and use it to evaluate the // gradient for that vertex. double flocal; PatchData subpatch; for( size_t i = 0; i < pd.num_free_vertices(); ++i ) { pd.get_subpatch( i, layers, subpatch, err ); MSQ_ERRZERO( err ); b = of->compute_subpatch_numerical_gradient( SAVE, TEMPORARY, subpatch, flocal, grad[i], err ); if( MSQ_CHKERR( err ) || !b ) { of->clear(); return false; } } of->clear(); return true; }
bool MBMesquite::ObjectiveFunction::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.
type | Evaluation type. |
pd | The patch. |
value_out | The passed-back value of the objective fuction. |
grad_out | The gradient of the OF wrt the coordinates of each *free* vertex in the patch. |
Hessian_out | The Hessian of the OF wrt the coordinates of each *free* vertex in the patch. |
Reimplemented in MBMesquite::PMeanPTemplate, MBMesquite::LPtoPTemplate, MBMesquite::CompositeOFAdd, MBMesquite::CompositeOFMultiply, MBMesquite::CompositeOFScalarMultiply, MBMesquite::CompositeOFScalarAdd, and MBMesquite::PatchPowerMeanP.
Definition at line 262 of file ObjectiveFunction.cpp.
References MBMesquite::MsqError::INVALID_STATE, and MSQ_SETERR.
Referenced by ObjectiveFunctionTests::compare_hessian_diagonal(), ObjectiveFunctionTests::compare_hessian_gradient(), ObjectiveFunctionTests::compare_numerical_hessian(), MBMesquite::OFEvaluator::evaluate(), ObjectiveFunctionTests::evaluate_internal(), MBMesquite::CompositeOFScalarAdd::evaluate_with_Hessian(), MBMesquite::CompositeOFScalarMultiply::evaluate_with_Hessian(), MBMesquite::CompositeOFAdd::evaluate_with_Hessian(), evaluate_with_Hessian_diagonal(), CompositeOFTest::get_hessians(), ObjectiveFunctionTests::test_handles_invalid_qm(), ObjectiveFunctionTests::test_handles_qm_error(), NumericalOFTest::test_Hessian_fails(), StdDevTemplateTest::test_hessian_fails(), StdDevTemplateTest::test_hessian_fails_sqr(), ObjectiveFunctionTests::test_negate_flag(), and MBMesquite::OFEvaluator::update().
{ MSQ_SETERR( err ) ( "No Hessian available for this objective function.\n" "Choose either a different objective function or a " "different solver.\n", MsqError::INVALID_STATE ); return false; }
bool MBMesquite::ObjectiveFunction::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.
type | Evaluation type. |
pd | The patch. |
value_out | The passed-back value of the objective fuction. |
grad_out | The gradient of the OF wrt the coordinates of each *free* vertex in the patch. |
hess_diag_out | The diagonal blocks of a Hessian. I.e. Decompose the Hessian into 3x3 submatrices and return only the submatrices (blocks) along the diagonal. |
Reimplemented in MBMesquite::PMeanPTemplate, MBMesquite::LPtoPTemplate, MBMesquite::CompositeOFAdd, MBMesquite::CompositeOFMultiply, MBMesquite::VarianceTemplate, MBMesquite::CompositeOFScalarMultiply, MBMesquite::CompositeOFScalarAdd, and MBMesquite::StdDevTemplate.
Definition at line 244 of file ObjectiveFunction.cpp.
References evaluate_with_Hessian(), MBMesquite::MsqHessian::get_block(), MBMesquite::hess(), MBMesquite::MsqHessian::initialize(), MSQ_ERRZERO, MBMesquite::MsqHessian::size(), and MBMesquite::Matrix3D::upper().
Referenced by ObjectiveFunctionTests::compare_diagonal_gradient(), ObjectiveFunctionTests::compare_hessian_diagonal(), ObjectiveFunctionTests::compare_numerical_hessian(), MBMesquite::OFEvaluator::evaluate(), ObjectiveFunctionTests::evaluate_internal(), MBMesquite::CompositeOFScalarAdd::evaluate_with_Hessian_diagonal(), MBMesquite::CompositeOFScalarMultiply::evaluate_with_Hessian_diagonal(), MBMesquite::CompositeOFMultiply::evaluate_with_Hessian_diagonal(), MBMesquite::CompositeOFAdd::evaluate_with_Hessian_diagonal(), ObjectiveFunctionTests::test_handles_invalid_qm(), ObjectiveFunctionTests::test_handles_qm_error(), ObjectiveFunctionTests::test_negate_flag(), and MBMesquite::OFEvaluator::update().
{ MsqHessian hess; hess.initialize( pd, err ); MSQ_ERRZERO( err ); bool val = evaluate_with_Hessian( type, pd, value_out, grad_out, hess, err ); MSQ_ERRZERO( err ); hess_diag_out.resize( hess.size() ); for( size_t i = 0; i < hess.size(); ++i ) hess_diag_out[i] = hess.get_block( i, i )->upper(); return val; }
double MBMesquite::ObjectiveFunction::get_eps | ( | PatchData & | pd, |
EvalType | type, | ||
double & | local_val, | ||
int | dim, | ||
size_t | vertex_index, | ||
MsqError & | err | ||
) | [protected] |
Returns eps used in the numerical gradient calculation.
Returns an appropiate value (eps) to use as a delta step for MsqVertex vertex in dimension k (i.e. k=0 -> x, k=1 -> y, k=2 -> z). The objective function value at the perturbed vertex position is given in local_val.
Definition at line 56 of file ObjectiveFunction.cpp.
References dim, eps, evaluate(), MBMesquite::PatchData::move_vertex(), MSQ_ERRZERO, MBMesquite::OF_FREE_EVALS_ONLY, MBMesquite::PatchData::set_vertex_coordinates(), and MBMesquite::PatchData::vertex_by_index().
Referenced by compute_patch_numerical_gradient(), and compute_subpatch_numerical_gradient().
{ double eps = 1.e-07; const double rho = 0.5; const int imax = 20; bool feasible = false; double tmp_var = 0.0; Vector3D delta( 0, 0, 0 ); for( int i = 0; i < imax; ++i ) { // perturb kth coord val and check feas if needed tmp_var = pd.vertex_by_index( vertex_index )[dim]; delta[dim] = eps; pd.move_vertex( delta, vertex_index, err ); feasible = evaluate( type, pd, local_val, OF_FREE_EVALS_ONLY, err ); MSQ_ERRZERO( err ); // revert kth coord val delta = pd.vertex_by_index( vertex_index ); delta[dim] = tmp_var; pd.set_vertex_coordinates( delta, vertex_index, err ); // if step was too big, shorten it and go again if( feasible ) return eps; else eps *= rho; } // end while looking for feasible eps return 0.0; } // end function get_eps
virtual bool MBMesquite::ObjectiveFunction::initialize_block_coordinate_descent | ( | MeshDomainAssoc * | mesh_and_domain, |
const Settings * | settings, | ||
PatchSet * | user_set, | ||
MsqError & | err | ||
) | [pure virtual] |
Initial accumulated value for block coordinate descent algorithms.
Set accumulated value of objective function to the value for the entire, unmodified mesh. This is the initial state for a block coordinate descent algorithm. The ObjectiveFunction will asked to add or remove values for a specific patch of the mesh during the optimization.
mesh | The Mesh |
domain | The MeshDomain |
user_set | User-defined patch set - not relevant for most OF templates. |
Implemented in MBMesquite::CompositeOFAdd, DummyOF, MBMesquite::CompositeOFMultiply, MBMesquite::CompositeOFScalarMultiply, MBMesquite::CompositeOFScalarAdd, FauxObjectiveFunction, and MBMesquite::ObjectiveFunctionTemplate.
Referenced by MBMesquite::OFEvaluator::initialize(), MBMesquite::CompositeOFScalarAdd::initialize_block_coordinate_descent(), MBMesquite::CompositeOFScalarMultiply::initialize_block_coordinate_descent(), MBMesquite::CompositeOFMultiply::initialize_block_coordinate_descent(), and MBMesquite::CompositeOFAdd::initialize_block_coordinate_descent().
virtual void MBMesquite::ObjectiveFunction::initialize_queue | ( | MeshDomainAssoc * | mesh_and_domain, |
const Settings * | settings, | ||
MsqError & | err | ||
) | [pure virtual] |
Called at start of instruction queue processing.
Do any preliminary global initialization, consistency checking, etc. This function is pure-virtual (abstract) in this class because every practical OF implementation at this time should have an implementation that at least recursively calls the same function on the underlying QualityMetric or ObjectiveFunction(s).
Implemented in FauxObjectiveFunction, DummyOF, MBMesquite::ObjectiveFunctionTemplate, MBMesquite::CompositeOFAdd, MBMesquite::CompositeOFMultiply, MBMesquite::CompositeOFScalarMultiply, and MBMesquite::CompositeOFScalarAdd.
Referenced by MBMesquite::CompositeOFScalarAdd::initialize_queue(), MBMesquite::CompositeOFScalarMultiply::initialize_queue(), MBMesquite::CompositeOFMultiply::initialize_queue(), MBMesquite::OFEvaluator::initialize_queue(), and MBMesquite::CompositeOFAdd::initialize_queue().
virtual int MBMesquite::ObjectiveFunction::min_patch_layers | ( | ) | const [pure virtual] |
Get the minimum number of layers of adjacent elements required in a patch to evaluate the objective function for a single free vertex.
Implemented in MBMesquite::CompositeOFAdd, MBMesquite::CompositeOFMultiply, MBMesquite::CompositeOFScalarMultiply, MBMesquite::CompositeOFScalarAdd, DummyOF, FauxObjectiveFunction, and MBMesquite::ObjectiveFunctionTemplate.
Referenced by evaluate_with_gradient(), MBMesquite::CompositeOFScalarAdd::min_patch_layers(), MBMesquite::CompositeOFScalarMultiply::min_patch_layers(), MBMesquite::CompositeOFMultiply::min_patch_layers(), and MBMesquite::CompositeOFAdd::min_patch_layers().