51 #include "FeasibleNewton.hpp"
52 #include "MsqFreeVertexIndexIterator.hpp"
53 #include "MsqDebug.hpp"
55 using namespace Mesquite;
89 MSQ_DBGOUT(2) <<
"\no Performing Feasible Newton optimization.\n";
91 const double sigma = 1e-4;
92 const double beta0 = 0.25;
93 const double beta1 = 0.80;
94 const double tol1 = 1e-8;
95 const double tol2 = 1e-12;
96 const double epsilon = 1e-10;
97 double original_value, new_value;
101 msq_std::vector<Vector3D> grad_vect(nv), d_vect(nv);
118 MSQ_DBGOUT(3) <<
" o gradient norm: " <<
length(grad, nv) << msq_stdio::endl;
136 <<
")original coordinates:\n ";
141 while (ind1.next()) {
142 MSQ_DBGOUT(3) <<
"\t\t\t" << toto1[ind1.value()];
156 double alpha =
inner(grad, d, nv);
165 if (alpha > -epsilon) {
166 MSQ_PRINT(1)(
"Newton direction not guaranteed descent. Ensure preconditioner is positive definite.\n");
212 if ((fn_bool && (original_value - new_value >= -alpha*beta - epsilon)) ||
230 while (beta >= tol1) {
241 else if (original_value - new_value >= -alpha*beta - epsilon ) {
272 MSQ_PRINT(1)(
"Sufficient decrease not obtained in linesearch; switching to gradient.\n");
274 alpha =
inner(grad, grad, nv);
275 if (alpha < 1) alpha = 1;
276 for (i = 0; i < nv; ++
i) {
277 d[
i] = -grad[
i] / alpha;
279 alpha =
inner(grad, d, nv);
284 while (beta >= tol2) {
295 else if (original_value - new_value >= -alpha*beta - epsilon ) {
316 MSQ_PRINT(1)(
"Sufficient decrease not obtained with gradient; critical point likely found.\n");
327 MSQ_DBGOUT(3) <<
" o Free vertices new coordinates: \n";
332 MSQ_DBGOUT(3) <<
"\t\t\t" << toto1[ind.value()];
#define MSQ_DBG(flag)
Check if a debug flag is activated - evaluates to a bool.
bool compute_hessian(PatchData &patch, MsqHessian &hessian, Vector3D *const &grad, double &OF_val, MsqError &err)
Calls compute_analytical_hessian. Function returns 'false' if the patch is not within a required feas...
virtual void set_patch_type(PatchData::PatchType patch_type, MsqError &err, int param1=0, int param2=0)
Sets the Patch Type.
PatchDataVerticesMemento * coordsMem
Used to hold the error state and return it to the application.
void accumulate_patch(PatchData &pd, MsqError &err)
Common code for both inner and outer termination criteria during inner iteration. ...
virtual void optimize_vertex_positions(PatchData &pd, MsqError &err)
void recreate_vertices_memento(PatchDataVerticesMemento *memento, MsqError &err, bool include_higher_order=false)
reinstantiates a memento to holds the current state of the PatchData coordinates. ...
bool evaluate(PatchData &patch, double &fval, MsqError &err)
The TerminationCriterion class contains functionality to terminate the VertexMover's optimization...
void set_name(msq_std::string name)
provides a name to the QualityImprover (use it in constructor).
Vector3D is the object that effeciently stores information about about three-deminsional vectors...
checks the gradient of objective function against a double and stops when
virtual void terminate_mesh_iteration(PatchData &pd, MsqError &err)
bool terminate()
Check if termination criterion has been met.
void set_to_vertices_memento(PatchDataVerticesMemento *memento, MsqError &err)
Restore the PatchData coordinates to the state contained in the memento.
int num_free_vertices(MsqError &err) const
Returns the number of elements in the current patch who are free to move.
PatchDataVerticesMemento * create_vertices_memento(MsqError &err, bool include_higher_order=false)
Creates a memento that holds the current state of the PatchData coordinates.
double length(Vector3D *const v, int n)
#define MSQ_CHKERR(err)
Mesquite's Error Checking macro.
void cg_solver(Vector3D x[], Vector3D b[], MsqError &err)
bool compute_gradient(PatchData &patch, Vector3D *const &grad, double &OF_val, MsqError &err, size_t array_size=0)
Calls either compute_numerical_gradient or compute_analytical_gradient depending on the value of grad...
#define MSQ_SETERR(err)
Macro to set error - use err.clear() to clear.
ObjectiveFunction * objFunc
virtual void initialize(PatchData &pd, MsqError &err)
void add_criterion_type_with_double(TCType tc_type, double eps, MsqError &err)
Sets the criterion by specifing the TCType and the eps value.
size_t num_vertices() const
number of vertices in the patch.
void accumulate_inner(PatchData &pd, MsqError &err)
Accumulate data during inner iteration.
TerminationCriterion * get_inner_termination_criterion()
return the inner termination criterion pointer
void move_free_vertices_constrained(Vector3D dk[], size_t nb_vtx, double step_size, MsqError &err)
Moves free vertices and then snaps the free vertices to the domain.
const MsqVertex * get_vertex_array(MsqError &err) const
Returns a pointer to the start of the vertex array.
FeasibleNewton(ObjectiveFunction *of)
#define MSQ_PRINT(flag)
Check debug flag and print printf-style formatted output.
#define MSQ_FUNCTION_TIMER(NAME)
void reorder()
Reorders the mesh data.
void initialize(PatchData &pd, MsqError &err)
creates a sparse structure for a Hessian, based on the connectivity information contained in the Patc...
virtual void initialize_mesh_iteration(PatchData &pd, MsqError &err)
#define MSQ_DBGOUT(flag)
Check debug flag and return ostream associated with flag.
iterates over indexes of free vetices in a PatchData.
MsqVertex is the Mesquite object that stores information about the vertices in the mesh...
double inner(const Vector3D lhs[], const Vector3D rhs[], int n)
#define MSQ_ERRRTN(err)
If passed error is true, return from a void function.
Base class for concrete Objective Functions ObjectiveFunction contains a pointer to a QualityMetric...