Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2019 Jakub Rydzewski (jr@fizyka.umk.pl). All rights reserved. 3 : 4 : See http://www.maze-code.github.io for more information. 5 : 6 : This file is part of maze. 7 : 8 : maze is free software: you can redistribute it and/or modify it under the 9 : terms of the GNU Lesser General Public License as published by the Free 10 : Software Foundation, either version 3 of the License, or (at your option) 11 : any later version. 12 : 13 : maze is distributed in the hope that it will be useful, but WITHOUT ANY 14 : WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS 15 : FOR A PARTICULAR PURPOSE. 16 : 17 : See the GNU Lesser General Public License for more details. 18 : 19 : You should have received a copy of the GNU Lesser General Public License 20 : along with maze. If not, see <https://www.gnu.org/licenses/>. 21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */ 22 : #ifndef __PLUMED_maze_Optimizer_h 23 : #define __PLUMED_maze_Optimizer_h 24 : 25 : /** 26 : * @file Optimizer.h 27 : * 28 : * @author J. Rydzewski (jr@fizyka.umk.pl) 29 : */ 30 : 31 : #include "colvar/Colvar.h" 32 : #include "tools/Communicator.h" 33 : #include "tools/OpenMP.h" 34 : #include "tools/NeighborList.h" 35 : #include "tools/Vector.h" 36 : 37 : #include "Core.h" 38 : #include "Loss.h" 39 : 40 : #define PLUMED_OPT_INIT(ao) Action(ao), Optimizer(ao) 41 : 42 : namespace PLMD { 43 : namespace maze { 44 : 45 : /** 46 : * @ingroup INHERIT 47 : * 48 : * @class Optimizer Optimizer.h "maze/Optimizer.h" 49 : * 50 : * @brief Base class for implementing optimizers for ligand unbinding. 51 : * 52 : * An optimizer is defined as a colvar that can be passed to Optimizer_Bias. 53 : */ 54 : class Optimizer: public colvar::Colvar { 55 : public: 56 : /** 57 : * PLMD constructor. 58 : * 59 : * @param[in] ao PLMD::ActionOptions& 60 : */ 61 : explicit Optimizer(const ActionOptions&); 62 : 63 : /** 64 : * Destructor. 65 : */ 66 14 : ~Optimizer() { /* Nothing to do. */ } 67 : 68 : /** 69 : * Registers PLMD keywords. 70 : * 71 : * @param[in] keys PLMD keywords 72 : */ 73 : static void registerKeywords(Keywords& keys); 74 : 75 : /** 76 : * The pairing function needs to be overridden by a specific optimizer. 77 : * 78 : * @param[in] distance distance between a pair of atoms 79 : */ 80 : virtual double pairing(double distance) const; 81 : 82 : /** 83 : * Optimal values needed for biasing are computed by methods overridding the 84 : * optimize function. 85 : */ 86 : virtual void optimize() = 0; 87 : 88 : /** 89 : * Calculate the optimal direction of pulling. 90 : */ 91 : void calculate(); 92 : 93 : /** 94 : * Prepare the neighbor list. 95 : */ 96 : void prepare(); 97 : 98 : /** 99 : * Score a ligand-protein configuration. 100 : * 101 : * @return score 102 : */ 103 : double score(); 104 : 105 : /** 106 : * Calculate sampling radius as the minimal distance between two groups in 107 : * neighbors list. 108 : * 109 : * @return minimal distance of ligand-protein atom pairs 110 : */ 111 : double sampling_radius(); 112 : 113 : /** 114 : * Load new positions of atoms in the neighbor list. 115 : */ 116 : void update_nl(); 117 : 118 : /** 119 : * Calculate the center of mass. 120 : * 121 : * @return center of mass 122 : */ 123 : Vector center_of_mass() const; 124 : 125 : public: 126 : /** 127 : * Getters and setters. 128 : */ 129 : 130 : std::string get_label() const; 131 : void set_label(const std::string&); 132 : 133 : // Start optimizer at time = 0. 134 : void start_step_0(); 135 : 136 : // Start optimizer at time = optimizer stride. 137 : void start_step_stride(); 138 : 139 : Vector get_opt() const; 140 : void set_opt(Vector); 141 : 142 : double get_opt_value() const; 143 : void set_opt_value(double); 144 : 145 : unsigned int get_optimizer_stride() const; 146 : void set_optimizer_stride(unsigned int); 147 : 148 : bool is_pbc_on() const; 149 : void pbc_on(); 150 : void pbc_off(); 151 : 152 : unsigned int get_n_iterations() const; 153 : void set_n_iterations(unsigned int); 154 : 155 : double get_sampling_radius() const; 156 : void set_sampling_radius(double); 157 : 158 : unsigned int get_rank_openmp() const; 159 : void set_rank_openmp(unsigned int); 160 : 161 : unsigned int get_stride_openmp() const; 162 : void set_stride_openmp(unsigned int); 163 : 164 : unsigned int get_n_threads_openmp() const; 165 : void set_n_threads_openmp(unsigned int); 166 : 167 : unsigned int get_nl_stride() const; 168 : void set_nl_stride(unsigned int); 169 : 170 : double get_nl_cutofff() const; 171 : void set_nl_cutoff(double); 172 : 173 : protected: 174 : //! Optimizer label. 175 : std::string label_; 176 : 177 : //! Start either at time = 0 or time = optimizer stride. 178 : bool first_step_; 179 : 180 : //! Biasing direction. 181 : Vector opt_; 182 : 183 : //! Current loss function value. 184 : double opt_value_; 185 : 186 : //! Optimizer stride. 187 : unsigned int optimizer_stride_; 188 : 189 : //! Periodic boundary conditions. 190 : bool pbc_; 191 : 192 : //! Number of global iterations. 193 : unsigned int n_iter_; 194 : 195 : //! Sampling radius. 196 : double sampling_r_; 197 : 198 : /** 199 : * OpenMP 200 : */ 201 : unsigned int rank_; 202 : unsigned int stride_; 203 : unsigned int n_threads_; 204 : 205 : //! Neighbor list of ligand-protein atom pairs. 206 : NeighborList *neighbor_list_; 207 : 208 : //! Neighbor list cut-off. 209 : double nl_cutoff_; 210 : 211 : //! Neighbor list stride. 212 : int nl_stride_; 213 : 214 : private: 215 : bool serial_; 216 : bool validate_list_; 217 : bool first_time_; 218 : 219 : //! Pointer to the loss function. 220 : Loss* loss_; 221 : std::vector<Loss*> vec_loss_; 222 : 223 : public: 224 : /* 225 : * Pointers to PLMD components. 226 : */ 227 : 228 : //! Biased cv. 229 : Value* value_x_; 230 : Value* value_y_; 231 : Value* value_z_; 232 : 233 : //! Loss value. 234 : Value* value_action_; 235 : //! Sampling radiues value. 236 : Value* value_sampling_radius_; 237 : }; 238 : 239 : /* 240 : * Getters and setters. 241 : */ 242 : 243 : inline void Optimizer::set_nl_cutoff(double nl_cutoff) { 244 : nl_cutoff_=nl_cutoff; 245 : } 246 : 247 : inline double Optimizer::get_nl_cutofff() const { 248 : return nl_cutoff_; 249 : } 250 : 251 : inline void Optimizer::set_nl_stride(unsigned int nl_stride) { 252 : nl_stride_=nl_stride; 253 : } 254 : 255 : inline unsigned int Optimizer::get_nl_stride() const { 256 : return nl_stride_; 257 : } 258 : 259 : inline void Optimizer::set_n_threads_openmp(unsigned int n_threads) { 260 : n_threads_=n_threads; 261 : } 262 : 263 : inline unsigned int Optimizer::get_n_threads_openmp() const { 264 300 : return n_threads_; 265 : } 266 : 267 : inline void Optimizer::set_stride_openmp(unsigned int stride) { 268 : stride_=stride; 269 : } 270 : 271 : inline unsigned int Optimizer::get_stride_openmp() const { 272 : return stride_; 273 : } 274 : 275 : inline void Optimizer::set_rank_openmp(unsigned int rank) { 276 : rank_=rank; 277 : } 278 : 279 : inline unsigned int Optimizer::get_rank_openmp() const { 280 : return rank_; 281 : } 282 : 283 : inline void Optimizer::set_sampling_radius(double sampling_r) { 284 : sampling_r_=sampling_r; 285 : } 286 : 287 : inline double Optimizer::get_sampling_radius() const { 288 : return sampling_r_; 289 : } 290 : 291 : inline void Optimizer::set_n_iterations(unsigned int n_iter) { 292 : n_iter_=n_iter; 293 : } 294 : 295 : inline unsigned int Optimizer::get_n_iterations() const { 296 33 : return n_iter_; 297 : } 298 : 299 : inline void Optimizer::pbc_off() { 300 : pbc_=false; 301 : } 302 : 303 : inline void Optimizer::pbc_on() { 304 : pbc_=true; 305 : } 306 : 307 : inline bool Optimizer::is_pbc_on() const { 308 : return pbc_==true; 309 : } 310 : 311 : inline void Optimizer::set_optimizer_stride(unsigned int optimizer_stride) { 312 : optimizer_stride_=optimizer_stride; 313 : } 314 : 315 : inline unsigned int Optimizer::get_optimizer_stride() const { 316 1 : return optimizer_stride_; 317 : } 318 : 319 : inline void Optimizer::set_opt_value(double opt_value) { 320 46 : opt_value_=opt_value; 321 30 : } 322 : 323 : inline double Optimizer::get_opt_value() const { 324 : return opt_value_; 325 : } 326 : 327 : inline void Optimizer::set_opt(Vector opt) { 328 44 : opt_=opt; 329 : } 330 : 331 : inline Vector Optimizer::get_opt() const { 332 3 : return opt_; 333 : } 334 : 335 : inline void Optimizer::set_label(const std::string& label) { 336 7 : label_=label; 337 7 : } 338 : 339 : inline std::string Optimizer::get_label() const { 340 : return label_; 341 : } 342 : 343 : inline void Optimizer::start_step_0() { 344 5 : first_step_=false; 345 : } 346 : 347 : inline void Optimizer::start_step_stride() { 348 2 : first_step_=true; 349 : } 350 : 351 : } // namespace maze 352 : } // namespace PLMD 353 : 354 : #endif // __PLUMED_maze_Optimizer_h