mcpele
1.0.0
The Monte Carlo Python Energy Landscape Explorer
|
00001 #ifndef _PELE_MODIFIED_FIRE_H__ 00002 #define _PELE_MODIFIED_FIRE_H__ 00003 00004 #include "base_potential.h" 00005 #include "array.h" 00006 #include "optimizer.h" 00007 00008 namespace pele{ 00044 class MODIFIED_FIRE : public GradientOptimizer{ 00045 private : 00046 double _dtstart, _dt, _dtmax, _maxstep, _Nmin, _finc, _fdec, _fa, _astart, _a, _fold, _ifnorm, _vnorm; 00047 pele::Array<double> _v, _dx, _xold, _gold; 00048 size_t _fire_iter_number, _N; 00049 bool _stepback; 00050 inline void _ForwardEuler_integration(); 00051 inline void _VelocityVerlet_integration(); 00052 public : 00053 00057 MODIFIED_FIRE(std::shared_ptr<pele::BasePotential> potential, pele::Array<double>& x0, 00058 double dtstart, double dtmax, double maxstep, size_t Nmin=5, 00059 double finc=1.1, double fdec=0.5, double fa=0.99, double astart=0.1, 00060 double tol=1e-4, bool stepback=true); 00064 virtual ~MODIFIED_FIRE() {} 00065 00070 void one_iteration(); 00071 00076 void initialize_func_gradient(); 00077 00078 void set_func_gradient(double f, Array<double> grad); 00079 00080 inline void reset(Array<double> &x0); 00081 00082 }; 00083 00084 inline void MODIFIED_FIRE::reset(Array<double> &x0) 00085 { 00086 //arrays are already wrapped by the integrator, must not wrap them again, just update their values, dont's use array.copy() 00087 // or will wrap a copy to the array 00088 if (!func_initialized_) { 00089 initialize_func_gradient(); 00090 } 00091 iter_number_ = 0; 00092 x_.assign(x0); 00093 f_ = potential_->get_energy_gradient(x_, g_); 00094 nfev_ = 1; 00095 //fire specific 00096 _fire_iter_number = 0; 00097 _dt = _dtstart; 00098 _a = _astart; 00099 _fold = f_; 00100 _xold.assign(x_); 00101 _gold.assign(g_); 00102 for(size_t k=0; k<x_.size();++k){ 00103 _v[k] = -g_[k]*_dt; 00104 } 00105 _ifnorm = 1./norm(g_); 00106 _vnorm = norm(_v); 00107 rms_ = 1. / (_ifnorm*sqrt(_N)); 00108 } 00109 00110 inline void MODIFIED_FIRE::_VelocityVerlet_integration() 00111 { 00112 /* the minuses in the following expressions are due to the fact that 00113 * the gradients rather than the forces appear in the expression 00114 */ 00115 for(size_t i=0; i<_N; ++i) { 00116 _v[i] -= 0.5 * _dt * (_gold[i] + g_[i]); //update velocity assumes all masses 1 00117 _dx[i] = _dt * (_v[i] - 0.5 * _dt * g_[i]); //build displacement vector, assumes all masses 1 00118 } 00119 _gold.assign(g_); //save gradient as old g 00120 double normdx = norm(_dx); 00121 00122 if(normdx > _maxstep){ 00123 _dx *= (_maxstep / normdx); //resize displacement vector is greater than _maxstep 00124 } 00125 00126 x_ += _dx; 00127 00128 f_ = potential_->get_energy_gradient(x_, g_); //update gradient 00129 } 00130 00131 inline void MODIFIED_FIRE::_ForwardEuler_integration() 00132 { 00133 /* the minuses in the following expressions are due to the fact that 00134 * the gradients rather than the forces appear in the expression 00135 */ 00136 for(size_t i=0; i<_N; ++i) { //this was after get_energy_gradient, moved for testing 00137 _v[i] -= _dt * g_[i]; //update velocity, assumes all masses are 1 00138 _dx[i] = _dt * _v[i]; //build displacement vector 00139 } 00140 00141 _gold.assign(g_); //save gradient as old g 00142 double normdx = norm(_dx); 00143 00144 if(normdx > _maxstep){ 00145 _dx *= (_maxstep / normdx); //resize displacement vector is greater than _maxstep 00146 } 00147 00148 x_ += _dx; 00149 00150 f_ = potential_->get_energy_gradient(x_, g_); //update gradient 00151 } 00152 00153 inline void MODIFIED_FIRE::one_iteration() 00154 { 00155 nfev_ += 1; 00156 iter_number_ += 1; 00157 _fire_iter_number += 1; //this is different from iter_number_ which does not get reset 00158 00159 //save old configuration in case next step has P < 0 00160 _fold = f_; //set f_ old before integration step 00161 _xold.assign(x_); //save x as xold, gold saved in integrator (because velocity verlet needs it, if vv is used) 00162 00163 /*equation written in this conditional statement _v = (1- _a)*_v + _a * funit * vnorm*/ 00164 00165 for (size_t i=0; i < _N; ++i) { 00166 _v[i] = (1. - _a) * _v[i] - _a * g_[i] * _ifnorm * _vnorm; 00167 } 00168 00169 /*run MD*/ 00170 this->_ForwardEuler_integration(); 00171 //this->_VelocityVerlet_integration(); 00172 00173 double P = -1 * dot(_v,g_); 00174 00175 if (P > 0) { 00176 if (_fire_iter_number > _Nmin) { 00177 _dt = std::min(_dt* _finc, _dtmax); 00178 _a *= _fa; 00179 } 00180 00181 _ifnorm = 1./norm(g_); 00182 _vnorm = norm(_v); 00183 rms_ = 1. / (_ifnorm * sqrt(_N)); //update rms 00184 } else { 00185 _dt *= _fdec; 00186 _a = _astart; 00187 _fire_iter_number = 0; 00188 _v.assign(0); 00189 00190 //reset position and gradient to the one before the step (core of modified fire) reset velocity to initial (0) 00191 if (_stepback == true) { 00192 f_ = _fold; 00193 //reset position and gradient to the one before the step (core of modified fire) reset velocity to initial (0) 00194 x_.assign(_xold); 00195 g_.assign(_gold); 00196 } 00197 if (verbosity_ > 1) { 00198 std::cout << "warning: step direction was uphill. inverting\n"; 00199 } 00200 } 00201 00202 // print some status information 00203 if ((iprint_ > 0) && (iter_number_ % iprint_ == 0)){ 00204 std::cout << "fire: " << iter_number_ 00205 << " fire_iter_number " << _fire_iter_number 00206 << " dt " << _dt 00207 << " a " << _a 00208 << " P " << P 00209 << " vnorm " << _vnorm 00210 << " E " << f_ 00211 << " rms " << rms_ 00212 << " nfev " << nfev_ << "\n"; 00213 } 00214 } 00215 } 00216 00217 #endif