c2992c309c6340e237e332014e24a30f44a73c7d
[namd.git] / colvars / src / colvarproxy.cpp
1 // -*- c++ -*-
2
3 // This file is part of the Collective Variables module (Colvars).
4 // The original version of Colvars and its updates are located at:
5 // https://github.com/colvars/colvars
6 // Please update all Colvars source files before making any changes.
7 // If you wish to distribute your changes, please submit them to the
8 // Colvars repository at GitHub.
9
10 #include <sstream>
11 #include <string.h>
12
13 #if defined(_OPENMP)
14 #include <omp.h>
15 #endif
16
17 #if defined(NAMD_TCL) || defined(VMDTCL)
18 #define COLVARS_TCL
19 #include <tcl.h>
20 #endif
21
22 #include "colvarmodule.h"
23 #include "colvarproxy.h"
24 #include "colvarscript.h"
25 #include "colvaratoms.h"
26
27
28
29 colvarproxy_system::colvarproxy_system()
30 {
31   reset_pbc_lattice();
32 }
33
34
35 colvarproxy_system::~colvarproxy_system() {}
36
37
38 void colvarproxy_system::add_energy(cvm::real energy) {}
39
40
41 void colvarproxy_system::request_total_force(bool yesno)
42 {
43   if (yesno == true)
44     cvm::error("Error: total forces are currently not implemented.\n",
45                COLVARS_NOT_IMPLEMENTED);
46 }
47
48
49 bool colvarproxy_system::total_forces_enabled() const
50 {
51   return false;
52 }
53
54
55 bool colvarproxy_system::total_forces_same_step() const
56 {
57   return false;
58 }
59
60
61 inline int round_to_integer(cvm::real x)
62 {
63   return std::floor(x+0.5);
64 }
65
66
67 void colvarproxy_system::update_pbc_lattice()
68 {
69   // Periodicity is assumed in all directions
70
71   if (boundaries_type == boundaries_unsupported ||
72       boundaries_type == boundaries_non_periodic) {
73     cvm::error("Error: setting PBC lattice with unsupported boundaries.\n",
74                BUG_ERROR);
75     return;
76   }
77
78   {
79     cvm::rvector const v = cvm::rvector::outer(unit_cell_y, unit_cell_z);
80     reciprocal_cell_x = v/(v*unit_cell_x);
81   }
82   {
83     cvm::rvector const v = cvm::rvector::outer(unit_cell_z, unit_cell_x);
84     reciprocal_cell_y = v/(v*unit_cell_y);
85   }
86   {
87     cvm::rvector const v = cvm::rvector::outer(unit_cell_x, unit_cell_y);
88     reciprocal_cell_z = v/(v*unit_cell_z);
89   }
90 }
91
92
93 void colvarproxy_system::reset_pbc_lattice()
94 {
95   unit_cell_x.reset();
96   unit_cell_y.reset();
97   unit_cell_z.reset();
98   reciprocal_cell_x.reset();
99   reciprocal_cell_y.reset();
100   reciprocal_cell_z.reset();
101 }
102
103
104 cvm::rvector colvarproxy_system::position_distance(cvm::atom_pos const &pos1,
105                                                    cvm::atom_pos const &pos2)
106   const
107 {
108   if (boundaries_type == boundaries_unsupported) {
109     cvm::error("Error: unsupported boundary conditions.\n", INPUT_ERROR);
110   }
111
112   cvm::rvector diff = (pos2 - pos1);
113
114   if (boundaries_type == boundaries_non_periodic) return diff;
115
116   cvm::real const x_shift = round_to_integer(reciprocal_cell_x*diff);
117   cvm::real const y_shift = round_to_integer(reciprocal_cell_y*diff);
118   cvm::real const z_shift = round_to_integer(reciprocal_cell_z*diff);
119
120   diff.x -= x_shift*unit_cell_x.x + y_shift*unit_cell_y.x +
121     z_shift*unit_cell_z.x;
122   diff.y -= x_shift*unit_cell_x.y + y_shift*unit_cell_y.y +
123     z_shift*unit_cell_z.y;
124   diff.z -= x_shift*unit_cell_x.z + y_shift*unit_cell_y.z +
125     z_shift*unit_cell_z.z;
126
127   return diff;
128 }
129
130
131
132 colvarproxy_atoms::colvarproxy_atoms() {}
133
134
135 colvarproxy_atoms::~colvarproxy_atoms()
136 {
137   reset();
138 }
139
140
141 int colvarproxy_atoms::reset()
142 {
143   atoms_ids.clear();
144   atoms_ncopies.clear();
145   atoms_masses.clear();
146   atoms_charges.clear();
147   atoms_positions.clear();
148   atoms_total_forces.clear();
149   atoms_new_colvar_forces.clear();
150   return COLVARS_OK;
151 }
152
153
154 int colvarproxy_atoms::add_atom_slot(int atom_id)
155 {
156   atoms_ids.push_back(atom_id);
157   atoms_ncopies.push_back(1);
158   atoms_masses.push_back(1.0);
159   atoms_charges.push_back(0.0);
160   atoms_positions.push_back(cvm::rvector(0.0, 0.0, 0.0));
161   atoms_total_forces.push_back(cvm::rvector(0.0, 0.0, 0.0));
162   atoms_new_colvar_forces.push_back(cvm::rvector(0.0, 0.0, 0.0));
163   return (atoms_ids.size() - 1);
164 }
165
166
167 int colvarproxy_atoms::init_atom(cvm::residue_id const &residue,
168                                  std::string const     &atom_name,
169                                  std::string const     &segment_id)
170 {
171   cvm::error("Error: initializing an atom by name and residue number is currently not supported.\n",
172              COLVARS_NOT_IMPLEMENTED);
173   return COLVARS_NOT_IMPLEMENTED;
174 }
175
176
177 int colvarproxy_atoms::check_atom_id(cvm::residue_id const &residue,
178                                      std::string const     &atom_name,
179                                      std::string const     &segment_id)
180 {
181   colvarproxy_atoms::init_atom(residue, atom_name, segment_id);
182   return COLVARS_NOT_IMPLEMENTED;
183 }
184
185
186 void colvarproxy_atoms::clear_atom(int index)
187 {
188   if (((size_t) index) >= atoms_ids.size()) {
189     cvm::error("Error: trying to disable an atom that was not previously requested.\n",
190                INPUT_ERROR);
191   }
192   if (atoms_ncopies[index] > 0) {
193     atoms_ncopies[index] -= 1;
194   }
195 }
196
197
198 int colvarproxy_atoms::load_atoms(char const *filename,
199                                   cvm::atom_group &atoms,
200                                   std::string const &pdb_field,
201                                   double)
202 {
203   return cvm::error("Error: loading atom identifiers from a file "
204                     "is currently not implemented.\n",
205                     COLVARS_NOT_IMPLEMENTED);
206 }
207
208
209 int colvarproxy_atoms::load_coords(char const *filename,
210                                    std::vector<cvm::atom_pos> &pos,
211                                    std::vector<int> const &sorted_ids,
212                                    std::string const &pdb_field,
213                                    double)
214 {
215   return cvm::error("Error: loading atomic coordinates from a file "
216                     "is currently not implemented.\n",
217                     COLVARS_NOT_IMPLEMENTED);
218 }
219
220
221
222 colvarproxy_atom_groups::colvarproxy_atom_groups() {}
223
224
225 colvarproxy_atom_groups::~colvarproxy_atom_groups()
226 {
227   reset();
228 }
229
230
231 int colvarproxy_atom_groups::reset()
232 {
233   atom_groups_ids.clear();
234   atom_groups_ncopies.clear();
235   atom_groups_masses.clear();
236   atom_groups_charges.clear();
237   atom_groups_coms.clear();
238   atom_groups_total_forces.clear();
239   atom_groups_new_colvar_forces.clear();
240   return COLVARS_OK;
241 }
242
243
244 int colvarproxy_atom_groups::add_atom_group_slot(int atom_group_id)
245 {
246   atom_groups_ids.push_back(atom_group_id);
247   atom_groups_ncopies.push_back(1);
248   atom_groups_masses.push_back(1.0);
249   atom_groups_charges.push_back(0.0);
250   atom_groups_coms.push_back(cvm::rvector(0.0, 0.0, 0.0));
251   atom_groups_total_forces.push_back(cvm::rvector(0.0, 0.0, 0.0));
252   atom_groups_new_colvar_forces.push_back(cvm::rvector(0.0, 0.0, 0.0));
253   return (atom_groups_ids.size() - 1);
254 }
255
256
257 int colvarproxy_atom_groups::scalable_group_coms()
258 {
259   return COLVARS_NOT_IMPLEMENTED;
260 }
261
262
263 int colvarproxy_atom_groups::init_atom_group(std::vector<int> const &atoms_ids)
264 {
265   cvm::error("Error: initializing a group outside of the Colvars module "
266              "is currently not supported.\n",
267              COLVARS_NOT_IMPLEMENTED);
268   return COLVARS_NOT_IMPLEMENTED;
269 }
270
271
272 void colvarproxy_atom_groups::clear_atom_group(int index)
273 {
274   if (((size_t) index) >= atom_groups_ids.size()) {
275     cvm::error("Error: trying to disable an atom group "
276                "that was not previously requested.\n",
277                INPUT_ERROR);
278   }
279   if (atom_groups_ncopies[index] > 0) {
280     atom_groups_ncopies[index] -= 1;
281   }
282 }
283
284
285
286 colvarproxy_smp::colvarproxy_smp()
287 {
288   b_smp_active = true; // May be disabled by user option
289   omp_lock_state = NULL;
290 #if defined(_OPENMP)
291   if (smp_thread_id() == 0) {
292     omp_init_lock(reinterpret_cast<omp_lock_t *>(omp_lock_state));
293   }
294 #endif
295 }
296
297
298 colvarproxy_smp::~colvarproxy_smp() {}
299
300
301 int colvarproxy_smp::smp_enabled()
302 {
303 #if defined(_OPENMP)
304   if (b_smp_active) {
305     return COLVARS_OK;
306   }
307   return COLVARS_ERROR;
308 #else
309   return COLVARS_NOT_IMPLEMENTED;
310 #endif
311 }
312
313
314 int colvarproxy_smp::smp_colvars_loop()
315 {
316 #if defined(_OPENMP)
317   colvarmodule *cv = cvm::main();
318   colvarproxy *proxy = cv->proxy;
319 #pragma omp parallel for
320   for (size_t i = 0; i < cv->variables_active_smp()->size(); i++) {
321     colvar *x = (*(cv->variables_active_smp()))[i];
322     int x_item = (*(cv->variables_active_smp_items()))[i];
323     if (cvm::debug()) {
324       cvm::log("["+cvm::to_str(proxy->smp_thread_id())+"/"+
325                cvm::to_str(proxy->smp_num_threads())+
326                "]: calc_colvars_items_smp(), i = "+cvm::to_str(i)+", cv = "+
327                x->name+", cvc = "+cvm::to_str(x_item)+"\n");
328     }
329     x->calc_cvcs(x_item, 1);
330   }
331   return cvm::get_error();
332 #else
333   return COLVARS_NOT_IMPLEMENTED;
334 #endif
335 }
336
337
338 int colvarproxy_smp::smp_biases_loop()
339 {
340 #if defined(_OPENMP)
341   colvarmodule *cv = cvm::main();
342 #pragma omp parallel
343   {
344 #pragma omp for
345     for (size_t i = 0; i < cv->biases_active()->size(); i++) {
346       colvarbias *b = (*(cv->biases_active()))[i];
347       if (cvm::debug()) {
348         cvm::log("Calculating bias \""+b->name+"\" on thread "+
349                  cvm::to_str(smp_thread_id())+"\n");
350       }
351       b->update();
352     }
353   }
354   return cvm::get_error();
355 #else
356   return COLVARS_NOT_IMPLEMENTED;
357 #endif
358 }
359
360
361 int colvarproxy_smp::smp_biases_script_loop()
362 {
363 #if defined(_OPENMP)
364   colvarmodule *cv = cvm::main();
365 #pragma omp parallel
366   {
367 #pragma omp single nowait
368     {
369       cv->calc_scripted_forces();
370     }
371 #pragma omp for
372     for (size_t i = 0; i < cv->biases_active()->size(); i++) {
373       colvarbias *b = (*(cv->biases_active()))[i];
374       if (cvm::debug()) {
375         cvm::log("Calculating bias \""+b->name+"\" on thread "+
376                  cvm::to_str(smp_thread_id())+"\n");
377       }
378       b->update();
379     }
380   }
381   return cvm::get_error();
382 #else
383   return COLVARS_NOT_IMPLEMENTED;
384 #endif
385 }
386
387
388
389
390 int colvarproxy_smp::smp_thread_id()
391 {
392 #if defined(_OPENMP)
393   return omp_get_thread_num();
394 #else
395   return COLVARS_NOT_IMPLEMENTED;
396 #endif
397 }
398
399
400 int colvarproxy_smp::smp_num_threads()
401 {
402 #if defined(_OPENMP)
403   return omp_get_max_threads();
404 #else
405   return COLVARS_NOT_IMPLEMENTED;
406 #endif
407 }
408
409
410 int colvarproxy_smp::smp_lock()
411 {
412 #if defined(_OPENMP)
413   omp_set_lock(reinterpret_cast<omp_lock_t *>(omp_lock_state));
414 #endif
415   return COLVARS_OK;
416 }
417
418
419 int colvarproxy_smp::smp_trylock()
420 {
421 #if defined(_OPENMP)
422   return omp_test_lock(reinterpret_cast<omp_lock_t *>(omp_lock_state)) ?
423     COLVARS_OK : COLVARS_ERROR;
424 #else
425   return COLVARS_OK;
426 #endif
427 }
428
429
430 int colvarproxy_smp::smp_unlock()
431 {
432 #if defined(_OPENMP)
433   omp_unset_lock(reinterpret_cast<omp_lock_t *>(omp_lock_state));
434 #endif
435   return COLVARS_OK;
436 }
437
438
439
440 colvarproxy_replicas::colvarproxy_replicas() {}
441
442
443 colvarproxy_replicas::~colvarproxy_replicas() {}
444
445
446 bool colvarproxy_replicas::replica_enabled()
447 {
448   return false;
449 }
450
451
452 int colvarproxy_replicas::replica_index()
453 {
454   return 0;
455 }
456
457
458 int colvarproxy_replicas::replica_num()
459 {
460   return 1;
461 }
462
463
464 void colvarproxy_replicas::replica_comm_barrier() {}
465
466
467 int colvarproxy_replicas::replica_comm_recv(char* msg_data,
468                                             int buf_len,
469                                             int src_rep)
470 {
471   return COLVARS_NOT_IMPLEMENTED;
472 }
473
474
475 int colvarproxy_replicas::replica_comm_send(char* msg_data,
476                                             int msg_len,
477                                             int dest_rep)
478 {
479   return COLVARS_NOT_IMPLEMENTED;
480 }
481
482
483
484
485 colvarproxy_script::colvarproxy_script()
486 {
487   script = NULL;
488 }
489
490
491 colvarproxy_script::~colvarproxy_script() {}
492
493
494 char const *colvarproxy_script::script_obj_to_str(unsigned char *obj)
495 {
496   cvm::error("Error: trying to print a script object without a scripting "
497              "language interface.\n", BUG_ERROR);
498   return reinterpret_cast<char *>(obj);
499 }
500
501
502 int colvarproxy_script::run_force_callback()
503 {
504   return COLVARS_NOT_IMPLEMENTED;
505 }
506
507
508 int colvarproxy_script::run_colvar_callback(
509       std::string const &name,
510       std::vector<const colvarvalue *> const &cvcs,
511       colvarvalue &value)
512 {
513   return COLVARS_NOT_IMPLEMENTED;
514 }
515
516
517 int colvarproxy_script::run_colvar_gradient_callback(
518       std::string const &name,
519       std::vector<const colvarvalue *> const &cvcs,
520       std::vector<cvm::matrix2d<cvm::real> > &gradient)
521 {
522   return COLVARS_NOT_IMPLEMENTED;
523 }
524
525
526
527 colvarproxy_tcl::colvarproxy_tcl()
528 {
529   _tcl_interp = NULL;
530 }
531
532
533 colvarproxy_tcl::~colvarproxy_tcl()
534 {
535 }
536
537
538 void colvarproxy_tcl::init_tcl_pointers()
539 {
540   cvm::error("Error: Tcl support is currently unavailable "
541              "outside NAMD or VMD.\n", COLVARS_NOT_IMPLEMENTED);
542 }
543
544
545 char const *colvarproxy_tcl::tcl_obj_to_str(unsigned char *obj)
546 {
547 #if defined(COLVARS_TCL)
548   return Tcl_GetString(reinterpret_cast<Tcl_Obj *>(obj));
549 #else
550   return NULL;
551 #endif
552 }
553
554
555 int colvarproxy_tcl::tcl_run_force_callback()
556 {
557 #if defined(COLVARS_TCL)
558   Tcl_Interp *const tcl_interp = reinterpret_cast<Tcl_Interp *>(_tcl_interp);
559   std::string cmd = std::string("calc_colvar_forces ")
560     + cvm::to_str(cvm::step_absolute());
561   int err = Tcl_Eval(tcl_interp, cmd.c_str());
562   if (err != TCL_OK) {
563     cvm::log(std::string("Error while executing calc_colvar_forces:\n"));
564     cvm::error(Tcl_GetStringResult(tcl_interp));
565     return COLVARS_ERROR;
566   }
567   return cvm::get_error();
568 #else
569   return COLVARS_NOT_IMPLEMENTED;
570 #endif
571 }
572
573
574 int colvarproxy_tcl::tcl_run_colvar_callback(
575                          std::string const &name,
576                          std::vector<const colvarvalue *> const &cvc_values,
577                          colvarvalue &value)
578 {
579 #if defined(COLVARS_TCL)
580
581   Tcl_Interp *const tcl_interp = reinterpret_cast<Tcl_Interp *>(_tcl_interp);
582   size_t i;
583   std::string cmd = std::string("calc_") + name;
584   for (i = 0; i < cvc_values.size(); i++) {
585     cmd += std::string(" {") + (*(cvc_values[i])).to_simple_string() +
586       std::string("}");
587   }
588   int err = Tcl_Eval(tcl_interp, cmd.c_str());
589   const char *result = Tcl_GetStringResult(tcl_interp);
590   if (err != TCL_OK) {
591     return cvm::error(std::string("Error while executing ")
592                       + cmd + std::string(":\n") +
593                       std::string(Tcl_GetStringResult(tcl_interp)), COLVARS_ERROR);
594   }
595   std::istringstream is(result);
596   if (value.from_simple_string(is.str()) != COLVARS_OK) {
597     cvm::log("Error parsing colvar value from script:");
598     cvm::error(result);
599     return COLVARS_ERROR;
600   }
601   return cvm::get_error();
602
603 #else
604
605   return COLVARS_NOT_IMPLEMENTED;
606
607 #endif
608 }
609
610
611 int colvarproxy_tcl::tcl_run_colvar_gradient_callback(
612                          std::string const &name,
613                          std::vector<const colvarvalue *> const &cvc_values,
614                          std::vector<cvm::matrix2d<cvm::real> > &gradient)
615 {
616 #if defined(COLVARS_TCL)
617
618   Tcl_Interp *const tcl_interp = reinterpret_cast<Tcl_Interp *>(_tcl_interp);
619   size_t i;
620   std::string cmd = std::string("calc_") + name + "_gradient";
621   for (i = 0; i < cvc_values.size(); i++) {
622     cmd += std::string(" {") + (*(cvc_values[i])).to_simple_string() +
623       std::string("}");
624   }
625   int err = Tcl_Eval(tcl_interp, cmd.c_str());
626   if (err != TCL_OK) {
627     return cvm::error(std::string("Error while executing ")
628                       + cmd + std::string(":\n") +
629                       std::string(Tcl_GetStringResult(tcl_interp)), COLVARS_ERROR);
630   }
631   Tcl_Obj **list;
632   int n;
633   Tcl_ListObjGetElements(tcl_interp, Tcl_GetObjResult(tcl_interp),
634                          &n, &list);
635   if (n != int(gradient.size())) {
636     cvm::error("Error parsing list of gradient values from script: found "
637                + cvm::to_str(n) + " values instead of " +
638                cvm::to_str(gradient.size()));
639     return COLVARS_ERROR;
640   }
641   for (i = 0; i < gradient.size(); i++) {
642     std::istringstream is(Tcl_GetString(list[i]));
643     if (gradient[i].from_simple_string(is.str()) != COLVARS_OK) {
644       cvm::log("Gradient matrix size: " + cvm::to_str(gradient[i].size()));
645       cvm::log("Gradient string: " + cvm::to_str(Tcl_GetString(list[i])));
646       cvm::error("Error parsing gradient value from script", COLVARS_ERROR);
647       return COLVARS_ERROR;
648     }
649   }
650
651   return cvm::get_error();
652
653 #else
654
655   return COLVARS_NOT_IMPLEMENTED;
656
657 #endif
658 }
659
660
661
662 colvarproxy_io::colvarproxy_io() {}
663
664
665 colvarproxy_io::~colvarproxy_io() {}
666
667
668 int colvarproxy_io::get_frame(long int&)
669 {
670   return COLVARS_NOT_IMPLEMENTED;
671 }
672
673
674 int colvarproxy_io::set_frame(long int)
675 {
676   return COLVARS_NOT_IMPLEMENTED;
677 }
678
679
680 std::ostream * colvarproxy_io::output_stream(std::string const &output_name,
681                                              std::ios_base::openmode mode)
682 {
683   if (cvm::debug()) {
684     cvm::log("Using colvarproxy::output_stream()\n");
685   }
686   std::list<std::ostream *>::iterator osi  = output_files.begin();
687   std::list<std::string>::iterator    osni = output_stream_names.begin();
688   for ( ; osi != output_files.end(); osi++, osni++) {
689     if (*osni == output_name) {
690       return *osi;
691     }
692   }
693   if (!(mode & (std::ios_base::app | std::ios_base::ate))) {
694     backup_file(output_name);
695   }
696   std::ofstream *os = new std::ofstream(output_name.c_str(), mode);
697   if (!os->is_open()) {
698     cvm::error("Error: cannot write to file/channel \""+output_name+"\".\n",
699                FILE_ERROR);
700     return NULL;
701   }
702   output_stream_names.push_back(output_name);
703   output_files.push_back(os);
704   return os;
705 }
706
707
708 int colvarproxy_io::flush_output_stream(std::ostream *os)
709 {
710   std::list<std::ostream *>::iterator osi  = output_files.begin();
711   std::list<std::string>::iterator    osni = output_stream_names.begin();
712   for ( ; osi != output_files.end(); osi++, osni++) {
713     if (*osi == os) {
714       ((std::ofstream *) (*osi))->flush();
715       return COLVARS_OK;
716     }
717   }
718   return cvm::error("Error: trying to flush an output file/channel "
719                     "that wasn't open.\n", BUG_ERROR);
720 }
721
722
723 int colvarproxy_io::close_output_stream(std::string const &output_name)
724 {
725   std::list<std::ostream *>::iterator osi  = output_files.begin();
726   std::list<std::string>::iterator    osni = output_stream_names.begin();
727   for ( ; osi != output_files.end(); osi++, osni++) {
728     if (*osni == output_name) {
729       ((std::ofstream *) (*osi))->close();
730       delete *osi;
731       output_files.erase(osi);
732       output_stream_names.erase(osni);
733       return COLVARS_OK;
734     }
735   }
736   return cvm::error("Error: trying to close an output file/channel "
737                     "that wasn't open.\n", BUG_ERROR);
738 }
739
740
741 int colvarproxy_io::backup_file(char const *filename)
742 {
743   return COLVARS_NOT_IMPLEMENTED;
744 }
745
746
747
748 colvarproxy::colvarproxy()
749 {
750   colvars = NULL;
751   b_simulation_running = true;
752   b_delete_requested = false;
753 }
754
755
756 colvarproxy::~colvarproxy() {}
757
758
759 int colvarproxy::reset()
760 {
761   int error_code = COLVARS_OK;
762   error_code |= colvarproxy_atoms::reset();
763   error_code |= colvarproxy_atom_groups::reset();
764   return error_code;
765 }
766
767
768 int colvarproxy::request_deletion()
769 {
770   return cvm::error("Error: \"delete\" command is only available in VMD; "
771                     "please use \"reset\" instead.\n",
772                     COLVARS_NOT_IMPLEMENTED);
773 }
774
775
776 int colvarproxy::setup()
777 {
778   return COLVARS_OK;
779 }
780
781
782 int colvarproxy::update_input()
783 {
784   return COLVARS_OK;
785 }
786
787
788 int colvarproxy::update_output()
789 {
790   return COLVARS_OK;
791 }
792
793
794 size_t colvarproxy::restart_frequency()
795 {
796   return 0;
797 }
798
799
800 int colvarproxy::get_version_from_string(char const *version_string)
801 {
802   std::string const v(version_string);
803   std::istringstream is(v.substr(0, 4) + v.substr(5, 2) + v.substr(8, 2));
804   int newint;
805   is >> newint;
806   return newint;
807 }
808