/usr/include/rdkit/Numerics/Optimizer/BFGSOpt.h is in librdkit-dev 201503-3.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 | //
// Copyright (C) 2004-2008 Greg Landrum and Rational Discovery LLC
//
// @@ All Rights Reserved @@
// This file is part of the RDKit.
// The contents are covered by the terms of the BSD license
// which is included in the file license.txt, found at the root
// of the RDKit source tree.
//
#include <math.h>
#include <RDGeneral/Invariant.h>
#include <cstring>
#include <algorithm>
namespace BFGSOpt {
const double FUNCTOL=1e-4; //!< Default tolerance for function convergence in the minimizer
const double MOVETOL=1e-7; //!< Default tolerance for x changes in the minimizer
const int MAXITS=200; //!< Default maximum number of iterations
const double EPS=3e-8; //!< Default gradient tolerance in the minimizer
const double TOLX=4.*EPS; //!< Default direction vector tolerance in the minimizer
const double MAXSTEP=100.0; //!< Default maximim step size in the minimizer
//! Do a Quasi-Newton minimization along a line.
/*!
See Numerical Recipes in C, Section 9.7 for a description of the algorithm.
\param dim the dimensionality of the space.
\param oldPt the current position, as an array.
\param oldVal the current function value.
\param grad the value of the function gradient at oldPt
\param dir the minimization direction
\param newPt used to return the final position
\param newVal used to return the final function value
\param func the function to minimize
\param maxStep the maximum allowable step size
\param resCode used to return the results of the search.
Possible values for resCode are on return are:
- 0: success
- 1: the stepsize got too small. This probably indicates success.
- -1: the direction is bad (orthogonal to the gradient)
*/
template <typename EnergyFunctor>
void linearSearch(unsigned int dim,double *oldPt,double oldVal,
double *grad,double *dir,double *newPt,
double &newVal,
EnergyFunctor func,
double maxStep,int &resCode){
PRECONDITION(oldPt,"bad input array");
PRECONDITION(grad,"bad input array");
PRECONDITION(dir,"bad input array");
PRECONDITION(newPt,"bad input array");
const unsigned int MAX_ITER_LINEAR_SEARCH = 1000;
double sum=0.0,slope=0.0,test=0.0,lambda=0.0;
double lambda2=0.0,lambdaMin=0.0,tmpLambda=0.0,val2=0.0;
resCode=-1;
// get the length of the direction vector:
sum=0.0;
for(unsigned int i=0;i<dim;i++)
sum +=dir[i]*dir[i];
sum=sqrt(sum);
// rescale if we're trying to move too far:
if(sum>maxStep){
for(unsigned int i=0;i<dim;i++)
dir[i] *= maxStep/sum;
}
// make sure our direction has at least some component along
// -grad
slope=0.0;
for(unsigned int i=0;i<dim;i++){
slope += dir[i]*grad[i];
}
if(slope>=0.0){
return;
}
test=0.0;
for(unsigned int i=0;i<dim;i++){
double temp=fabs(dir[i])/std::max(fabs(oldPt[i]),1.0);
if(temp>test) test=temp;
}
lambdaMin = MOVETOL/test;
lambda = 1.0;
unsigned int it = 0;
while(it < MAX_ITER_LINEAR_SEARCH){
//std::cerr << "\t" << it<<" : "<<lambda << " " << lambdaMin << std::endl;
if(lambda<lambdaMin){
// the position change is too small
resCode=1;
break;
}
for(unsigned int i=0;i<dim;i++){
newPt[i]=oldPt[i]+lambda*dir[i];
}
newVal = func(newPt);
if( newVal-oldVal <= FUNCTOL*lambda*slope ){
// we're converged on the function:
resCode=0;
return;
}
// if we made it this far, we need to backtrack:
if(it == 0){
// it's the first step:
tmpLambda = -slope / (2.0*(newVal-oldVal-slope));
} else {
double rhs1 = newVal-oldVal-lambda*slope;
double rhs2 = val2-oldVal-lambda2*slope;
double a = (rhs1/(lambda*lambda) - rhs2/(lambda2*lambda2))/
(lambda-lambda2);
double b = (-lambda2*rhs1/(lambda*lambda)+lambda*rhs2/(lambda2*lambda2))/
(lambda-lambda2);
if( a==0.0 ){
tmpLambda = -slope/(2.0*b);
} else {
double disc=b*b-3*a*slope;
if(disc<0.0){
tmpLambda = 0.5*lambda;
} else if(b<=0.0) {
tmpLambda = (-b+sqrt(disc))/(3.0*a);
} else {
tmpLambda = -slope/(b+sqrt(disc));
}
}
if( tmpLambda > 0.5*lambda ){
tmpLambda = 0.5*lambda;
}
}
lambda2 = lambda;
val2 = newVal;
lambda = std::max(tmpLambda,0.1*lambda);
++it;
}
// nothing was done
//std::cerr<<" RETURN AT END: "<<it<<" "<<resCode<<std::endl;
for(unsigned int i=0;i<dim;i++){
newPt[i]=oldPt[i];
}
}
#define CLEANUP() { delete [] grad; delete [] dGrad; delete [] hessDGrad;\
delete [] newPos; delete [] xi; delete [] invHessian; }
//! Do a BFGS minimization of a function.
/*!
See Numerical Recipes in C, Section 10.7 for a description of the algorithm.
\param dim the dimensionality of the space.
\param pos the starting position, as an array.
\param gradTol tolerance for gradient convergence
\param numIters used to return the number of iterations required
\param funcVal used to return the final function value
\param func the function to minimize
\param gradFunc calculates the gradient of func
\param funcTol tolerance for changes in the function value for convergence.
\param maxIts maximum number of iterations allowed
\return a flag indicating success (or type of failure). Possible values are:
- 0: success
- 1: too many iterations were required
*/
template <typename EnergyFunctor,typename GradientFunctor>
int minimize(unsigned int dim,double *pos,
double gradTol,
unsigned int &numIters,
double &funcVal,
EnergyFunctor func,
GradientFunctor gradFunc,
double funcTol=TOLX,
unsigned int maxIts=MAXITS){
PRECONDITION(pos,"bad input array");
PRECONDITION(gradTol>0,"bad tolerance");
double sum,maxStep,fp;
double *grad,*dGrad,*hessDGrad;
double *newPos,*xi;
double *invHessian;
grad = new double[dim];
dGrad = new double[dim];
hessDGrad = new double[dim];
newPos = new double[dim];
xi = new double[dim];
invHessian = new double[dim*dim];
// evaluate the function and gradient in our current position:
fp=func(pos);
gradFunc(pos,grad);
sum = 0.0;
memset(invHessian,0,dim*dim*sizeof(double));
for(unsigned int i=0;i<dim;i++){
unsigned int itab=i*dim;
// initialize the inverse hessian to be identity:
invHessian[itab+i]=1.0;
// the first line dir is -grad:
xi[i] = -grad[i];
sum += pos[i]*pos[i];
}
// pick a max step size:
maxStep = MAXSTEP * std::max(sqrt(sum),static_cast<double>(dim));
for(unsigned int iter=1;iter<=maxIts;iter++){
numIters=iter;
int status;
// do the line search:
linearSearch(dim,pos,fp,grad,xi,newPos,funcVal,func,maxStep,status);
CHECK_INVARIANT(status>=0,"bad direction in linearSearch");
// save the function value for the next search:
fp = funcVal;
// set the direction of this line and save the gradient:
double test=0.0;
for(unsigned int i=0;i<dim;i++){
xi[i] = newPos[i]-pos[i];
pos[i] = newPos[i];
double temp=fabs(xi[i])/std::max(fabs(pos[i]),1.0);
if(temp>test) test=temp;
dGrad[i] = grad[i];
}
//std::cerr<<" iter: "<<iter<<" "<<fp<<" "<<test<<" "<<TOLX<<std::endl;
if(test<TOLX) {
CLEANUP();
return 0;
}
// update the gradient:
double gradScale=gradFunc(pos,grad);
// is the gradient converged?
test=0.0;
double term=std::max(funcVal*gradScale,1.0);
for(unsigned int i=0;i<dim;i++){
double temp=fabs(grad[i])*std::max(fabs(pos[i]),1.0);
test=std::max(test,temp);
dGrad[i] = grad[i]-dGrad[i];
}
test /= term;
//std::cerr<<" "<<gradScale<<" "<<test<<" "<<gradTol<<std::endl;
if(test<gradTol){
CLEANUP();
return 0;
}
//for(unsigned int i=0;i<dim;i++){
// figure out how much the gradient changed:
//}
// compute hessian*dGrad:
double fac=0,fae=0,sumDGrad=0,sumXi=0;
for(unsigned int i=0;i<dim;i++){
unsigned int itab=i*dim;
hessDGrad[i] = 0.0;
for(unsigned int j=0;j<dim;j++){
hessDGrad[i] += invHessian[itab+j]*dGrad[j];
}
fac += dGrad[i]*xi[i];
fae += dGrad[i]*hessDGrad[i];
sumDGrad += dGrad[i]*dGrad[i];
sumXi += xi[i]*xi[i];
}
if(fac > sqrt(EPS*sumDGrad*sumXi)){
fac = 1.0/fac;
double fad = 1.0/fae;
for(unsigned int i=0;i<dim;i++){
dGrad[i] = fac*xi[i] - fad*hessDGrad[i];
}
for(unsigned int i=0;i<dim;i++){
unsigned int itab=i*dim;
for(unsigned int j=i;j<dim;j++){
invHessian[itab+j] += fac*xi[i]*xi[j] -
fad*hessDGrad[i]*hessDGrad[j] +
fae*dGrad[i]*dGrad[j];
invHessian[j*dim+i] = invHessian[itab+j];
}
}
}
// generate the next direction to move:
for(unsigned int i=0;i<dim;i++){
unsigned int itab=i*dim;
xi[i] = 0.0;
for(unsigned int j=0;j<dim;j++){
xi[i] -= invHessian[itab+j]*grad[j];
}
}
}
CLEANUP();
return 1;
}
}
|