Sketcher: New Features: SparseQR decomposition and Solver advanced control TaskBox

==================================================================================

The solver has been adapted to use Eigen's SparseQR QR decomposition algorithm. The original
Dense QR implementation is maintained and can be selected using the Advanced Control TaskBox (see below).

The use of SparseQR provides over an order of magnitude improvement in solving time in complex sketches due to
the Sparse nature of the Jacobian matrix of the system of equations.

The solver advanced control is a new TaskBox in the Sketcher that allows to select which algorithms are to be used for
the different solving operations and tweak its parameters. It is not intended to be a user control, but means to debug
solving problems and improve the algorithms and their configuration.

This commit also introduces multithread support for Eigen. Currently it is only limited to products and does not provide
a substantial speed improvement. It is expected to have more multithreaded operations in Eigen in the future.

As a bonus, the TaskBoxes in the Taskbar of the Sketcher remember the last state (collapsed or deployed).
This commit is contained in:
Abdullah Tahiri 2015-06-14 19:30:29 +02:00 committed by wmayer
parent 02df1acb5f
commit 796c9d79d4
11 changed files with 1555 additions and 139 deletions

View File

@ -4,6 +4,14 @@ else(MSVC)
add_definitions(-DHAVE_LIMITS_H -DHAVE_CONFIG_H)
endif(MSVC)
find_package(OpenMP)
if (OPENMP_FOUND)
set (CMAKE_C_FLAGS "${CMAKE_C_FLAGS} ${OpenMP_C_FLAGS}")
set (CMAKE_CXX_FLAGS "${CMAKE_CXX_FLAGS} ${OpenMP_CXX_FLAGS}")
endif()
include_directories(
${CMAKE_BINARY_DIR}
${CMAKE_SOURCE_DIR}/src

View File

@ -55,7 +55,6 @@
#include <iostream>
using namespace Sketcher;
using namespace Base;
using namespace Part;
@ -63,7 +62,8 @@ using namespace Part;
TYPESYSTEM_SOURCE(Sketcher::Sketch, Base::Persistence)
Sketch::Sketch()
: GCSsys(), ConstraintsCounter(0), isInitMove(false)
: GCSsys(), ConstraintsCounter(0), isInitMove(false),
defaultSolver(GCS::DogLeg),defaultSolverRedundant(GCS::DogLeg),debugMode(GCS::Minimal)
{
}
@ -110,6 +110,9 @@ int Sketch::setUpSketch(const std::vector<Part::Geometry *> &GeoList,
const std::vector<Constraint *> &ConstraintList,
int extGeoCount)
{
Base::TimeInfo start_time;
clear();
std::vector<Part::Geometry *> intGeoList, extGeoList;
@ -131,9 +134,16 @@ int Sketch::setUpSketch(const std::vector<Part::Geometry *> &GeoList,
}
GCSsys.clearByTag(-1);
GCSsys.declareUnknowns(Parameters);
GCSsys.initSolution();
GCSsys.initSolution(defaultSolverRedundant);
GCSsys.getConflicting(Conflicting);
GCSsys.getRedundant(Redundant);
if(debugMode==GCS::Minimal || debugMode==GCS::IterationLevel) {
Base::TimeInfo end_time;
Base::Console().Log("Sketcher::setUpSketch()-T:%s\n",Base::TimeInfo::diffTime(start_time,end_time).c_str());
}
return GCSsys.dofsNumber();
}
@ -2086,77 +2096,137 @@ int Sketch::solve(void)
GCSsys.clearByTag(-1);
isFine = true;
}
int ret;
bool valid_solution;
for (int soltype=0; soltype < (isInitMove ? 1 : 4); soltype++) {
std::string solvername;
switch (soltype) {
case 0: // solving with the default DogLeg solver
// (or with SQP if we are in moving mode)
solvername = isInitMove ? "SQP" : "DogLeg";
ret = GCSsys.solve(isFine, GCS::DogLeg);
break;
case 1: // solving with the LevenbergMarquardt solver
solvername = "LevenbergMarquardt";
ret = GCSsys.solve(isFine, GCS::LevenbergMarquardt);
break;
case 2: // solving with the BFGS solver
solvername = "BFGS";
ret = GCSsys.solve(isFine, GCS::BFGS);
break;
case 3: // last resort: augment the system with a second subsystem and use the SQP solver
solvername = "SQP(augmented system)";
InitParameters.resize(Parameters.size());
int i=0;
for (std::vector<double*>::iterator it = Parameters.begin(); it != Parameters.end(); ++it, i++) {
InitParameters[i] = **it;
GCSsys.addConstraintEqual(*it, &InitParameters[i], -1);
}
GCSsys.initSolution();
ret = GCSsys.solve(isFine);
break;
std::string solvername;
int defaultsoltype;
if(isInitMove){
solvername = "DogLeg"; // DogLeg is used for dragging (same as before)
ret = GCSsys.solve(isFine, GCS::DogLeg);
}
else{
switch (defaultSolver) {
case 0:
solvername = "BFGS";
ret = GCSsys.solve(isFine, GCS::BFGS);
defaultsoltype=2;
break;
case 1: // solving with the LevenbergMarquardt solver
solvername = "LevenbergMarquardt";
ret = GCSsys.solve(isFine, GCS::LevenbergMarquardt);
defaultsoltype=1;
break;
case 2: // solving with the BFGS solver
solvername = "DogLeg";
ret = GCSsys.solve(isFine, GCS::DogLeg);
defaultsoltype=0;
break;
}
}
// if successfully solved try to write the parameters back
if (ret == GCS::Success) {
GCSsys.applySolution();
valid_solution = updateGeometry();
if (!valid_solution) {
GCSsys.undoSolution();
updateGeometry();
Base::Console().Warning("Invalid solution from %s solver.\n", solvername.c_str());
}else
{
updateNonDrivingConstraints();
}
} else {
valid_solution = false;
if(debugMode==GCS::Minimal || debugMode==GCS::IterationLevel){
Base::Console().Log("Sketcher::Solve()-%s- Failed!! Falling back...\n",solvername.c_str());
}
}
// if successfully solved try to write the parameters back
if (ret == GCS::Success) {
GCSsys.applySolution();
valid_solution = updateGeometry();
if (!valid_solution) {
GCSsys.undoSolution();
updateGeometry();
Base::Console().Warning("Invalid solution from %s solver.\n", solvername.c_str());
}else
{
updateNonDrivingConstraints();
if(!valid_solution && !isInitMove) { // Fall back to other solvers
for (int soltype=0; soltype < 4; soltype++) {
if(soltype==defaultsoltype){
continue; // skip default solver
}
} else {
valid_solution = false;
//Base::Console().Log("NotSolved ");
}
if (soltype == 3) // cleanup temporary constraints of the augmented system
GCSsys.clearByTag(-1);
if (valid_solution) {
if (soltype == 1)
Base::Console().Log("Important: the LevenbergMarquardt solver succeeded where the DogLeg solver had failed.\n");
else if (soltype == 2)
Base::Console().Log("Important: the BFGS solver succeeded where the DogLeg and LevenbergMarquardt solvers have failed.\n");
else if (soltype == 3)
Base::Console().Log("Important: the SQP solver succeeded where all single subsystem solvers have failed.\n");
if (soltype > 0) {
Base::Console().Log("If you see this message please report a way of reproducing this result at\n");
Base::Console().Log("http://www.freecadweb.org/tracker/main_page.php\n");
switch (soltype) {
case 0:
solvername = "DogLeg";
ret = GCSsys.solve(isFine, GCS::DogLeg);
break;
case 1: // solving with the LevenbergMarquardt solver
solvername = "LevenbergMarquardt";
ret = GCSsys.solve(isFine, GCS::LevenbergMarquardt);
break;
case 2: // solving with the BFGS solver
solvername = "BFGS";
ret = GCSsys.solve(isFine, GCS::BFGS);
break;
case 3: // last resort: augment the system with a second subsystem and use the SQP solver
solvername = "SQP(augmented system)";
InitParameters.resize(Parameters.size());
int i=0;
for (std::vector<double*>::iterator it = Parameters.begin(); it != Parameters.end(); ++it, i++) {
InitParameters[i] = **it;
GCSsys.addConstraintEqual(*it, &InitParameters[i], -1);
}
GCSsys.initSolution();
ret = GCSsys.solve(isFine);
break;
}
break;
}
} // soltype
// if successfully solved try to write the parameters back
if (ret == GCS::Success) {
GCSsys.applySolution();
valid_solution = updateGeometry();
if (!valid_solution) {
GCSsys.undoSolution();
updateGeometry();
Base::Console().Warning("Invalid solution from %s solver.\n", solvername.c_str());
}else
{
updateNonDrivingConstraints();
}
} else {
valid_solution = false;
if(debugMode==GCS::Minimal || debugMode==GCS::IterationLevel){
Base::Console().Log("Sketcher::Solve()-%s- Failed!! Falling back...\n",solvername.c_str());
}
}
if (soltype == 3) // cleanup temporary constraints of the augmented system
GCSsys.clearByTag(-1);
if (valid_solution) {
if (soltype == 1)
Base::Console().Log("Important: the LevenbergMarquardt solver succeeded where the DogLeg solver had failed.\n");
else if (soltype == 2)
Base::Console().Log("Important: the BFGS solver succeeded where the DogLeg and LevenbergMarquardt solvers have failed.\n");
else if (soltype == 3)
Base::Console().Log("Important: the SQP solver succeeded where all single subsystem solvers have failed.\n");
if (soltype > 0) {
Base::Console().Log("If you see this message please report a way of reproducing this result at\n");
Base::Console().Log("http://www.freecadweb.org/tracker/main_page.php\n");
}
break;
}
} // soltype
}
Base::TimeInfo end_time;
//Base::Console().Log("T:%s\n",Base::TimeInfo::diffTime(start_time,end_time).c_str());
if(debugMode==GCS::Minimal || debugMode==GCS::IterationLevel){
Base::Console().Log("Sketcher::Solve()-%s-T:%s\n",solvername.c_str(),Base::TimeInfo::diffTime(start_time,end_time).c_str());
}
SolveTime = Base::TimeInfo::diffTimeF(start_time,end_time);
return ret;
}

View File

@ -374,6 +374,33 @@ protected:
bool isInitMove;
bool isFine;
public:
GCS::Algorithm defaultSolver;
GCS::Algorithm defaultSolverRedundant;
inline void setDebugMode(GCS::DebugMode mode) {debugMode=mode;GCSsys.debugMode=mode;}
inline GCS::DebugMode getDebugMode(void) {return debugMode;}
inline void setMaxIter(int maxiter){GCSsys.maxIter=maxiter;}
inline void setMaxIterRedundant(int maxiter){GCSsys.maxIterRedundant=maxiter;}
inline void setSketchSizeMultiplier(bool mult){GCSsys.sketchSizeMultiplier=mult;}
inline void setSketchSizeMultiplierRedundant(bool mult){GCSsys.sketchSizeMultiplierRedundant=mult;}
inline void setConvergence(double conv){GCSsys.convergence=conv;}
inline void setConvergenceRedundant(double conv){GCSsys.convergenceRedundant=conv;}
inline void setQRAlgorithm(GCS::QRAlgorithm alg){GCSsys.qrAlgorithm=alg;}
inline void setLM_eps(double val){GCSsys.LM_eps=val;}
inline void setLM_eps1(double val){GCSsys.LM_eps1=val;}
inline void setLM_tau(double val){GCSsys.LM_tau=val;}
inline void setDL_tolg(double val){GCSsys.DL_tolg=val;}
inline void setDL_tolx(double val){GCSsys.DL_tolx=val;}
inline void setDL_tolf(double val){GCSsys.DL_tolf=val;}
inline void setLM_epsRedundant(double val){GCSsys.LM_epsRedundant=val;}
inline void setLM_eps1Redundant(double val){GCSsys.LM_eps1Redundant=val;}
inline void setLM_tauRedundant(double val){GCSsys.LM_tauRedundant=val;}
inline void setDL_tolgRedundant(double val){GCSsys.DL_tolgRedundant=val;}
inline void setDL_tolxRedundant(double val){GCSsys.DL_tolxRedundant=val;}
inline void setDL_tolfRedundant(double val){GCSsys.DL_tolfRedundant=val;}
protected:
GCS::DebugMode debugMode;
private:

View File

@ -26,15 +26,16 @@
#include "GCS.h"
#include "qp_eq.h"
#include <Eigen/QR>
#include <Eigen/Sparse>
#include <Eigen/OrderingMethods>
#undef _GCS_DEBUG
#undef _GCS_DEBUG_SOLVER_JACOBIAN_QR_DECOMPOSITION_TRIANGULAR_MATRIX
#if defined(_GCS_DEBUG) || defined(_GCS_DEBUG_SOLVER_JACOBIAN_QR_DECOMPOSITION_TRIANGULAR_MATRIX)
#include <FCConfig.h>
#include <Base/Console.h>
#endif // _GCS_DEBUG
#include <boost/graph/adjacency_list.hpp>
#include <boost/graph/connected_components.hpp>
@ -153,8 +154,19 @@ System::System()
c2p(), p2c(),
subSystems(0), subSystemsAux(0),
reference(0),
hasUnknowns(false), hasDiagnosis(false), isInit(false)
hasUnknowns(false), hasDiagnosis(false), isInit(false),
maxIter(100), maxIterRedundant(100),
sketchSizeMultiplier(true), sketchSizeMultiplierRedundant(true),
convergence(1e-10), convergenceRedundant(1e-10),
qrAlgorithm(EigenSparseQR), debugMode(Minimal),
LM_eps(1E-10), LM_eps1(1E-80), LM_tau(1E-3),
DL_tolg(1E-80), DL_tolx(1E-80), DL_tolf(1E-10),
LM_epsRedundant(1E-10), LM_eps1Redundant(1E-80), LM_tauRedundant(1E-3),
DL_tolgRedundant(1E-80), DL_tolxRedundant(1E-80), DL_tolfRedundant(1E-10)
{
// currently Eigen only supports multithreading for multiplications
// There is no appreciable gain from using more threads
Eigen::setNbThreads(1);
}
/*DeepSOIC: seriously outdated, needs redesign
@ -860,7 +872,7 @@ void System::declareUnknowns(VEC_pD &params)
hasUnknowns = true;
}
void System::initSolution()
void System::initSolution(Algorithm alg)
{
// - Stores the current parameters values in the vector "reference"
// - identifies any decoupled subsystems and partitions the original
@ -881,7 +893,7 @@ void System::initSolution()
// diagnose conflicting or redundant constraints
if (!hasDiagnosis) {
diagnose();
diagnose(alg);
if (!hasDiagnosis)
return;
}
@ -1007,14 +1019,14 @@ void System::resetToReference()
}
}
int System::solve(VEC_pD &params, bool isFine, Algorithm alg)
int System::solve(VEC_pD &params, bool isFine, Algorithm alg, bool isRedundantsolving)
{
declareUnknowns(params);
initSolution();
return solve(isFine, alg);
return solve(isFine, alg, isRedundantsolving);
}
int System::solve(bool isFine, Algorithm alg)
int System::solve(bool isFine, Algorithm alg, bool isRedundantsolving)
{
if (!isInit)
return Failed;
@ -1029,11 +1041,11 @@ int System::solve(bool isFine, Algorithm alg)
isReset = true;
}
if (subSystems[cid] && subSystemsAux[cid])
res = std::max(res, solve(subSystems[cid], subSystemsAux[cid], isFine));
res = std::max(res, solve(subSystems[cid], subSystemsAux[cid], isFine, isRedundantsolving));
else if (subSystems[cid])
res = std::max(res, solve(subSystems[cid], isFine, alg));
res = std::max(res, solve(subSystems[cid], isFine, alg, isRedundantsolving));
else if (subSystemsAux[cid])
res = std::max(res, solve(subSystemsAux[cid], isFine, alg));
res = std::max(res, solve(subSystemsAux[cid], isFine, alg, isRedundantsolving));
}
if (res == Success) {
for (std::set<Constraint *>::const_iterator constr=redundant.begin();
@ -1042,7 +1054,7 @@ int System::solve(bool isFine, Algorithm alg)
//convergence, which makes no sense. Potentially I fixed bug, and
//chances are low I've broken anything.
double err = (*constr)->error();
if (err*err > XconvergenceFine) {
if (err*err > (isRedundantsolving?convergenceRedundant:convergence)) {
res = Converged;
return res;
}
@ -1051,19 +1063,19 @@ int System::solve(bool isFine, Algorithm alg)
return res;
}
int System::solve(SubSystem *subsys, bool isFine, Algorithm alg)
int System::solve(SubSystem *subsys, bool isFine, Algorithm alg, bool isRedundantsolving)
{
if (alg == BFGS)
return solve_BFGS(subsys, isFine);
return solve_BFGS(subsys, isFine, isRedundantsolving);
else if (alg == LevenbergMarquardt)
return solve_LM(subsys);
return solve_LM(subsys, isRedundantsolving);
else if (alg == DogLeg)
return solve_DL(subsys);
return solve_DL(subsys, isRedundantsolving);
else
return Failed;
}
int System::solve_BFGS(SubSystem *subsys, bool isFine)
int System::solve_BFGS(SubSystem *subsys, bool isFine, bool isRedundantsolving)
{
int xsize = subsys->pSize();
if (xsize == 0)
@ -1092,13 +1104,17 @@ int System::solve_BFGS(SubSystem *subsys, bool isFine)
subsys->getParams(x);
h = x - h; // = x - xold
double convergence = isFine ? XconvergenceFine : XconvergenceRough;
int maxIterNumber = MaxIterations * xsize;
//double convergence = isFine ? convergence : XconvergenceRough;
int maxIterNumber = isRedundantsolving?
(sketchSizeMultiplierRedundant?maxIterRedundant * xsize:maxIterRedundant):
(sketchSizeMultiplier?maxIter * xsize:maxIter);
double divergingLim = 1e6*err + 1e12;
double h_norm;
for (int iter=1; iter < maxIterNumber; iter++) {
if (h.norm() <= convergence || err <= smallF)
h_norm = h.norm();
if (h_norm <= isRedundantsolving?convergenceRedundant:convergence || err <= smallF)
break;
if (err > divergingLim || err != err) // check for diverging and NaN
break;
@ -1127,18 +1143,30 @@ int System::solve_BFGS(SubSystem *subsys, bool isFine)
h = x;
subsys->getParams(x);
h = x - h; // = x - xold
if(debugMode==IterationLevel) {
std::stringstream stream;
// Iteration: 1, residual: 1e-3, tolg: 1e-5, tolx: 1e-3
stream << "BFGS, Iteration: " << iter
<< ", err(eps): " << err
<< ", h_norm: " << h_norm << "\n";
const std::string tmp = stream.str();
Base::Console().Log(tmp.c_str());
}
}
subsys->revertParams();
if (err <= smallF)
return Success;
if (h.norm() <= convergence)
if (h.norm() <= isRedundantsolving?convergenceRedundant:convergence)
return Converged;
return Failed;
}
int System::solve_LM(SubSystem* subsys)
int System::solve_LM(SubSystem* subsys, bool isRedundantsolving)
{
int xsize = subsys->pSize();
int csize = subsys->cSize();
@ -1157,11 +1185,15 @@ int System::solve_LM(SubSystem* subsys)
subsys->calcResidual(e);
e*=-1;
int maxIterNumber = MaxIterations * xsize;
int maxIterNumber = isRedundantsolving?
(sketchSizeMultiplierRedundant?maxIterRedundant * xsize:maxIterRedundant):
(sketchSizeMultiplier?maxIter * xsize:maxIter);
double divergingLim = 1e6*e.squaredNorm() + 1e12;
double eps=1e-10, eps1=1e-80;
double tau=1e-3;
double eps=isRedundantsolving?LM_epsRedundant:LM_eps;
double eps1=isRedundantsolving?LM_eps1Redundant:LM_eps1;
double tau=isRedundantsolving?LM_tauRedundant:LM_tau;
double nu=2, mu=0;
int iter=0, stop=0;
for (iter=0; iter < maxIterNumber && !stop; ++iter) {
@ -1197,6 +1229,7 @@ int System::solve_LM(SubSystem* subsys)
if (iter == 0)
mu = tau * diag_A.lpNorm<Eigen::Infinity>();
double h_norm;
// determine increment using adaptive damping
int k=0;
while (k < 50) {
@ -1218,7 +1251,7 @@ int System::solve_LM(SubSystem* subsys)
// compute par's new estimate and ||d_par||^2
x_new = x + h;
double h_norm = h.squaredNorm();
h_norm = h.squaredNorm();
if (h_norm <= eps1*eps1*x.norm()) { // relative change in p is small, stop
stop = 3;
@ -1262,6 +1295,19 @@ int System::solve_LM(SubSystem* subsys)
stop = 7;
break;
}
if(debugMode==IterationLevel) {
std::stringstream stream;
// Iteration: 1, residual: 1e-3, tolg: 1e-5, tolx: 1e-3
stream << "LM, Iteration: " << iter
<< ", err(eps): " << err
<< ", g_inf(eps1): " << g_inf
<< ", h_norm: " << h_norm << "\n";
const std::string tmp = stream.str();
Base::Console().Log(tmp.c_str());
}
}
if (iter >= maxIterNumber)
@ -1273,9 +1319,11 @@ int System::solve_LM(SubSystem* subsys)
}
int System::solve_DL(SubSystem* subsys)
int System::solve_DL(SubSystem* subsys, bool isRedundantsolving)
{
double tolg=1e-80, tolx=1e-80, tolf=1e-10;
double tolg=isRedundantsolving?DL_tolgRedundant:DL_tolg;
double tolx=isRedundantsolving?DL_tolxRedundant:DL_tolx;
double tolf=isRedundantsolving?DL_tolfRedundant:DL_tolf;
int xsize = subsys->pSize();
int csize = subsys->cSize();
@ -1301,7 +1349,10 @@ int System::solve_DL(SubSystem* subsys)
double g_inf = g.lpNorm<Eigen::Infinity>();
double fx_inf = fx.lpNorm<Eigen::Infinity>();
int maxIterNumber = MaxIterations * xsize;
int maxIterNumber = isRedundantsolving?
(sketchSizeMultiplierRedundant?maxIterRedundant * xsize:maxIterRedundant):
(sketchSizeMultiplier?maxIter * xsize:maxIter);
double divergingLim = 1e6*err + 1e12;
double delta=0.1;
@ -1412,7 +1463,21 @@ int System::solve_DL(SubSystem* subsys)
}
else
reduce--;
if(debugMode==IterationLevel) {
std::stringstream stream;
// Iteration: 1, residual: 1e-3, tolg: 1e-5, tolx: 1e-3
stream << "DL, Iteration: " << iter
<< ", fx_inf(tolf): " << fx_inf
<< ", g_inf(tolg): " << g_inf
<< ", delta(f(tolx)): " << delta
<< ", err(divergingLim): " << err << "\n";
const std::string tmp = stream.str();
Base::Console().Log(tmp.c_str());
}
// count this iteration and start again
iter++;
}
@ -1424,7 +1489,7 @@ int System::solve_DL(SubSystem* subsys)
// The following solver variant solves a system compound of two subsystems
// treating the first of them as of higher priority than the second
int System::solve(SubSystem *subsysA, SubSystem *subsysB, bool isFine)
int System::solve(SubSystem *subsysA, SubSystem *subsysB, bool isFine, bool isRedundantsolving)
{
int xsizeA = subsysA->pSize();
int xsizeB = subsysB->pSize();
@ -1470,8 +1535,11 @@ int System::solve(SubSystem *subsysA, SubSystem *subsysB, bool isFine)
subsysA->calcJacobi(plistAB,JA);
subsysA->calcResidual(resA);
double convergence = isFine ? XconvergenceFine : XconvergenceRough;
int maxIterNumber = MaxIterations;
//double convergence = isFine ? XconvergenceFine : XconvergenceRough;
int maxIterNumber = isRedundantsolving?
(sketchSizeMultiplierRedundant?maxIterRedundant * xsize:maxIterRedundant):
(sketchSizeMultiplier?maxIter * xsize:maxIter);
double divergingLim = 1e6*subsysA->error() + 1e12;
double mu = 0;
@ -1565,7 +1633,7 @@ int System::solve(SubSystem *subsysA, SubSystem *subsysB, bool isFine)
}
double err = subsysA->error();
if (h.norm() <= convergence && err <= smallF)
if (h.norm() <= isRedundantsolving?convergenceRedundant:convergence && err <= smallF)
break;
if (err > divergingLim || err != err) // check for diverging and NaN
break;
@ -1574,7 +1642,7 @@ int System::solve(SubSystem *subsysA, SubSystem *subsysB, bool isFine)
int ret;
if (subsysA->error() <= smallF)
ret = Success;
else if (h.norm() <= convergence)
else if (h.norm() <= isRedundantsolving?convergenceRedundant:convergence)
ret = Converged;
else
ret = Failed;
@ -1603,7 +1671,7 @@ void System::undoSolution()
resetToReference();
}
int System::diagnose()
int System::diagnose(Algorithm alg)
{
// Analyses the constrainess grad of the system and provides feedback
// The vector "conflictingTags" will hold a group of conflicting constraints
@ -1636,6 +1704,16 @@ int System::diagnose()
}
}
Eigen::SparseMatrix<double> SJ;
if(qrAlgorithm==EigenSparseQR){
// this creation is not optimized (done using triplets)
// however the time this takes is negligible compared to the
// time the QR decomposition itself takes
SJ = J.sparseView();
SJ.makeCompressed();
}
#ifdef _GCS_DEBUG
// Debug code starts
std::stringstream stream;
@ -1646,26 +1724,74 @@ int System::diagnose()
const std::string tmp = stream.str();
Base::Console().Warning(tmp.c_str());
Base::Console().Log(tmp.c_str());
// Debug code ends
#endif
if (J.rows() > 0) {
Eigen::FullPivHouseholderQR<Eigen::MatrixXd> qrJT(J.topRows(count).transpose());
Eigen::MatrixXd Q = qrJT.matrixQ ();
int paramsNum = qrJT.rows();
int constrNum = qrJT.cols();
qrJT.setThreshold(1e-13);
int rank = qrJT.rank();
Eigen::MatrixXd R;
int paramsNum;
int constrNum;
int rank;
Eigen::FullPivHouseholderQR<Eigen::MatrixXd> qrJT;
Eigen::SparseQR<Eigen::SparseMatrix<double>, Eigen::COLAMDOrdering<int> > SqrJT;
if(qrAlgorithm==EigenDenseQR){
if (J.rows() > 0) {
qrJT=Eigen::FullPivHouseholderQR<Eigen::MatrixXd>(J.topRows(count).transpose());
Eigen::MatrixXd Q = qrJT.matrixQ ();
paramsNum = qrJT.rows();
constrNum = qrJT.cols();
qrJT.setThreshold(1e-13);
rank = qrJT.rank();
Eigen::MatrixXd R;
if (constrNum >= paramsNum)
R = qrJT.matrixQR().triangularView<Eigen::Upper>();
else
R = qrJT.matrixQR().topRows(constrNum)
.triangularView<Eigen::Upper>();
if (constrNum >= paramsNum)
R = qrJT.matrixQR().triangularView<Eigen::Upper>();
else
R = qrJT.matrixQR().topRows(constrNum)
.triangularView<Eigen::Upper>();
}
}
else if(qrAlgorithm==EigenSparseQR){
if (SJ.rows() > 0) {
SqrJT=Eigen::SparseQR<Eigen::SparseMatrix<double>, Eigen::COLAMDOrdering<int> >(SJ.topRows(count).transpose());
// Do not ask for Q Matrix!!
// At Eigen 3.2 still has a bug that this only works for square matrices
// if enabled it will crash
//Eigen::SparseMatrix<double> Q = qrJT.matrixQ();
//qrJT.matrixQ().evalTo(Q);
paramsNum = SqrJT.rows();
constrNum = SqrJT.cols();
SqrJT.setPivotThreshold(1e-13);
rank = SqrJT.rank();
if (constrNum >= paramsNum)
R = SqrJT.matrixR().triangularView<Eigen::Upper>();
else
R = SqrJT.matrixR().topRows(constrNum)
.triangularView<Eigen::Upper>();
}
}
if(debugMode==IterationLevel) {
std::stringstream stream;
stream << (qrAlgorithm==EigenSparseQR?"EigenSparseQR":(qrAlgorithm==EigenDenseQR?"DenseQR":""));
stream << ", Threads: " << Eigen::nbThreads()
#ifdef EIGEN_VECTORIZE
<< ", Vectorization: On"
#endif
<< ", Params: " << paramsNum
<< ", Constr: " << constrNum
<< ", Rank: " << rank << "\n";
const std::string tmp = stream.str();
Base::Console().Log(tmp.c_str());
}
if (J.rows() > 0) {
#ifdef _GCS_DEBUG_SOLVER_JACOBIAN_QR_DECOMPOSITION_TRIANGULAR_MATRIX
// Debug code starts
std::stringstream stream;
@ -1676,7 +1802,7 @@ int System::diagnose()
const std::string tmp = stream.str();
Base::Console().Warning(tmp.c_str());
Base::Console().Log(tmp.c_str());
// Debug code ends
#endif
@ -1696,11 +1822,23 @@ int System::diagnose()
for (int j=rank; j < constrNum; j++) {
for (int row=0; row < rank; row++) {
if (fabs(R(row,j)) > 1e-10) {
int origCol = qrJT.colsPermutation().indices()[row];
int origCol;
if(qrAlgorithm==EigenDenseQR)
origCol=qrJT.colsPermutation().indices()[row];
else if(qrAlgorithm==EigenSparseQR)
origCol=SqrJT.colsPermutation().indices()[row];
conflictGroups[j-rank].push_back(clist[origCol]);
}
}
int origCol = qrJT.colsPermutation().indices()[j];
int origCol;
if(qrAlgorithm==EigenDenseQR)
origCol=qrJT.colsPermutation().indices()[j];
else if(qrAlgorithm==EigenSparseQR)
origCol=SqrJT.colsPermutation().indices()[j];
conflictGroups[j-rank].push_back(clist[origCol]);
}
@ -1750,13 +1888,32 @@ int System::diagnose()
clistTmp.push_back(*constr);
SubSystem *subSysTmp = new SubSystem(clistTmp, plist);
int res = solve(subSysTmp);
int res = solve(subSysTmp,true,alg,true);
if(debugMode==Minimal) {
std::string solvername;
switch (alg) {
case 0:
solvername = "BFGS";
break;
case 1: // solving with the LevenbergMarquardt solver
solvername = "LevenbergMarquardt";
break;
case 2: // solving with the BFGS solver
solvername = "DogLeg";
break;
}
Base::Console().Log("Sketcher::RedundantSolving-%s-\n",solvername.c_str());
}
if (res == Success) {
subSysTmp->applySolution();
for (std::set<Constraint *>::const_iterator constr=skipped.begin();
constr != skipped.end(); constr++) {
double err = (*constr)->error();
if (err * err < XconvergenceFine)
if (err * err < convergenceRedundant)
redundant.insert(*constr);
}
resetToReference();

View File

@ -24,17 +24,15 @@
#define PLANEGCS_GCS_H
#include "SubSystem.h"
#include <boost/concept_check.hpp>
namespace GCS
{
///////////////////////////////////////
// BFGS Solver parameters
// Other BFGS Solver parameters
///////////////////////////////////////
#define XconvergenceRough 1e-8
#define XconvergenceFine 1e-10
#define smallF 1e-20
#define MaxIterations 100 //Note that the total number of iterations allowed is MaxIterations *xLength
///////////////////////////////////////
// Solver
@ -51,6 +49,17 @@ namespace GCS
LevenbergMarquardt = 1,
DogLeg = 2
};
enum QRAlgorithm {
EigenDenseQR = 0,
EigenSparseQR = 1
};
enum DebugMode {
NoDebug = 0,
Minimal = 1,
IterationLevel = 2
};
class System
{
@ -83,9 +92,31 @@ namespace GCS
bool hasDiagnosis; // if dofs, conflictingTags, redundantTags are up to date
bool isInit; // if plists, clists, reductionmaps are up to date
int solve_BFGS(SubSystem *subsys, bool isFine);
int solve_LM(SubSystem *subsys);
int solve_DL(SubSystem *subsys);
int solve_BFGS(SubSystem *subsys, bool isFine=true, bool isRedundantsolving=false);
int solve_LM(SubSystem *subsys, bool isRedundantsolving=false);
int solve_DL(SubSystem *subsys, bool isRedundantsolving=false);
public:
int maxIter;
int maxIterRedundant;
bool sketchSizeMultiplier; // if true note that the total number of iterations allowed is MaxIterations *xLength
bool sketchSizeMultiplierRedundant;
double convergence;
double convergenceRedundant;
QRAlgorithm qrAlgorithm;
DebugMode debugMode;
double LM_eps;
double LM_eps1;
double LM_tau;
double DL_tolg;
double DL_tolx;
double DL_tolf;
double LM_epsRedundant;
double LM_eps1Redundant;
double LM_tauRedundant;
double DL_tolgRedundant;
double DL_tolxRedundant;
double DL_tolfRedundant;
public:
System();
/*System(std::vector<Constraint *> clist_);*/
@ -192,19 +223,22 @@ namespace GCS
void rescaleConstraint(int id, double coeff);
void declareUnknowns(VEC_pD &params);
void initSolution();
void initSolution(Algorithm alg=DogLeg);
int solve(bool isFine=true, Algorithm alg=DogLeg);
int solve(VEC_pD &params, bool isFine=true, Algorithm alg=DogLeg);
int solve(SubSystem *subsys, bool isFine=true, Algorithm alg=DogLeg);
int solve(SubSystem *subsysA, SubSystem *subsysB, bool isFine=true);
int solve(bool isFine=true, Algorithm alg=DogLeg, bool isRedundantsolving=false);
int solve(VEC_pD &params, bool isFine=true, Algorithm alg=DogLeg, bool isRedundantsolving=false);
int solve(SubSystem *subsys, bool isFine=true, Algorithm alg=DogLeg, bool isRedundantsolving=false);
int solve(SubSystem *subsysA, SubSystem *subsysB, bool isFine=true, bool isRedundantsolving=false);
void applySolution();
void undoSolution();
//FIXME: looks like XconvergenceFine is not the solver precision, at least in DogLeg solver.
// Note: Yes, every solver has a different way of interpreting precision
// but one has to study what is this needed for in order to decide
// what to return (this is unchanged from previous versions)
double getFinePrecision(){ return convergence;}
double getFinePrecision(){ return XconvergenceFine;}//FIXME: looks like XconvergenceFine is not the solver precision, at least in DogLeg slover.
int diagnose();
int diagnose(Algorithm alg=DogLeg);
int dofsNumber() { return hasDiagnosis ? dofs : -1; }
void getConflicting(VEC_I &conflictingOut) const
{ conflictingOut = hasDiagnosis ? conflictingTags : VEC_I(0); }

View File

@ -32,6 +32,7 @@ set(SketcherGui_MOC_HDRS
TaskSketcherCreateCommands.h
TaskSketcherGeneral.h
TaskSketcherMessages.h
TaskSketcherSolverAdvanced.h
TaskSketcherValidation.h
TaskDlgEditSketch.h
SketchOrientationDialog.h
@ -48,6 +49,7 @@ set(SketcherGui_UIC_SRCS
TaskSketcherElements.ui
TaskSketcherGeneral.ui
TaskSketcherMessages.ui
TaskSketcherSolverAdvanced.ui
TaskSketcherValidation.ui
InsertDatum.ui
SketchOrientationDialog.ui
@ -90,6 +92,9 @@ SET(SketcherGui_SRCS
TaskSketcherMessages.ui
TaskSketcherMessages.cpp
TaskSketcherMessages.h
TaskSketcherSolverAdvanced.ui
TaskSketcherSolverAdvanced.h
TaskSketcherSolverAdvanced.cpp
TaskSketcherValidation.ui
TaskSketcherValidation.cpp
TaskSketcherValidation.h

View File

@ -47,11 +47,25 @@ TaskDlgEditSketch::TaskDlgEditSketch(ViewProviderSketch *sketchView)
Elements = new TaskSketcherElements(sketchView);
General = new TaskSketcherGeneral(sketchView);
Messages = new TaskSketcherMessages(sketchView);
SolverAdvanced = new TaskSketcherSolverAdvanced(sketchView);
Content.push_back(Messages);
Content.push_back(SolverAdvanced);
Content.push_back(General);
Content.push_back(Constraints);
Content.push_back(Elements);
ParameterGrp::handle hGrp = App::GetApplication().GetParameterGroupByPath("User parameter:BaseApp/Preferences/Mod/Sketcher");
if( !hGrp->GetBool("ShowMessagesWidget",true))
Messages->hideGroupBox();
if( !hGrp->GetBool("ShowSolverAdvancedWidget",false))
SolverAdvanced->hideGroupBox();
if( !hGrp->GetBool("ShowEditControlWidget",false))
General->hideGroupBox();
if( !hGrp->GetBool("ShowConstraintsWidget",true))
Constraints->hideGroupBox();
if( !hGrp->GetBool("ShowElementsWidget",true))
Elements->hideGroupBox();
App::Document* document = sketchView->getObject()->getDocument();
connectUndoDocument =
@ -90,12 +104,19 @@ void TaskDlgEditSketch::clicked(int)
}
bool TaskDlgEditSketch::accept()
{
{
return true;
}
bool TaskDlgEditSketch::reject()
{
ParameterGrp::handle hGrp = App::GetApplication().GetParameterGroupByPath("User parameter:BaseApp/Preferences/Mod/Sketcher");
hGrp->SetBool("ShowMessagesWidget",Messages->isGroupVisible());
hGrp->SetBool("ShowSolverAdvancedWidget",SolverAdvanced->isGroupVisible());
hGrp->SetBool("ShowEditControlWidget",General->isGroupVisible());
hGrp->SetBool("ShowConstraintsWidget",Constraints->isGroupVisible());
hGrp->SetBool("ShowElementsWidget",Elements->isGroupVisible());
std::string document = getDocumentName(); // needed because resetEdit() deletes this instance
Gui::Command::doCommand(Gui::Command::Gui,"Gui.getDocument('%s').resetEdit()", document.c_str());
Gui::Command::doCommand(Gui::Command::Doc,"App.getDocument('%s').recompute()", document.c_str());

View File

@ -31,6 +31,7 @@
#include "TaskSketcherElements.h"
#include "TaskSketcherGeneral.h"
#include "TaskSketcherMessages.h"
#include "TaskSketcherSolverAdvanced.h"
#include <boost/signals.hpp>
typedef boost::signals::connection Connection;
@ -75,6 +76,7 @@ protected:
TaskSketcherElements *Elements;
TaskSketcherGeneral *General;
TaskSketcherMessages *Messages;
TaskSketcherSolverAdvanced *SolverAdvanced;
Connection connectUndoDocument;
Connection connectRedoDocument;
};

View File

@ -0,0 +1,466 @@
/***************************************************************************
* Copyright (c) 2015 Abdullah Tahiri <abdullah.tahiri.yo@gmail.com *
* *
* This file is part of the FreeCAD CAx development system. *
* *
* This library is free software; you can redistribute it and/or *
* modify it under the terms of the GNU Library General Public *
* License as published by the Free Software Foundation; either *
* version 2 of the License, or (at your option) any later version. *
* *
* This library is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU Library General Public License for more details. *
* *
* You should have received a copy of the GNU Library General Public *
* License along with this library; see the file COPYING.LIB. If not, *
* write to the Free Software Foundation, Inc., 59 Temple Place, *
* Suite 330, Boston, MA 02111-1307, USA *
* *
***************************************************************************/
#include "PreCompiled.h"
#ifndef _PreComp_
#endif
#include "ui_TaskSketcherSolverAdvanced.h"
#include "TaskSketcherSolverAdvanced.h"
#include <Gui/Application.h>
#include <Gui/Document.h>
#include <Gui/BitmapFactory.h>
#include <Gui/ViewProvider.h>
#include <Gui/WaitCursor.h>
#include <Gui/Selection.h>
#include <Gui/Command.h>
#include <boost/bind.hpp>
#include <Mod/Sketcher/App/SketchObject.h>
#include <QString>
#include "ViewProviderSketch.h"
using namespace SketcherGui;
using namespace Gui::TaskView;
TaskSketcherSolverAdvanced::TaskSketcherSolverAdvanced(ViewProviderSketch *sketchView)
: TaskBox(Gui::BitmapFactory().pixmap("document-new"),tr("Advanced solver control"),true, 0)
, sketchView(sketchView)
{
// we need a separate container widget to add all controls to
proxy = new QWidget(this);
ui = new Ui_TaskSketcherSolverAdvanced();
ui->setupUi(proxy);
QMetaObject::connectSlotsByName(this);
this->groupLayout()->addWidget(proxy);
ui->comboBoxDefaultSolver->onRestore();
ui->spinBoxMaxIter->onRestore();
ui->checkBoxSketchSizeMultiplier->onRestore();
ui->lineEditCovergence->onRestore();
ui->comboBoxQRMethod->onRestore();
ui->comboBoxRedundantDefaultSolver->onRestore();
ui->spinBoxRedundantSolverMaxIterations->onRestore();
ui->checkBoxRedundantSketchSizeMultiplier->onRestore();
ui->lineEditRedundantConvergence->onRestore();
ui->comboBoxDebugMode->onRestore();
updateSketchObject();
}
TaskSketcherSolverAdvanced::~TaskSketcherSolverAdvanced()
{
delete ui;
}
void TaskSketcherSolverAdvanced::updateDefaultMethodParameters(void)
{
ParameterGrp::handle hGrp = App::GetApplication().GetParameterGroupByPath("User parameter:BaseApp/Preferences/Mod/Sketcher/SolverAdvanced");
switch(ui->comboBoxDefaultSolver->currentIndex())
{
case 0: // BFGS
ui->labelSolverParam1->setText(QString::fromLatin1(""));
ui->labelSolverParam2->setText(QString::fromLatin1(""));
ui->labelSolverParam3->setText(QString::fromLatin1(""));
ui->lineEditSolverParam1->setDisabled(true);
ui->lineEditSolverParam2->setDisabled(true);
ui->lineEditSolverParam3->setDisabled(true);
break;
case 1: // LM
{
ui->labelSolverParam1->setText(QString::fromLatin1("Eps"));
ui->labelSolverParam2->setText(QString::fromLatin1("Eps1"));
ui->labelSolverParam3->setText(QString::fromLatin1("Tau"));
ui->lineEditSolverParam1->setEnabled(true);
ui->lineEditSolverParam2->setEnabled(true);
ui->lineEditSolverParam3->setEnabled(true);
double eps = hGrp->GetFloat("LM_eps",1E-10);
double eps1 = hGrp->GetFloat("LM_eps1",1E-80);
double tau = hGrp->GetFloat("LM_tau",1E-3);
ui->lineEditSolverParam1->setText(QString::number(eps).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
ui->lineEditSolverParam2->setText(QString::number(eps1).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
ui->lineEditSolverParam3->setText(QString::number(tau).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
sketchView->getSketchObject()->getSolvedSketch().setLM_eps(eps);
sketchView->getSketchObject()->getSolvedSketch().setLM_eps1(eps1);
sketchView->getSketchObject()->getSolvedSketch().setLM_tau(eps1);
break;
}
case 2: // DogLeg
{
ui->labelSolverParam1->setText(QString::fromLatin1("Tolg"));
ui->labelSolverParam2->setText(QString::fromLatin1("Tolx"));
ui->labelSolverParam3->setText(QString::fromLatin1("Tolf"));
ui->lineEditSolverParam1->setEnabled(true);
ui->lineEditSolverParam2->setEnabled(true);
ui->lineEditSolverParam3->setEnabled(true);
double tolg = hGrp->GetFloat("DL_tolg",1E-80);
double tolx = hGrp->GetFloat("DL_tolx",1E-80);
double tolf = hGrp->GetFloat("DL_tolf",1E-10);
ui->lineEditSolverParam1->setText(QString::number(tolg).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
ui->lineEditSolverParam2->setText(QString::number(tolx).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
ui->lineEditSolverParam3->setText(QString::number(tolf).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
sketchView->getSketchObject()->getSolvedSketch().setDL_tolg(tolg);
sketchView->getSketchObject()->getSolvedSketch().setDL_tolf(tolf);
sketchView->getSketchObject()->getSolvedSketch().setDL_tolx(tolx);
break;
}
}
}
void TaskSketcherSolverAdvanced::updateRedundantMethodParameters(void)
{
ParameterGrp::handle hGrp = App::GetApplication().GetParameterGroupByPath("User parameter:BaseApp/Preferences/Mod/Sketcher/SolverAdvanced");
switch(ui->comboBoxRedundantDefaultSolver->currentIndex())
{
case 0: // BFGS
ui->labelRedundantSolverParam1->setText(QString::fromLatin1(""));
ui->labelRedundantSolverParam2->setText(QString::fromLatin1(""));
ui->labelRedundantSolverParam3->setText(QString::fromLatin1(""));
ui->lineEditRedundantSolverParam1->setDisabled(true);
ui->lineEditRedundantSolverParam2->setDisabled(true);
ui->lineEditRedundantSolverParam3->setDisabled(true);
break;
case 1: // LM
{
ui->labelRedundantSolverParam1->setText(QString::fromLatin1("R.Eps"));
ui->labelRedundantSolverParam2->setText(QString::fromLatin1("R.Eps1"));
ui->labelRedundantSolverParam3->setText(QString::fromLatin1("R.Tau"));
ui->lineEditRedundantSolverParam1->setEnabled(true);
ui->lineEditRedundantSolverParam2->setEnabled(true);
ui->lineEditRedundantSolverParam3->setEnabled(true);
double eps = hGrp->GetFloat("Redundant_LM_eps",1E-10);
double eps1 = hGrp->GetFloat("Redundant_LM_eps1",1E-80);
double tau = hGrp->GetFloat("Redundant_LM_tau",1E-3);
ui->lineEditRedundantSolverParam1->setText(QString::number(eps).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
ui->lineEditRedundantSolverParam2->setText(QString::number(eps1).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
ui->lineEditRedundantSolverParam3->setText(QString::number(tau).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
sketchView->getSketchObject()->getSolvedSketch().setLM_epsRedundant(eps);
sketchView->getSketchObject()->getSolvedSketch().setLM_eps1Redundant(eps1);
sketchView->getSketchObject()->getSolvedSketch().setLM_tauRedundant(eps1);
break;
}
case 2: // DogLeg
{
ui->labelRedundantSolverParam1->setText(QString::fromLatin1("R.Tolg"));
ui->labelRedundantSolverParam2->setText(QString::fromLatin1("R.Tolx"));
ui->labelRedundantSolverParam3->setText(QString::fromLatin1("R.Tolf"));
ui->lineEditRedundantSolverParam1->setEnabled(true);
ui->lineEditRedundantSolverParam2->setEnabled(true);
ui->lineEditRedundantSolverParam3->setEnabled(true);
double tolg = hGrp->GetFloat("Redundant_DL_tolg",1E-80);
double tolx = hGrp->GetFloat("Redundant_DL_tolx",1E-80);
double tolf = hGrp->GetFloat("Redundant_DL_tolf",1E-10);
ui->lineEditRedundantSolverParam1->setText(QString::number(tolg).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
ui->lineEditRedundantSolverParam2->setText(QString::number(tolx).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
ui->lineEditRedundantSolverParam3->setText(QString::number(tolf).remove(QString::fromLatin1("+").replace(QString::fromLatin1("e0"),QString::fromLatin1("E")).toUpper()));
sketchView->getSketchObject()->getSolvedSketch().setDL_tolgRedundant(tolg);
sketchView->getSketchObject()->getSolvedSketch().setDL_tolfRedundant(tolf);
sketchView->getSketchObject()->getSolvedSketch().setDL_tolxRedundant(tolx);
break;
}
}
}
void TaskSketcherSolverAdvanced::on_lineEditSolverParam1_editingFinished()
{
QString text = ui->lineEditSolverParam1->text();
double val = text.toDouble();
QString sci = QString::number(val);
sci.remove(QString::fromLatin1("+"));
sci.replace(QString::fromLatin1("e0"),QString::fromLatin1("E"));
ui->lineEditSolverParam1->setText(sci.toUpper());
switch(ui->comboBoxDefaultSolver->currentIndex())
{
case 1: // LM
{
sketchView->getSketchObject()->getSolvedSketch().setLM_eps(val);
ui->lineEditSolverParam1->setEntryName("LM_eps");
ui->lineEditSolverParam1->onSave();
break;
}
case 2: // DogLeg
{
sketchView->getSketchObject()->getSolvedSketch().setDL_tolg(val);
ui->lineEditSolverParam1->setEntryName("DL_tolg");
ui->lineEditSolverParam1->onSave();
break;
}
}
}
void TaskSketcherSolverAdvanced::on_lineEditRedundantSolverParam1_editingFinished()
{
QString text = ui->lineEditRedundantSolverParam1->text();
double val = text.toDouble();
QString sci = QString::number(val);
sci.remove(QString::fromLatin1("+"));
sci.replace(QString::fromLatin1("e0"),QString::fromLatin1("E"));
ui->lineEditRedundantSolverParam1->setText(sci.toUpper());
switch(ui->comboBoxDefaultSolver->currentIndex())
{
case 1: // LM
{
sketchView->getSketchObject()->getSolvedSketch().setLM_epsRedundant(val);
ui->lineEditRedundantSolverParam1->setEntryName("Redundant_LM_eps");
ui->lineEditRedundantSolverParam1->onSave();
break;
}
case 2: // DogLeg
{
sketchView->getSketchObject()->getSolvedSketch().setDL_tolgRedundant(val);
ui->lineEditRedundantSolverParam1->setEntryName("Redundant_DL_tolg");
ui->lineEditRedundantSolverParam1->onSave();
break;
}
}
}
void TaskSketcherSolverAdvanced::on_lineEditSolverParam2_editingFinished()
{
QString text = ui->lineEditSolverParam2->text();
double val = text.toDouble();
QString sci = QString::number(val);
sci.remove(QString::fromLatin1("+"));
sci.replace(QString::fromLatin1("e0"),QString::fromLatin1("E"));
ui->lineEditSolverParam2->setText(sci.toUpper());
switch(ui->comboBoxDefaultSolver->currentIndex())
{
case 1: // LM
{
sketchView->getSketchObject()->getSolvedSketch().setLM_eps1(val);
ui->lineEditSolverParam2->setEntryName("LM_eps1");
ui->lineEditSolverParam2->onSave();
break;
}
case 2: // DogLeg
{
sketchView->getSketchObject()->getSolvedSketch().setDL_tolx(val);
ui->lineEditSolverParam2->setEntryName("DL_tolx");
ui->lineEditSolverParam2->onSave();
break;
}
}
}
void TaskSketcherSolverAdvanced::on_lineEditRedundantSolverParam2_editingFinished()
{
QString text = ui->lineEditRedundantSolverParam2->text();
double val = text.toDouble();
QString sci = QString::number(val);
sci.remove(QString::fromLatin1("+"));
sci.replace(QString::fromLatin1("e0"),QString::fromLatin1("E"));
ui->lineEditRedundantSolverParam2->setText(sci.toUpper());
switch(ui->comboBoxDefaultSolver->currentIndex())
{
case 1: // LM
{
sketchView->getSketchObject()->getSolvedSketch().setLM_eps1Redundant(val);
ui->lineEditRedundantSolverParam2->setEntryName("Redundant_LM_eps1");
ui->lineEditRedundantSolverParam2->onSave();
break;
}
case 2: // DogLeg
{
sketchView->getSketchObject()->getSolvedSketch().setDL_tolxRedundant(val);
ui->lineEditRedundantSolverParam2->setEntryName("Redundant_DL_tolx");
ui->lineEditRedundantSolverParam2->onSave();
break;
}
}
}
void TaskSketcherSolverAdvanced::on_lineEditSolverParam3_editingFinished()
{
QString text = ui->lineEditSolverParam3->text();
double val = text.toDouble();
QString sci = QString::number(val);
sci.remove(QString::fromLatin1("+"));
sci.replace(QString::fromLatin1("e0"),QString::fromLatin1("E"));
ui->lineEditSolverParam3->setText(sci.toUpper());
switch(ui->comboBoxDefaultSolver->currentIndex())
{
case 1: // LM
{
sketchView->getSketchObject()->getSolvedSketch().setLM_tau(val);
ui->lineEditSolverParam3->setEntryName("LM_tau");
ui->lineEditSolverParam3->onSave();
break;
}
case 2: // DogLeg
{
sketchView->getSketchObject()->getSolvedSketch().setDL_tolf(val);
ui->lineEditSolverParam3->setEntryName("DL_tolf");
ui->lineEditSolverParam3->onSave();
break;
}
}
}
void TaskSketcherSolverAdvanced::on_lineEditRedundantSolverParam3_editingFinished()
{
QString text = ui->lineEditRedundantSolverParam3->text();
double val = text.toDouble();
QString sci = QString::number(val);
sci.remove(QString::fromLatin1("+"));
sci.replace(QString::fromLatin1("e0"),QString::fromLatin1("E"));
ui->lineEditRedundantSolverParam3->setText(sci.toUpper());
switch(ui->comboBoxDefaultSolver->currentIndex())
{
case 1: // LM
{
sketchView->getSketchObject()->getSolvedSketch().setLM_tauRedundant(val);
ui->lineEditRedundantSolverParam3->setEntryName("Redundant_LM_tau");
ui->lineEditRedundantSolverParam3->onSave();
break;
}
case 2: // DogLeg
{
sketchView->getSketchObject()->getSolvedSketch().setDL_tolfRedundant(val);
ui->lineEditRedundantSolverParam3->setEntryName("Redundant_DL_tolf");
ui->lineEditRedundantSolverParam3->onSave();
break;
}
}
}
void TaskSketcherSolverAdvanced::on_comboBoxDefaultSolver_currentIndexChanged(int index)
{
ui->comboBoxDefaultSolver->onSave();
sketchView->getSketchObject()->getSolvedSketch().defaultSolver=(GCS::Algorithm) index;
updateDefaultMethodParameters();
}
void TaskSketcherSolverAdvanced::on_spinBoxMaxIter_valueChanged(int i)
{
ui->spinBoxMaxIter->onSave();
sketchView->getSketchObject()->getSolvedSketch().setMaxIter(i);
}
void TaskSketcherSolverAdvanced::on_checkBoxSketchSizeMultiplier_stateChanged(int state)
{
if(state==Qt::Checked) {
ui->checkBoxSketchSizeMultiplier->onSave();
sketchView->getSketchObject()->getSolvedSketch().setSketchSizeMultiplier(true);
}
else if (state==Qt::Unchecked) {
ui->checkBoxSketchSizeMultiplier->onSave();
sketchView->getSketchObject()->getSolvedSketch().setSketchSizeMultiplier(false);
}
}
void TaskSketcherSolverAdvanced::on_lineEditCovergence_editingFinished()
{
QString text = ui->lineEditCovergence->text();
double val = text.toDouble();
QString sci = QString::number(val);
sci.remove(QString::fromLatin1("+"));
sci.replace(QString::fromLatin1("e0"),QString::fromLatin1("E"));
ui->lineEditCovergence->setText(sci.toUpper());
ui->lineEditCovergence->onSave();
sketchView->getSketchObject()->getSolvedSketch().setConvergence(val);
}
void TaskSketcherSolverAdvanced::on_lineEditRedundantConvergence_editingFinished()
{
QString text = ui->lineEditRedundantConvergence->text();
double val = text.toDouble();
QString sci = QString::number(val);
sci.remove(QString::fromLatin1("+"));
sci.replace(QString::fromLatin1("e0"),QString::fromLatin1("E"));
ui->lineEditRedundantConvergence->setText(sci.toUpper());
ui->lineEditRedundantConvergence->onSave();
sketchView->getSketchObject()->getSolvedSketch().setConvergenceRedundant(val);
}
void TaskSketcherSolverAdvanced::on_comboBoxQRMethod_currentIndexChanged(int index)
{
sketchView->getSketchObject()->getSolvedSketch().setQRAlgorithm((GCS::QRAlgorithm) index);
ui->comboBoxQRMethod->onSave();
}
void TaskSketcherSolverAdvanced::on_comboBoxRedundantDefaultSolver_currentIndexChanged(int index)
{
ui->comboBoxRedundantDefaultSolver->onSave();
sketchView->getSketchObject()->getSolvedSketch().defaultSolverRedundant=(GCS::Algorithm) index;
updateRedundantMethodParameters();
}
void TaskSketcherSolverAdvanced::on_spinBoxRedundantSolverMaxIterations_valueChanged(int i)
{
ui->spinBoxRedundantSolverMaxIterations->onSave();
sketchView->getSketchObject()->getSolvedSketch().setMaxIterRedundant(i);
}
void TaskSketcherSolverAdvanced::on_checkBoxRedundantSketchSizeMultiplier_stateChanged(int state)
{
if(state==Qt::Checked) {
ui->checkBoxRedundantSketchSizeMultiplier->onSave();
sketchView->getSketchObject()->getSolvedSketch().setSketchSizeMultiplierRedundant(true);
}
else if (state==Qt::Unchecked) {
ui->checkBoxRedundantSketchSizeMultiplier->onSave();
sketchView->getSketchObject()->getSolvedSketch().setSketchSizeMultiplierRedundant(true);
}
}
void TaskSketcherSolverAdvanced::on_comboBoxDebugMode_currentIndexChanged(int index)
{
ui->comboBoxDebugMode->onSave();
sketchView->getSketchObject()->getSolvedSketch().setDebugMode((GCS::DebugMode) index);
}
void TaskSketcherSolverAdvanced::updateSketchObject(void)
{
sketchView->getSketchObject()->getSolvedSketch().setDebugMode((GCS::DebugMode) ui->comboBoxDebugMode->currentIndex());
sketchView->getSketchObject()->getSolvedSketch().setSketchSizeMultiplierRedundant(ui->checkBoxRedundantSketchSizeMultiplier->isChecked());
sketchView->getSketchObject()->getSolvedSketch().setMaxIterRedundant(ui->spinBoxRedundantSolverMaxIterations->value());
sketchView->getSketchObject()->getSolvedSketch().defaultSolverRedundant=(GCS::Algorithm) ui->comboBoxRedundantDefaultSolver->currentIndex();
sketchView->getSketchObject()->getSolvedSketch().setQRAlgorithm((GCS::QRAlgorithm) ui->comboBoxQRMethod->currentIndex());
sketchView->getSketchObject()->getSolvedSketch().setConvergenceRedundant(ui->lineEditRedundantConvergence->text().toDouble());
sketchView->getSketchObject()->getSolvedSketch().setConvergence(ui->lineEditCovergence->text().toDouble());
sketchView->getSketchObject()->getSolvedSketch().setSketchSizeMultiplier(ui->checkBoxSketchSizeMultiplier->isChecked());
sketchView->getSketchObject()->getSolvedSketch().setMaxIter(ui->spinBoxMaxIter->value());
sketchView->getSketchObject()->getSolvedSketch().defaultSolver=(GCS::Algorithm) ui->comboBoxDefaultSolver->currentIndex();
updateDefaultMethodParameters();
updateRedundantMethodParameters();
}
#include "moc_TaskSketcherSolverAdvanced.cpp"

View File

@ -0,0 +1,81 @@
/***************************************************************************
* Copyright (c) 2015 Abdullah Tahiri <abdullah.tahiri.yo@gmail.com *
* *
* This file is part of the FreeCAD CAx development system. *
* *
* This library is free software; you can redistribute it and/or *
* modify it under the terms of the GNU Library General Public *
* License as published by the Free Software Foundation; either *
* version 2 of the License, or (at your option) any later version. *
* *
* This library is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU Library General Public License for more details. *
* *
* You should have received a copy of the GNU Library General Public *
* License along with this library; see the file COPYING.LIB. If not, *
* write to the Free Software Foundation, Inc., 59 Temple Place, *
* Suite 330, Boston, MA 02111-1307, USA *
* *
***************************************************************************/
#ifndef GUI_TASKVIEW_TaskSketcherSolverAdvanced_H
#define GUI_TASKVIEW_TaskSketcherSolverAdvanced_H
#include <Gui/TaskView/TaskView.h>
#include <Gui/Selection.h>
#include <boost/signals.hpp>
class Ui_TaskSketcherSolverAdvanced;
namespace App {
class Property;
}
namespace SketcherGui {
class ViewProviderSketch;
class TaskSketcherSolverAdvanced : public Gui::TaskView::TaskBox
{
Q_OBJECT
public:
TaskSketcherSolverAdvanced(ViewProviderSketch *sketchView);
~TaskSketcherSolverAdvanced();
private Q_SLOTS:
void on_comboBoxDefaultSolver_currentIndexChanged(int index);
void on_spinBoxMaxIter_valueChanged(int i);
void on_checkBoxSketchSizeMultiplier_stateChanged(int state);
void on_lineEditCovergence_editingFinished();
void on_comboBoxQRMethod_currentIndexChanged(int index);
void on_comboBoxRedundantDefaultSolver_currentIndexChanged(int index);
void on_lineEditRedundantConvergence_editingFinished();
void on_spinBoxRedundantSolverMaxIterations_valueChanged(int i);
void on_checkBoxRedundantSketchSizeMultiplier_stateChanged(int state);
void on_comboBoxDebugMode_currentIndexChanged(int index);
void on_lineEditSolverParam1_editingFinished();
void on_lineEditRedundantSolverParam1_editingFinished();
void on_lineEditSolverParam2_editingFinished();
void on_lineEditRedundantSolverParam2_editingFinished();
void on_lineEditSolverParam3_editingFinished();
void on_lineEditRedundantSolverParam3_editingFinished();
protected:
void updateDefaultMethodParameters(void);
void updateRedundantMethodParameters(void);
void updateSketchObject(void);
protected:
ViewProviderSketch *sketchView;
private:
QWidget* proxy;
Ui_TaskSketcherSolverAdvanced* ui;
};
} //namespace SketcherGui
#endif // GUI_TASKVIEW_TaskSketcherSolverAdvanced_H

View File

@ -0,0 +1,545 @@
<?xml version="1.0" encoding="UTF-8"?>
<ui version="4.0">
<class>TaskSketcherSolverAdvanced</class>
<widget class="QWidget" name="TaskSketcherSolverAdvanced">
<property name="geometry">
<rect>
<x>0</x>
<y>0</y>
<width>326</width>
<height>619</height>
</rect>
</property>
<property name="windowTitle">
<string>Form</string>
</property>
<layout class="QVBoxLayout" name="verticalLayout">
<item>
<layout class="QHBoxLayout" name="horizontalLayout_4">
<item>
<widget class="QLabel" name="labelDefaultSolver">
<property name="toolTip">
<string>Default algorithm used for Sketch solving</string>
</property>
<property name="text">
<string>Default Solver:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefComboBox" name="comboBoxDefaultSolver">
<property name="currentIndex">
<number>2</number>
</property>
<property name="prefEntry" stdset="0">
<cstring>DefaultSolver</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
<item>
<property name="text">
<string>BFGS</string>
</property>
</item>
<item>
<property name="text">
<string>LevenbergMarquardt</string>
</property>
</item>
<item>
<property name="text">
<string>DogLeg</string>
</property>
</item>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_2">
<item>
<widget class="QLabel" name="labelMaxIter">
<property name="toolTip">
<string>Maximum number of iterations of the default algorithm</string>
</property>
<property name="text">
<string>Maximum Iterations:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefSpinBox" name="spinBoxMaxIter">
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="maximum">
<number>999</number>
</property>
<property name="value">
<number>100</number>
</property>
<property name="prefEntry" stdset="0">
<cstring>MaxIter</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_3">
<item>
<widget class="QLabel" name="labelSketchSizeMultiplier">
<property name="toolTip">
<string>If selected, the Maximum iterations value is multiplied by the sketch size</string>
</property>
<property name="text">
<string>Sketch size multiplier:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefCheckBox" name="checkBoxSketchSizeMultiplier">
<property name="sizePolicy">
<sizepolicy hsizetype="Minimum" vsizetype="Fixed">
<horstretch>0</horstretch>
<verstretch>0</verstretch>
</sizepolicy>
</property>
<property name="layoutDirection">
<enum>Qt::RightToLeft</enum>
</property>
<property name="text">
<string/>
</property>
<property name="prefEntry" stdset="0">
<cstring>SketchSizeMultiplier</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_9">
<item>
<widget class="QLabel" name="labelSolverConvergence">
<property name="toolTip">
<string>Error threshold under which convergence is reached</string>
</property>
<property name="text">
<string>Convergence:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefLineEdit" name="lineEditCovergence">
<property name="layoutDirection">
<enum>Qt::LeftToRight</enum>
</property>
<property name="text">
<string notr="true">1E-10</string>
</property>
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="prefEntry" stdset="0">
<cstring>Covergence</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_10">
<item>
<widget class="QLabel" name="labelSolverParam1">
<property name="text">
<string>Param1</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefLineEdit" name="lineEditSolverParam1">
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="prefEntry" stdset="0">
<cstring>param</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_11">
<item>
<widget class="QLabel" name="labelSolverParam2">
<property name="text">
<string>Param2</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefLineEdit" name="lineEditSolverParam2">
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="prefEntry" stdset="0">
<cstring>param</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_12">
<item>
<widget class="QLabel" name="labelSolverParam3">
<property name="text">
<string>Param3</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefLineEdit" name="lineEditSolverParam3">
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="prefEntry" stdset="0">
<cstring>param</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout">
<item>
<widget class="QLabel" name="labelQRAlgorithm">
<property name="toolTip">
<string>Algorithm used for the rank revealing QR decomposition</string>
</property>
<property name="text">
<string>QR Algorithm:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefComboBox" name="comboBoxQRMethod">
<property name="currentIndex">
<number>1</number>
</property>
<property name="prefEntry" stdset="0">
<cstring>QRMethod</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
<item>
<property name="text">
<string>Eigen Dense QR</string>
</property>
</item>
<item>
<property name="text">
<string>Eigen Sparse QR</string>
</property>
</item>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_5">
<item>
<widget class="QLabel" name="labelRedundantSolver">
<property name="toolTip">
<string>Solving algorithm used for determination of Redundant constraints</string>
</property>
<property name="text">
<string>Redundant Solver:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefComboBox" name="comboBoxRedundantDefaultSolver">
<property name="prefEntry" stdset="0">
<cstring>RedundantDefaultSolver</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
<item>
<property name="text">
<string>BFGS</string>
</property>
</item>
<item>
<property name="text">
<string>LevenbergMarquardt</string>
</property>
</item>
<item>
<property name="text">
<string>DogLeg</string>
</property>
</item>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_6">
<item>
<widget class="QLabel" name="labelRedundantSolverMaxIterations">
<property name="toolTip">
<string>Maximum number of iterations of the solver used for determination of Redundant constraints</string>
</property>
<property name="text">
<string>Red. Max Iterations:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefSpinBox" name="spinBoxRedundantSolverMaxIterations">
<property name="toolTip">
<string/>
</property>
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="maximum">
<number>999</number>
</property>
<property name="value">
<number>100</number>
</property>
<property name="prefEntry" stdset="0">
<cstring>RedundantSolverMaxIterations</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_7">
<item>
<widget class="QLabel" name="labelRedundantSketchSizeMultiplier">
<property name="toolTip">
<string>If selected, the Maximum iterations value for the redundant algorithm is multiplied by the sketch size</string>
</property>
<property name="text">
<string>Red. Sketch size multiplier:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefCheckBox" name="checkBoxRedundantSketchSizeMultiplier">
<property name="layoutDirection">
<enum>Qt::RightToLeft</enum>
</property>
<property name="text">
<string/>
</property>
<property name="prefEntry" stdset="0">
<cstring>RedundantSketchSizeMultiplier</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_13">
<item>
<widget class="QLabel" name="labelRedundantConvergence">
<property name="toolTip">
<string>Error threshold under which convergence is reached for the solving of redundant constraints</string>
</property>
<property name="text">
<string>Red. Convergence</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefLineEdit" name="lineEditRedundantConvergence">
<property name="text">
<string>1E-10</string>
</property>
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="prefEntry" stdset="0">
<cstring>RedundantConvergence</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_14">
<item>
<widget class="QLabel" name="labelRedundantSolverParam1">
<property name="text">
<string>Red. Param1</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefLineEdit" name="lineEditRedundantSolverParam1">
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="prefEntry" stdset="0">
<cstring>param</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_15">
<item>
<widget class="QLabel" name="labelRedundantSolverParam2">
<property name="text">
<string>Red. Param2</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefLineEdit" name="lineEditRedundantSolverParam2">
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="prefEntry" stdset="0">
<cstring>param</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_16">
<item>
<widget class="QLabel" name="labelRedundantSolverParam3">
<property name="text">
<string>Red. Param3</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefLineEdit" name="lineEditRedundantSolverParam3">
<property name="alignment">
<set>Qt::AlignRight|Qt::AlignTrailing|Qt::AlignVCenter</set>
</property>
<property name="prefEntry" stdset="0">
<cstring>param</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
</widget>
</item>
</layout>
</item>
<item>
<layout class="QHBoxLayout" name="horizontalLayout_8">
<item>
<widget class="QLabel" name="labelDebugMode">
<property name="toolTip">
<string>Degree of verbosity of the debug output to the console</string>
</property>
<property name="text">
<string>Console Debug mode:</string>
</property>
</widget>
</item>
<item>
<widget class="Gui::PrefComboBox" name="comboBoxDebugMode">
<property name="currentIndex">
<number>1</number>
</property>
<property name="prefEntry" stdset="0">
<cstring>DebugMode</cstring>
</property>
<property name="prefPath" stdset="0">
<cstring>Mod/Sketcher/SolverAdvanced</cstring>
</property>
<item>
<property name="text">
<string>None</string>
</property>
</item>
<item>
<property name="text">
<string>Minimum</string>
</property>
</item>
<item>
<property name="text">
<string>Iteration Level</string>
</property>
</item>
</widget>
</item>
</layout>
</item>
</layout>
</widget>
<customwidgets>
<customwidget>
<class>Gui::PrefLineEdit</class>
<extends>QLineEdit</extends>
<header>Gui/PrefWidgets.h</header>
</customwidget>
<customwidget>
<class>Gui::PrefComboBox</class>
<extends>QComboBox</extends>
<header>Gui/PrefWidgets.h</header>
</customwidget>
<customwidget>
<class>Gui::PrefCheckBox</class>
<extends>QCheckBox</extends>
<header>Gui/PrefWidgets.h</header>
</customwidget>
<customwidget>
<class>Gui::PrefSpinBox</class>
<extends>QSpinBox</extends>
<header>Gui/PrefWidgets.h</header>
</customwidget>
</customwidgets>
<resources/>
<connections/>
</ui>