NAMD
colvarproxy_namd.C
Go to the documentation of this file.
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 <errno.h>
11 
12 #include "common.h"
13 #include "fstream_namd.h"
14 #include "Debug.h"
15 #include "BackEnd.h"
16 #include "InfoStream.h"
17 #include "Node.h"
18 #include "Molecule.h"
19 #include "GridForceGrid.h"
20 #include "GridForceGrid.inl"
21 #include "PDB.h"
22 #include "PDBData.h"
23 #include "ReductionMgr.h"
24 #include "ScriptTcl.h"
25 #include "NamdState.h"
26 #include "Controller.h"
27 #include "PatchData.h"
28 
29 #ifdef NAMD_TCL
30 #include <tcl.h>
31 #endif
32 
33 // For replica exchange
34 #include "converse.h"
35 #include "DataExchanger.h"
36 
37 #include "colvarmodule.h"
38 #include "colvar.h"
39 #include "colvarbias.h"
40 #include "colvaratoms.h"
41 #include "colvarproxy.h"
42 #include "colvarproxy_namd.h"
44 #include "colvarscript.h"
45 
46 
48 {
49  engine_name_ = "NAMD";
50 #if CMK_SMP && USE_CKLOOP
51  charm_lock_state = CmiCreateLock();
52 #endif
53 
54  version_int = get_version_from_string(COLVARPROXY_VERSION);
55 #if CMK_TRACE_ENABLED
56  if ( 0 == CkMyPe() ) {
57  traceRegisterUserEvent("GM COLVAR item", GLOBAL_MASTER_CKLOOP_CALC_ITEM);
58  traceRegisterUserEvent("GM COLVAR bias", GLOBAL_MASTER_CKLOOP_CALC_BIASES );
59  traceRegisterUserEvent("GM COLVAR scripted bias", GLOBAL_MASTER_CKLOOP_CALC_SCRIPTED_BIASES );
60  }
61 #endif
62  first_timestep = true;
63  requestTotalForce(total_force_requested);
64 
65  boltzmann_ = 0.001987191;
66 
67  angstrom_value_ = 1.;
68 
69  // initialize pointers to NAMD configuration data
71 
72  if (cvm::debug())
73  iout << "Info: initializing the colvars proxy object.\n" << endi;
74 
75  // find the configuration file, if provided
76  StringList *config = Node::Object()->configList->find("colvarsConfig");
77 
78  // find the input state file
79  StringList *input_restart = Node::Object()->configList->find("colvarsInput");
80  colvarproxy_io::set_input_prefix(input_restart ? input_restart->data : "");
81 
83  set_integration_timestep(simparams->dt);
84  set_time_step_factor(simparams->globalMasterFrequency);
85 
87 
88  // both fields are taken from data structures already available
89  updated_masses_ = updated_charges_ = true;
90 
91  // Take the output prefixes from the NAMD input
92  colvarproxy_io::set_output_prefix(std::string(simparams->outputFilename));
93  colvarproxy_io::set_restart_output_prefix(std::string(simparams->restartFilename));
94  colvarproxy_io::set_default_restart_frequency(simparams->restartFrequency);
95 
96  if (simparams->accelMDOn) {
97  accelMDOn = true;
98  } else {
99  accelMDOn = false;
100  }
101  amd_weight_factor = 1.0;
102 
103  // check if it is possible to save output configuration
104  if ((!output_prefix_str.size()) && (!restart_output_prefix_str.size())) {
105  error("Error: neither the final output state file or "
106  "the output restart file could be defined, exiting.\n");
107  }
108 
109  init_atoms_map();
110 
111  // initialize module: this object will be the communication proxy
112  colvars = new colvarmodule(this);
113 
114  cvm::log("Using NAMD interface, version "+
115  cvm::to_str(COLVARPROXY_VERSION)+".\n");
116  colvars->cite_feature("NAMD engine");
117  colvars->cite_feature("Colvars-NAMD interface");
118 
119  errno = 0;
120  for ( ; config; config = config->next ) {
121  add_config("configfile", config->data);
122  }
123 
124  // Trigger immediate initialization of the module
125  colvarproxy::parse_module_config();
127  colvars->update_engine_parameters();
128  colvars->setup_input();
129  colvars->setup_output();
130 
131  // save to Node for Tcl script access
132  Node::Object()->colvars = colvars;
133 
134 #ifdef NAMD_TCL
135  have_scripts = true;
136 #else
137  have_scripts = false;
138 #endif
139 
140  if (simparams->firstTimestep != 0) {
141  colvars->set_initial_step(static_cast<cvm::step_number>(simparams->firstTimestep));
142  }
143 
144 #if !defined (NAMD_UNIFIED_REDUCTION)
146 #endif
147 
148  #if defined(NODEGROUP_FORCE_REGISTER) && !defined(NAMD_UNIFIED_REDUCTION)
149  CProxy_PatchData cpdata(CkpvAccess(BOCclass_group).patchData);
150  PatchData *patchData = cpdata.ckLocalBranch();
151  nodeReduction = patchData->reduction;
152  #endif
153 
154  if (cvm::debug())
155  iout << "Info: done initializing the colvars proxy object.\n" << endi;
156 }
157 
158 
160 {
161 #if CMK_SMP && USE_CKLOOP
162  CmiDestroyLock(charm_lock_state);
163 #endif
164 #if !defined (NAMD_UNIFIED_REDUCTION)
165  delete reduction;
166 #endif
167 }
168 
169 
171 {
172  int error_code = COLVARS_OK;
173  if (simparams->rescaleFreq > 0) {
174  error_code |= set_target_temperature(simparams->rescaleTemp);
175  } else if (simparams->reassignFreq > 0) {
176  error_code |= set_target_temperature(simparams->reassignTemp);
177  } else if (simparams->langevinOn) {
178  error_code |= set_target_temperature(simparams->langevinTemp);
179  } else if (simparams->tCoupleOn) {
180  error_code |= set_target_temperature(simparams->tCoupleTemp);
181  } else if (simparams->loweAndersenOn) {
182  error_code |= set_target_temperature(simparams->loweAndersenTemp);
183  } else if (simparams->stochRescaleOn) {
184  error_code |= set_target_temperature(simparams->stochRescaleTemp);
185  } else {
186  error_code |= set_target_temperature(0.0);
187  }
188  return error_code;
189 }
190 
191 
192 
194 {
195  size_t const n_all_atoms = Node::Object()->molecule->numAtoms;
196  atoms_map.assign(n_all_atoms, -1);
197 }
198 
199 
202 {
203  if (atoms_map.size() != Node::Object()->molecule->numAtoms) {
204  init_atoms_map();
205  }
206 
207  for (AtomIDList::const_iterator a_i = begin; a_i != end; a_i++) {
208 
209  if (atoms_map[*a_i] >= 0) continue;
210 
211  for (size_t i = 0; i < atoms_ids.size(); i++) {
212  if (atoms_ids[i] == *a_i) {
213  atoms_map[*a_i] = i;
214  break;
215  }
216  }
217 
218  if (atoms_map[*a_i] < 0) {
219  // this atom is probably managed by another GlobalMaster:
220  // add it here anyway to avoid having to test for array boundaries at each step
221  int const index = add_atom_slot(*a_i);
222  atoms_map[*a_i] = index;
223  modifyRequestedAtoms().add(*a_i);
224  update_atom_properties(index);
225  }
226  }
227 
228  if (cvm::debug()) {
229  cvm::log("atoms_map = "+cvm::to_str(atoms_map)+".\n");
230  }
231 
232  return COLVARS_OK;
233 }
234 
235 
237 {
238  int error_code = colvarproxy::setup();
239 
240  if (colvars->size() == 0) {
241  // Module is empty, nothing to do
242  return COLVARS_OK;
243  }
244 
245  log("Updating NAMD interface:\n");
246 
247  errno = 0;
248 
249  if (simparams->wrapAll) {
250  log("Warning: enabling wrapAll can lead to inconsistent results "
251  "for Colvars calculations: please disable wrapAll, "
252  "as is the default option in NAMD.\n");
253  }
254 
255  log("updating atomic data ("+cvm::to_str(atoms_ids.size())+" atoms).\n");
256 
257  size_t i;
258  for (i = 0; i < atoms_ids.size(); i++) {
260 
261  // zero out mutable arrays
262  atoms_positions[i] = cvm::rvector(0.0, 0.0, 0.0);
263  atoms_total_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
264  atoms_new_colvar_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
265  }
266 
267  size_t n_group_atoms = 0;
268  for (int ig = 0; ig < modifyRequestedGroups().size(); ig++) {
269  n_group_atoms += modifyRequestedGroups()[ig].size();
270  }
271 
272  log("updating group data ("+cvm::to_str(atom_groups_ids.size())+
273  " scalable groups, "+
274  cvm::to_str(n_group_atoms)+" atoms in total).\n");
275 
276  // Note: groupMassBegin, groupMassEnd may be used here, but they won't work for charges
277  for (int ig = 0; ig < modifyRequestedGroups().size(); ig++) {
278 
279  // update mass and charge
281 
282  atom_groups_coms[ig] = cvm::rvector(0.0, 0.0, 0.0);
283  atom_groups_total_forces[ig] = cvm::rvector(0.0, 0.0, 0.0);
284  atom_groups_new_colvar_forces[ig] = cvm::rvector(0.0, 0.0, 0.0);
285  }
286 
287 #if NAMD_VERSION_NUMBER >= 34471681
288  log("updating grid object data ("+cvm::to_str(volmaps_ids.size())+
289  " grid objects in total).\n");
290  for (int imap = 0; imap < modifyGridObjForces().size(); imap++) {
291  volmaps_new_colvar_forces[imap] = 0.0;
292  }
293 #endif
294 
295  size_t const new_features_hash =
296  std::hash<std::string>{}(colvars->feature_report(0));
297  if (new_features_hash != features_hash) {
298  // Nag only once, there may be many run commands
299  log(std::string("\n")+colvars->feature_report(0)+std::string("\n"));
300  features_hash = new_features_hash;
301  }
302 
304  log("updating target temperature (T = "+
305  cvm::to_str(target_temperature())+" K).\n");
306 
307  // Note: not needed currently, but may be in the future if NAMD allows
308  // redefining the timestep
309  set_integration_timestep(simparams->dt);
310 
311  return error_code;
312 }
313 
314 
316 {
317  if (cvm::debug()) {
318  cvm::log("colvarproxy_namd::reset()\n");
319  }
320 
321  int error_code = COLVARS_OK;
322 
323  // Unrequest all positions, total forces, etc from NAMD
327 
328  modifyRequestedGroups().clear();
330 
331 #if NAMD_VERSION_NUMBER >= 34471681
334 #endif
335 
336  requestTotalForce(false);
337 
338  atoms_map.clear();
339 
340  // Clear internal atomic data
341  error_code |= colvarproxy::reset();
342 
343  return error_code;
344 }
345 
346 
348 {
349  errno = 0;
350 
351  if (first_timestep) {
352 
353  // First run after the proxy is constructed
354 
356  colvars->update_engine_parameters();
357  colvars->setup_input();
358  colvars->setup_output();
359  // Controller is only available after full startup phase, so now
361 
362  first_timestep = false;
363 
364  } else {
365 
366  // Use the time step number inherited from GlobalMaster
367  if ( step - previous_NAMD_step == time_step_factor() ) {
368  colvars->it += time_step_factor();
369  b_simulation_continuing = false;
370  } else {
371 
372  // Cases covered by this condition:
373  // - run 0
374  // - beginning of a new run statement
375  // The internal counter is not incremented, and the objects are made
376  // aware of this via the following flag
377  b_simulation_continuing = true;
378 
379  // Update NAMD output and restart prefixes
380  colvarproxy_io::set_output_prefix(std::string(simparams->outputFilename));
381  colvarproxy_io::set_restart_output_prefix(std::string(simparams->restartFilename));
382  colvarproxy_io::set_default_restart_frequency(simparams->restartFrequency);
383  colvars->setup_output();
384 
385  }
386  }
387 
390 
391  {
392  Vector const a = lattice->a();
393  Vector const b = lattice->b();
394  Vector const c = lattice->c();
395  unit_cell_x.set(a.x, a.y, a.z);
396  unit_cell_y.set(b.x, b.y, c.z);
397  unit_cell_z.set(c.x, c.y, c.z);
398  }
399 
400  if (!lattice->a_p() && !lattice->b_p() && !lattice->c_p()) {
401  boundaries_type = boundaries_non_periodic;
402  reset_pbc_lattice();
403  } else if (lattice->a_p() && lattice->b_p() && lattice->c_p()) {
404  if (lattice->orthogonal()) {
405  boundaries_type = boundaries_pbc_ortho;
406  } else {
407  boundaries_type = boundaries_pbc_triclinic;
408  }
409  colvarproxy_system::update_pbc_lattice();
410  } else {
411  boundaries_type = boundaries_unsupported;
412  }
413 
414  if (cvm::debug()) {
415  cvm::log(std::string(cvm::line_marker)+
416  "colvarproxy_namd, step no. "+cvm::to_str(colvars->it)+"\n"+
417  "Updating atomic data arrays.\n");
418  }
419 
420  // must delete the forces applied at the previous step: we can do
421  // that because they have already been used and copied to other
422  // memory locations
425 
426  // If new atomic positions or forces have been requested by other
427  // GlobalMaster objects, add these to the atom map as well
428  size_t const n_all_atoms = Node::Object()->molecule->numAtoms;
429  if ( (atoms_map.size() != n_all_atoms) ||
430  (int(atoms_ids.size()) < (getAtomIdEnd() - getAtomIdBegin())) ||
431  (int(atoms_ids.size()) < (getForceIdEnd() - getForceIdBegin())) ) {
434  }
435 
436  // prepare local arrays
437  for (size_t i = 0; i < atoms_ids.size(); i++) {
438  atoms_positions[i] = cvm::rvector(0.0, 0.0, 0.0);
439  atoms_total_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
440  atoms_new_colvar_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
441  }
442 
443  for (size_t i = 0; i < atom_groups_ids.size(); i++) {
444  atom_groups_total_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
445  atom_groups_new_colvar_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
446  }
447 
448 #if NAMD_VERSION_NUMBER >= 34471681
449  for (int imap = 0; imap < volmaps_ids.size(); imap++) {
450  volmaps_new_colvar_forces[imap] = 0.0;
451  }
452 #endif
453 
454  {
455  if (cvm::debug()) {
456  cvm::log("Updating positions arrays.\n");
457  }
458  size_t n_positions = 0;
462 
463  for ( ; a_i != a_e; ++a_i, ++p_i ) {
464  atoms_positions[atoms_map[*a_i]] = cvm::rvector((*p_i).x, (*p_i).y, (*p_i).z);
465  n_positions++;
466  }
467 
468  // The following had to be relaxed because some atoms may be forced without their position being requested
469  // if (n_positions < atoms_ids.size()) {
470  // cvm::error("Error: did not receive the positions of all atoms.\n", COLVARS_BUG_ERROR);
471  // }
472  }
473 
474  if (total_force_requested && cvm::step_relative() > 0) {
475 
476  // sort the force arrays the previous step
477  // (but only do so if there *is* a previous step!)
478 
479  {
480  if (cvm::debug()) {
481  cvm::log("Updating total forces arrays.\n");
482  }
483  size_t n_total_forces = 0;
487 
488  for ( ; a_i != a_e; ++a_i, ++f_i ) {
489  atoms_total_forces[atoms_map[*a_i]] = cvm::rvector((*f_i).x, (*f_i).y, (*f_i).z);
490  n_total_forces++;
491  }
492 
493  if ( (! b_simulation_continuing) &&
494  (n_total_forces < atoms_ids.size()) ) {
495  cvm::error("Error: total forces were requested, but total forces "
496  "were not received for all atoms.\n"
497  "The most probable cause is combination of energy "
498  "minimization with a biasing method that requires MD (e.g. ABF).\n"
499  "Always run minimization and ABF separately.", COLVARS_INPUT_ERROR);
500  }
501  }
502 
503  {
504  if (cvm::debug()) {
505  cvm::log("Updating group total forces arrays.\n");
506  }
509  size_t i = 0;
510  if ( (! b_simulation_continuing) &&
511  ((f_e - f_i) != ((int) atom_groups_ids.size())) ) {
512  cvm::error("Error: total forces were requested for scalable groups, "
513  "but they are not in the same number from the number of groups.\n"
514  "The most probable cause is combination of energy "
515  "minimization with a biasing method that requires MD (e.g. ABF).\n"
516  "Always run minimization and ABF separately.", COLVARS_INPUT_ERROR);
517  }
518  for ( ; f_i != f_e; f_i++, i++) {
519  atom_groups_total_forces[i] = cvm::rvector((*f_i).x, (*f_i).y, (*f_i).z);
520  }
521  }
522  }
523 
524  {
525  if (cvm::debug()) {
526  cvm::log("Updating group positions arrays.\n");
527  }
528  // update group data (only coms available so far)
529  size_t ig;
530  // note: getGroupMassBegin() could be used here, but masses and charges
531  // have already been calculated from the last call to setup()
533  for (ig = 0; gp_i != getGroupPositionEnd(); gp_i++, ig++) {
534  atom_groups_coms[ig] = cvm::rvector(gp_i->x, gp_i->y, gp_i->z);
535  }
536  }
537 
538 #if NAMD_VERSION_NUMBER >= 34471681
539  {
540  if (cvm::debug()) {
541  log("Updating grid objects.\n");
542  }
543  // Using a simple nested loop: there probably won't be so many maps that
544  // this becomes performance-limiting
547  for ( ; gov_i != getGridObjValueEnd(); goi_i++, gov_i++) {
548  for (size_t imap = 0; imap < volmaps_ids.size(); imap++) {
549  if (volmaps_ids[imap] == *goi_i) {
550  volmaps_values[imap] = *gov_i;
551  break;
552  }
553  }
554  }
555  }
556 #endif
557 
558  if (cvm::debug()) {
559  print_input_atomic_data();
560  }
561 
562  // call the collective variable module
563  if (colvars->calc() != COLVARS_OK) {
564  cvm::error("Error in the collective variables module.\n", COLVARS_ERROR);
565  }
566 
567  if (cvm::debug()) {
568  print_output_atomic_data();
569  }
570 
571  // communicate all forces to the MD integrator
572  for (size_t i = 0; i < atoms_ids.size(); i++) {
573  cvm::rvector const &f = atoms_new_colvar_forces[i];
574  modifyForcedAtoms().add(atoms_ids[i]);
575  modifyAppliedForces().add(Vector(f.x, f.y, f.z));
576  }
577 
578  if (atom_groups_new_colvar_forces.size() > 0) {
581  for (int ig = 0; gf_i != modifyGroupForces().end(); gf_i++, ig++) {
582  cvm::rvector const &f = atom_groups_new_colvar_forces[ig];
583  *gf_i = Vector(f.x, f.y, f.z);
584  }
585  }
586 
587 #if NAMD_VERSION_NUMBER >= 34471681
588  if (volmaps_new_colvar_forces.size() > 0) {
593  for ( ; goi_i != getGridObjIndexEnd(); goi_i++, gof_i++) {
594  for (size_t imap = 0; imap < volmaps_ids.size(); imap++) {
595  if (volmaps_ids[imap] == *goi_i) {
596  *gof_i = volmaps_new_colvar_forces[imap];
597  break;
598  }
599  }
600  }
601  }
602 #endif
603 
604  // send MISC energy
605  #if defined(NODEGROUP_FORCE_REGISTER) && !defined(NAMD_UNIFIED_REDUCTION)
607  reduction->submit();
608  }
609  #else
610  #if !defined(NAMD_UNIFIED_REDUCTION)
611  reduction->submit();
612  #else
613  // submitReduction();
614  #endif
615  #endif
616 
617  // NAMD does not destruct GlobalMaster objects, so we must remember
618  // to write all output files at the end of a run
619  if (step == simparams->N) {
620  post_run();
621  }
622 }
623 
625  // This aMD factor is from previous step!
626  amd_weight_factor = std::exp(controller->accelMDdV / (target_temperature() * boltzmann()));
627 }
628 
629 
630 // Callback functions
631 
633 {
634 #ifdef NAMD_TCL
635  // Store pointer to NAMD's Tcl interpreter
636  set_tcl_interp(Node::Object()->getScript()->interp);
637 #else
638  colvarproxy::init_tcl_pointers(); // Create dedicated interpreter
639 #endif
640 }
641 
643 {
644  return colvarproxy::tcl_run_force_callback();
645 }
646 
648  std::string const &name,
649  std::vector<const colvarvalue *> const &cvc_values,
650  colvarvalue &value)
651 {
652  return colvarproxy::tcl_run_colvar_callback(name, cvc_values, value);
653 }
654 
656  std::string const &name,
657  std::vector<const colvarvalue *> const &cvc_values,
658  std::vector<cvm::matrix2d<cvm::real> > &gradient)
659 {
660  return colvarproxy::tcl_run_colvar_gradient_callback(name, cvc_values,
661  gradient);
662 }
663 
664 
665 void colvarproxy_namd::add_energy(cvm::real energy)
666 {
667  #if defined(NODEGROUP_FORCE_REGISTER) && !defined(NAMD_UNIFIED_REDUCTION)
669  nodeReduction->item(REDUCTION_MISC_ENERGY) += energy;
670  } else {
672  }
673  #else
674  #if !defined(NAMD_UNIFIED_REDUCTION)
676  #else
678  #endif
679  #endif
680 }
681 
683 {
684  if (cvm::debug()) {
685  cvm::log("colvarproxy_namd::request_total_force()\n");
686  }
687  total_force_requested = yesno;
688  requestTotalForce(total_force_requested);
689  if (cvm::debug()) {
690  cvm::log("colvarproxy_namd::request_total_force() end\n");
691  }
692 }
693 
694 
695 void colvarproxy_namd::log(std::string const &message)
696 {
697  std::istringstream is(message);
698  std::string line;
699  while (std::getline(is, line))
700  iout << "colvars: " << line << "\n";
701  iout << endi;
702 }
703 
704 
705 void colvarproxy_namd::error(std::string const &message)
706 {
707  log(message);
708  switch (cvm::get_error()) {
709  case COLVARS_FILE_ERROR:
710  errno = EIO; break;
711  case COLVARS_NOT_IMPLEMENTED:
712  errno = ENOSYS; break;
713  case COLVARS_MEMORY_ERROR:
714  errno = ENOMEM; break;
715  }
716  char const *msg = "Error in the collective variables module "
717  "(see above for details)";
718  if (errno) {
719  NAMD_err(msg);
720  } else {
721  NAMD_die(msg);
722  }
723 }
724 
725 
727 {
728  // NAMD's internal numbering starts from zero
729  int const aid = (atom_number-1);
730 
731  if (cvm::debug())
732  cvm::log("Adding atom "+cvm::to_str(atom_number)+
733  " for collective variables calculation.\n");
734 
735  if ( (aid < 0) || (aid >= Node::Object()->molecule->numAtoms) ) {
736  cvm::error("Error: invalid atom number specified, "+
737  cvm::to_str(atom_number)+"\n", COLVARS_INPUT_ERROR);
738  return COLVARS_INPUT_ERROR;
739  }
740 
741  return aid;
742 }
743 
744 
746 {
747  return COLVARS_OK;
748 }
749 
750 
751 int colvarproxy_namd::init_atom(int atom_number)
752 {
753  // save time by checking first whether this atom has been requested before
754  // (this is more common than a non-valid atom number)
755  int aid = (atom_number-1);
756 
757  for (size_t i = 0; i < atoms_ids.size(); i++) {
758  if (atoms_ids[i] == aid) {
759  // this atom id was already recorded
760  atoms_refcount[i] += 1;
761  return i;
762  }
763  }
764 
765  aid = check_atom_id(atom_number);
766 
767  if (aid < 0) {
768  return COLVARS_INPUT_ERROR;
769  }
770 
771  int const index = add_atom_slot(aid);
772  atoms_map[aid] = index;
773  modifyRequestedAtoms().add(aid);
774  update_atom_properties(index);
775  return index;
776 }
777 
778 
779 int colvarproxy_namd::check_atom_id(cvm::residue_id const &residue,
780  std::string const &atom_name,
781  std::string const &segment_id)
782 {
783  int const aid =
784  (segment_id.size() ?
785  Node::Object()->molecule->get_atom_from_name(segment_id.c_str(),
786  residue,
787  atom_name.c_str()) :
789  residue,
790  atom_name.c_str()));
791 
792  if (aid < 0) {
793  // get_atom_from_name() has returned an error value
794  cvm::error("Error: could not find atom \""+
795  atom_name+"\" in residue "+
796  cvm::to_str(residue)+
797  ( (segment_id != "MAIN") ?
798  (", segment \""+segment_id+"\"") :
799  ("") )+
800  "\n", COLVARS_INPUT_ERROR);
801  return COLVARS_INPUT_ERROR;
802  }
803 
804  return aid;
805 }
806 
807 
808 
812 int colvarproxy_namd::init_atom(cvm::residue_id const &residue,
813  std::string const &atom_name,
814  std::string const &segment_id)
815 {
816  int const aid = check_atom_id(residue, atom_name, segment_id);
817 
818  for (size_t i = 0; i < atoms_ids.size(); i++) {
819  if (atoms_ids[i] == aid) {
820  // this atom id was already recorded
821  atoms_refcount[i] += 1;
822  return i;
823  }
824  }
825 
826  if (cvm::debug())
827  cvm::log("Adding atom \""+
828  atom_name+"\" in residue "+
829  cvm::to_str(residue)+
830  " (index "+cvm::to_str(aid)+
831  ") for collective variables calculation.\n");
832 
833  int const index = add_atom_slot(aid);
834  atoms_map[aid] = index;
835  modifyRequestedAtoms().add(aid);
836  update_atom_properties(index);
837  return index;
838 }
839 
840 
842 {
843  colvarproxy::clear_atom(index);
844  // TODO remove it from GlobalMaster arrays?
845 }
846 
847 
849 {
850  Molecule *mol = Node::Object()->molecule;
851  // update mass
852  double const mass = mol->atommass(atoms_ids[index]);
853  if (mass <= 0.001) {
854  this->log("Warning: near-zero mass for atom "+
855  cvm::to_str(atoms_ids[index]+1)+
856  "; expect unstable dynamics if you apply forces to it.\n");
857  }
858  atoms_masses[index] = mass;
859  // update charge
860  atoms_charges[index] = mol->atomcharge(atoms_ids[index]);
861 }
862 
863 
864 cvm::rvector colvarproxy_namd::position_distance(cvm::atom_pos const &pos1,
865  cvm::atom_pos const &pos2)
866  const
867 {
868  Position const p1(pos1.x, pos1.y, pos1.z);
869  Position const p2(pos2.x, pos2.y, pos2.z);
870  // return p2 - p1
871  Vector const d = this->lattice->delta(p2, p1);
872  return cvm::rvector(d.x, d.y, d.z);
873 }
874 
875 
876 
885 };
886 
887 
888 e_pdb_field pdb_field_str2enum(std::string const &pdb_field_str)
889 {
890  e_pdb_field pdb_field = e_pdb_none;
891 
892  if (colvarparse::to_lower_cppstr(pdb_field_str) ==
893  colvarparse::to_lower_cppstr("O")) {
894  pdb_field = e_pdb_occ;
895  }
896 
897  if (colvarparse::to_lower_cppstr(pdb_field_str) ==
898  colvarparse::to_lower_cppstr("B")) {
899  pdb_field = e_pdb_beta;
900  }
901 
902  if (colvarparse::to_lower_cppstr(pdb_field_str) ==
903  colvarparse::to_lower_cppstr("X")) {
904  pdb_field = e_pdb_x;
905  }
906 
907  if (colvarparse::to_lower_cppstr(pdb_field_str) ==
908  colvarparse::to_lower_cppstr("Y")) {
909  pdb_field = e_pdb_y;
910  }
911 
912  if (colvarparse::to_lower_cppstr(pdb_field_str) ==
913  colvarparse::to_lower_cppstr("Z")) {
914  pdb_field = e_pdb_z;
915  }
916 
917  if (pdb_field == e_pdb_none) {
918  cvm::error("Error: unsupported PDB field, \""+
919  pdb_field_str+"\".\n", COLVARS_INPUT_ERROR);
920  }
921 
922  return pdb_field;
923 }
924 
925 
926 int colvarproxy_namd::load_coords_pdb(char const *pdb_filename,
927  std::vector<cvm::atom_pos> &pos,
928  const std::vector<int> &indices,
929  std::string const &pdb_field_str,
930  double const pdb_field_value)
931 {
932  if (pdb_field_str.size() == 0 && indices.size() == 0) {
933  cvm::error("Bug alert: either PDB field should be defined or list of "
934  "atom IDs should be available when loading atom coordinates!\n", COLVARS_BUG_ERROR);
935  }
936 
937  e_pdb_field pdb_field_index;
938  bool const use_pdb_field = (pdb_field_str.size() > 0);
939  if (use_pdb_field) {
940  pdb_field_index = pdb_field_str2enum(pdb_field_str);
941  }
942 
943  // next index to be looked up in PDB file (if list is supplied)
944  std::vector<int>::const_iterator current_index = indices.begin();
945 
946  PDB *pdb = new PDB(pdb_filename);
947  size_t const pdb_natoms = pdb->num_atoms();
948 
949  if (pos.size() != pdb_natoms) {
950 
951  bool const pos_allocated = (pos.size() > 0);
952 
953  size_t ipos = 0, ipdb = 0;
954  for ( ; ipdb < pdb_natoms; ipdb++) {
955 
956  if (use_pdb_field) {
957  // PDB field mode: skip atoms with wrong value in PDB field
958  double atom_pdb_field_value = 0.0;
959 
960  switch (pdb_field_index) {
961  case e_pdb_occ:
962  atom_pdb_field_value = (pdb->atom(ipdb))->occupancy();
963  break;
964  case e_pdb_beta:
965  atom_pdb_field_value = (pdb->atom(ipdb))->temperaturefactor();
966  break;
967  case e_pdb_x:
968  atom_pdb_field_value = (pdb->atom(ipdb))->xcoor();
969  break;
970  case e_pdb_y:
971  atom_pdb_field_value = (pdb->atom(ipdb))->ycoor();
972  break;
973  case e_pdb_z:
974  atom_pdb_field_value = (pdb->atom(ipdb))->zcoor();
975  break;
976  default:
977  break;
978  }
979 
980  if ( (pdb_field_value) &&
981  (atom_pdb_field_value != pdb_field_value) ) {
982  continue;
983  } else if (atom_pdb_field_value == 0.0) {
984  continue;
985  }
986 
987  } else {
988  // Atom ID mode: use predefined atom IDs from the atom group
989  if (((int) ipdb) != *current_index) {
990  // Skip atoms not in the list
991  continue;
992  } else {
993  current_index++;
994  }
995  }
996 
997  if (!pos_allocated) {
998  pos.push_back(cvm::atom_pos(0.0, 0.0, 0.0));
999  } else if (ipos >= pos.size()) {
1000  cvm::error("Error: the PDB file \""+
1001  std::string(pdb_filename)+
1002  "\" contains coordinates for "
1003  "more atoms than needed.\n", COLVARS_BUG_ERROR);
1004  }
1005 
1006  pos[ipos] = cvm::atom_pos((pdb->atom(ipdb))->xcoor(),
1007  (pdb->atom(ipdb))->ycoor(),
1008  (pdb->atom(ipdb))->zcoor());
1009  ipos++;
1010  if (!use_pdb_field && current_index == indices.end())
1011  break;
1012  }
1013 
1014  if (ipos < pos.size() || (!use_pdb_field && current_index != indices.end())) {
1015  size_t n_requested = use_pdb_field ? pos.size() : indices.size();
1016  cvm::error("Error: number of matching records in the PDB file \""+
1017  std::string(pdb_filename)+"\" ("+cvm::to_str(ipos)+
1018  ") does not match the number of requested coordinates ("+
1019  cvm::to_str(n_requested)+").\n", COLVARS_INPUT_ERROR);
1020  return COLVARS_ERROR;
1021  }
1022  } else {
1023 
1024  // when the PDB contains exactly the number of atoms of the array,
1025  // ignore the fields and just read coordinates
1026  for (size_t ia = 0; ia < pos.size(); ia++) {
1027  pos[ia] = cvm::atom_pos((pdb->atom(ia))->xcoor(),
1028  (pdb->atom(ia))->ycoor(),
1029  (pdb->atom(ia))->zcoor());
1030  }
1031  }
1032 
1033  delete pdb;
1034  return COLVARS_OK;
1035 }
1036 
1037 int colvarproxy_namd::load_atoms_pdb(char const *pdb_filename,
1038  cvm::atom_group& atoms,
1039  std::string const &pdb_field_str,
1040  double const pdb_field_value)
1041 {
1042  if (pdb_field_str.size() == 0)
1043  cvm::error("Error: must define which PDB field to use "
1044  "in order to define atoms from a PDB file.\n", COLVARS_INPUT_ERROR);
1045 
1046  PDB *pdb = new PDB(pdb_filename);
1047  size_t const pdb_natoms = pdb->num_atoms();
1048 
1049  e_pdb_field pdb_field_index = pdb_field_str2enum(pdb_field_str);
1050  auto modify_atoms = atoms.get_atom_modifier();
1051  for (size_t ipdb = 0; ipdb < pdb_natoms; ipdb++) {
1052 
1053  double atom_pdb_field_value = 0.0;
1054 
1055  switch (pdb_field_index) {
1056  case e_pdb_occ:
1057  atom_pdb_field_value = (pdb->atom(ipdb))->occupancy();
1058  break;
1059  case e_pdb_beta:
1060  atom_pdb_field_value = (pdb->atom(ipdb))->temperaturefactor();
1061  break;
1062  case e_pdb_x:
1063  atom_pdb_field_value = (pdb->atom(ipdb))->xcoor();
1064  break;
1065  case e_pdb_y:
1066  atom_pdb_field_value = (pdb->atom(ipdb))->ycoor();
1067  break;
1068  case e_pdb_z:
1069  atom_pdb_field_value = (pdb->atom(ipdb))->zcoor();
1070  break;
1071  default:
1072  break;
1073  }
1074 
1075  if ( (pdb_field_value) &&
1076  (atom_pdb_field_value != pdb_field_value) ) {
1077  continue;
1078  } else if (atom_pdb_field_value == 0.0) {
1079  continue;
1080  }
1081 
1082  if (atoms.is_enabled(colvardeps::f_ag_scalable)) {
1083  modify_atoms.add_atom_id(ipdb);
1084  } else {
1085  modify_atoms.add_atom(cvm::atom_group::init_atom_from_proxy(this, ipdb+1));
1086  }
1087  }
1088 
1089  delete pdb;
1090  return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
1091 }
1092 
1093 
1094 std::ostream & colvarproxy_namd::output_stream(std::string const &output_name,
1095  std::string const description)
1096 {
1097  if (cvm::debug()) {
1098  cvm::log("Using colvarproxy_namd::output_stream()\n");
1099  }
1100 
1101  if (!io_available()) {
1102  cvm::error("Error: trying to access an output file/channel "
1103  "from the wrong thread.\n", COLVARS_BUG_ERROR);
1104  return *output_stream_error_;
1105  }
1106 
1107  if (output_streams_.count(output_name) > 0) {
1108  return *(output_streams_[output_name]);
1109  }
1110 
1111  backup_file(output_name.c_str());
1112 
1113  output_streams_[output_name] = new ofstream_namd(output_name.c_str(), std::ios::binary);
1114  if (! output_streams_[output_name]->good()) {
1115  cvm::error("Error: cannot write to "+description+" \""+output_name+"\".\n",
1116  COLVARS_FILE_ERROR);
1117  }
1118 
1119  return *(output_streams_[output_name]);
1120 }
1121 
1122 
1123 int colvarproxy_namd::flush_output_stream(std::string const &output_name)
1124 {
1125  if (!io_available()) {
1126  return COLVARS_OK;
1127  }
1128 
1129  if (output_streams_.count(output_name) > 0) {
1130  (reinterpret_cast<ofstream_namd *>(output_streams_[output_name]))->flush();
1131  return COLVARS_OK;
1132  }
1133 
1134  return cvm::error("Error: trying to flush an output file/channel "
1135  "that wasn't open.\n", COLVARS_BUG_ERROR);
1136 }
1137 
1138 
1140 {
1141  if (!io_available()) {
1142  return COLVARS_OK;
1143  }
1144 
1145  for (std::map<std::string, std::ostream *>::iterator osi = output_streams_.begin();
1146  osi != output_streams_.end();
1147  osi++) {
1148  (reinterpret_cast<ofstream_namd *>(osi->second))->flush();
1149  }
1150 
1151  return COLVARS_OK;
1152 }
1153 
1154 
1155 int colvarproxy_namd::close_output_stream(std::string const &output_name)
1156 {
1157  if (!io_available()) {
1158  return cvm::error("Error: trying to access an output file/channel "
1159  "from the wrong thread.\n", COLVARS_BUG_ERROR);
1160  }
1161 
1162  if (output_streams_.count(output_name) > 0) {
1163  (reinterpret_cast<ofstream_namd *>(output_streams_[output_name]))->close();
1164  delete output_streams_[output_name];
1165  output_streams_.erase(output_name);
1166  }
1167 
1168  return COLVARS_OK;
1169 }
1170 
1171 
1173 {
1174  if (! io_available()) {
1175  return COLVARS_OK;
1176  }
1177 
1178  for (std::map<std::string, std::ostream *>::iterator osi = output_streams_.begin();
1179  osi != output_streams_.end();
1180  osi++) {
1181  (reinterpret_cast<ofstream_namd *>(osi->second))->close();
1182  }
1183  output_streams_.clear();
1184 
1185  return COLVARS_OK;
1186 }
1187 
1188 
1189 int colvarproxy_namd::backup_file(char const *filename)
1190 {
1191  if (std::string(filename).rfind(std::string(".colvars.state")) != std::string::npos) {
1192  NAMD_backup_file(filename, ".old");
1193  } else {
1194  NAMD_backup_file(filename, ".BAK");
1195  }
1196  return COLVARS_OK;
1197 }
1198 
1199 
1200 int colvarproxy_namd::init_atom_group(std::vector<int> const &atoms_ids)
1201 {
1202  if (cvm::debug())
1203  cvm::log("Requesting from NAMD a group of size "+cvm::to_str(atoms_ids.size())+
1204  " for collective variables calculation.\n");
1205 
1206  colvars->cite_feature("Scalable center-of-mass computation (NAMD)");
1207 
1208  // Note: modifyRequestedGroups is supposed to be in sync with the colvarproxy arrays,
1209  // and to stay that way during a simulation
1210 
1211  // compare this new group to those already allocated inside GlobalMaster
1212  int ig;
1213  for (ig = 0; ig < modifyRequestedGroups().size(); ig++) {
1214  AtomIDList const &namd_group = modifyRequestedGroups()[ig];
1215  bool b_match = true;
1216 
1217  if (namd_group.size() != ((int) atoms_ids.size())) {
1218  b_match = false;
1219  } else {
1220  int ia;
1221  for (ia = 0; ia < namd_group.size(); ia++) {
1222  int const aid = atoms_ids[ia];
1223  if (namd_group[ia] != aid) {
1224  b_match = false;
1225  break;
1226  }
1227  }
1228  }
1229 
1230  if (b_match) {
1231  if (cvm::debug())
1232  cvm::log("Group was already added.\n");
1233  // this group already exists
1234  atom_groups_refcount[ig] += 1;
1235  return ig;
1236  }
1237  }
1238 
1239  // add this group (note: the argument of add_atom_group_slot() is redundant for NAMD, and provided only for consistency)
1240  size_t const index = add_atom_group_slot(atom_groups_ids.size());
1241  modifyRequestedGroups().resize(atom_groups_ids.size());
1242  // the following is done in calculate()
1243  // modifyGroupForces().resize(atom_groups_ids.size());
1244  AtomIDList &namd_group = modifyRequestedGroups()[index];
1245  namd_group.resize(atoms_ids.size());
1246  int const n_all_atoms = Node::Object()->molecule->numAtoms;
1247  for (size_t ia = 0; ia < atoms_ids.size(); ia++) {
1248  int const aid = atoms_ids[ia];
1249  if (cvm::debug())
1250  cvm::log("Adding atom "+cvm::to_str(aid+1)+
1251  " for collective variables calculation.\n");
1252  if ( (aid < 0) || (aid >= n_all_atoms) ) {
1253  cvm::error("Error: invalid atom number specified, "+
1254  cvm::to_str(aid+1)+"\n", COLVARS_INPUT_ERROR);
1255  return -1;
1256  }
1257  namd_group[ia] = aid;
1258  }
1259 
1260  update_group_properties(index);
1261 
1262  if (cvm::debug()) {
1263  cvm::log("Group has index "+cvm::to_str(index)+"\n");
1264  cvm::log("modifyRequestedGroups length = "+cvm::to_str(modifyRequestedGroups().size())+
1265  ", modifyGroupForces length = "+cvm::to_str(modifyGroupForces().size())+"\n");
1266  }
1267 
1268  return index;
1269 }
1270 
1271 
1273 {
1274  // do nothing, keep the NAMD arrays in sync with the colvarproxy ones
1275  colvarproxy::clear_atom_group(index);
1276 }
1277 
1278 
1280 {
1281  AtomIDList const &namd_group = modifyRequestedGroups()[index];
1282  if (cvm::debug()) {
1283  cvm::log("Re-calculating total mass and charge for scalable group no. "+cvm::to_str(index+1)+" ("+
1284  cvm::to_str(namd_group.size())+" atoms).\n");
1285  }
1286 
1287  cvm::real total_mass = 0.0;
1288  cvm::real total_charge = 0.0;
1289  for (int i = 0; i < namd_group.size(); i++) {
1290  total_mass += Node::Object()->molecule->atommass(namd_group[i]);
1291  total_charge += Node::Object()->molecule->atomcharge(namd_group[i]);
1292  }
1293  atom_groups_masses[index] = total_mass;
1294  atom_groups_charges[index] = total_charge;
1295 
1296  if (cvm::debug()) {
1297  cvm::log("total mass = "+cvm::to_str(total_mass)+
1298  ", total charge = "+cvm::to_str(total_charge)+"\n");
1299  }
1300 
1301  return COLVARS_OK;
1302 }
1303 
1304 
1305 
1306 int colvarproxy_namd::set_unit_system(std::string const &units_in, bool /*check_only*/)
1307 {
1308  if (units_in != "real") {
1309  cvm::error("Error: Specified unit system \"" + units_in + "\" is unsupported in NAMD. Supported units are \"real\" (A, kcal/mol).\n");
1310  return COLVARS_ERROR;
1311  }
1312  return COLVARS_OK;
1313 }
1314 
1315 
1316 #if NAMD_VERSION_NUMBER >= 34471681
1317 
1318 
1320 {
1321  return COLVARS_OK;
1322 }
1323 
1324 
1326 {
1327  for (size_t i = 0; i < volmaps_ids.size(); i++) {
1328  if (volmaps_ids[i] == volmap_id) {
1329  // this map has already been requested
1330  volmaps_refcount[i] += 1;
1331  return i;
1332  }
1333  }
1334 
1335  int error_code = check_volmap_by_id(volmap_id);
1336  int index = -1;
1337  if (error_code == COLVARS_OK) {
1338  index = add_volmap_slot(volmap_id);
1339  modifyRequestedGridObjects().add(volmap_id);
1340  }
1341 
1342  return (error_code == COLVARS_OK) ? index : -1;
1343 }
1344 
1345 
1346 int colvarproxy_namd::init_volmap_by_name(char const *volmap_name)
1347 {
1348  if (volmap_name == NULL) {
1349  return cvm::error("Error: no grid object name provided.", COLVARS_INPUT_ERROR);
1350  }
1351 
1352  int error_code = COLVARS_OK;
1353 
1354  error_code |= check_volmap_by_name(volmap_name);
1355 
1356  int index = -1;
1357  if (error_code == COLVARS_OK) {
1358 
1359  int volmap_id = simparams->mgridforcelist.index_for_key(volmap_name);
1360 
1361  // Check that the scale factor is correctly set to zero
1362  Molecule *mol = Node::Object()->molecule;
1363  GridforceGrid const *grid = mol->get_gridfrc_grid(volmap_id);
1364  Vector const gfScale = grid->get_scale();
1365  if ((gfScale.x != 0.0) || (gfScale.y != 0.0) || (gfScale.z != 0.0)) {
1366  error_code |= cvm::error("Error: GridForce map \""+
1367  std::string(volmap_name)+
1368  "\" has non-zero scale factors.\n",
1369  COLVARS_INPUT_ERROR);
1370  }
1371 
1372  for (size_t i = 0; i < volmaps_ids.size(); i++) {
1373  if (volmaps_ids[i] == volmap_id) {
1374  // this map has already been requested
1375  volmaps_refcount[i] += 1;
1376  return i;
1377  }
1378  }
1379 
1380  index = add_volmap_slot(volmap_id);
1381  modifyRequestedGridObjects().add(volmap_id);
1382  }
1383 
1384  return (error_code == COLVARS_OK) ? index : -1;
1385 }
1386 
1387 
1389 {
1390  Molecule *mol = Node::Object()->molecule;
1391  if ((volmap_id < 0) || (volmap_id >= mol->numGridforceGrids)) {
1392  return cvm::error("Error: invalid numeric ID ("+cvm::to_str(volmap_id)+
1393  ") for map.\n", COLVARS_INPUT_ERROR);
1394  }
1395  colvars->cite_feature("GridForces volumetric map implementation for NAMD");
1396  return COLVARS_OK;
1397 }
1398 
1399 
1400 int colvarproxy_namd::check_volmap_by_name(char const *volmap_name)
1401 {
1402  if (volmap_name == NULL) {
1403  return cvm::error("Error: no grid object name provided.", COLVARS_INPUT_ERROR);
1404  }
1405  int volmap_id = simparams->mgridforcelist.index_for_key(volmap_name);
1406  if (volmap_id < 0) {
1407  return cvm::error("Error: invalid map name \""+std::string(volmap_name)+
1408  "\".\n", COLVARS_INPUT_ERROR);
1409  }
1410  colvars->cite_feature("GridForces volumetric map implementation for NAMD");
1411  return COLVARS_OK;
1412 }
1413 
1414 
1416 {
1417  // TODO remove from GlobalMaster
1418  colvarproxy::clear_volmap(index);
1419 }
1420 
1421 
1422 int colvarproxy_namd::get_volmap_id_from_name(char const *volmap_name)
1423 {
1424  int const volmap_id =
1425  simparams->mgridforcelist.index_for_key(volmap_name);
1426  if (volmap_id < 0) {
1427  // Print error
1428  check_volmap_by_name(volmap_name);
1429  }
1430  return volmap_id;
1431 }
1432 
1433 
1434 template<class T, int flags>
1436  cvm::atom_group* ag,
1437  cvm::real *value,
1438  cvm::real *atom_field)
1439 {
1440  float V = 0.0f;
1441  Vector dV(0.0);
1442  for (size_t i = 0; i < ag->size(); ++i) {
1443  if (g->compute_VdV(Position(ag->pos_x(i), ag->pos_y(i), ag->pos_z(i)), V, dV)) {
1444  // out-of-bounds atom
1445  V = 0.0f;
1446  dV = 0.0;
1447  } else {
1448  if (flags & volmap_flag_use_atom_field) {
1449  *value += V * atom_field[i];
1450  if (flags & volmap_flag_gradients) {
1451  const cvm::rvector grad = atom_field[i] * cvm::rvector(dV.x, dV.y, dV.z);
1452  ag->grad_x(i) += grad.x;
1453  ag->grad_y(i) += grad.y;
1454  ag->grad_z(i) += grad.z;
1455  }
1456  } else {
1457  *value += V;
1458  if (flags & volmap_flag_gradients) {
1459  ag->grad_x(i) += dV.x;
1460  ag->grad_y(i) += dV.y;
1461  ag->grad_z(i) += dV.z;
1462  }
1463  }
1464  }
1465  }
1466 }
1467 
1468 
1469 template<class T>
1471  T const *g,
1472  cvm::atom_group* ag,
1473  cvm::real *value,
1474  cvm::real *atom_field)
1475 {
1476  if (flags & volmap_flag_use_atom_field) {
1477  int const new_flags = volmap_flag_use_atom_field | volmap_flag_gradients;
1478  GridForceGridLoop<T, new_flags>(g, ag, value, atom_field);
1479  } else {
1480  int const new_flags = volmap_flag_gradients;
1481  GridForceGridLoop<T, new_flags>(g, ag, value, atom_field);
1482  }
1483 }
1484 
1486  int volmap_id,
1487  cvm::atom_group* ag,
1488  cvm::real *value,
1489  cvm::real *atom_field)
1490 {
1491  Molecule *mol = Node::Object()->molecule;
1492  GridforceGrid *grid = mol->get_gridfrc_grid(volmap_id);
1493  // Inheritance is not possible with GridForceGrid's design
1495  GridforceFullMainGrid *g = dynamic_cast<GridforceFullMainGrid *>(grid);
1496  getGridForceGridValue<GridforceFullMainGrid>(flags, g, ag,
1497  value, atom_field);
1498  } else if (grid->get_grid_type() == GridforceGrid::GridforceGridTypeLite) {
1499  GridforceLiteGrid *g = dynamic_cast<GridforceLiteGrid *>(grid);
1500  getGridForceGridValue<GridforceLiteGrid>(flags, g, ag,
1501  value, atom_field);
1502  }
1503  return COLVARS_OK;
1504 }
1505 
1506 #endif
1507 
1508 #if CMK_SMP && USE_CKLOOP // SMP only
1509 
1510 colvarproxy::smp_mode_t colvarproxy_namd::get_smp_mode() const {
1511  return smp_mode;
1512 }
1513 
1514 int colvarproxy_namd::set_smp_mode(smp_mode_t mode) {
1515  smp_mode = mode;
1516  return COLVARS_OK;
1517 }
1518 
1519 
1520 int colvarproxy_namd::smp_loop(int n_items, std::function<int (int)> const &worker)
1521 {
1522  auto cmkWorker = [&](int start, int end, void * /* result */) {
1523 #if CMK_TRACE_ENABLED
1524  double before = CmiWallTimer();
1525 #endif
1526  for (int i = start; i <= end; i++) {
1527  worker(i);
1528  }
1529 #if CMK_TRACE_ENABLED
1530  traceUserBracketEvent(GLOBAL_MASTER_CKLOOP_CALC_ITEM, before, CmiWallTimer());
1531 #endif
1532  };
1533  const int numChunks = smp_num_threads() > n_items ?
1534  n_items :
1535  smp_num_threads();
1536  cvm::increase_depth();
1537  CkLoop_Parallelize(numChunks, 0, n_items - 1, cmkWorker, nullptr, CKLOOP_NONE, nullptr);
1538  cvm::decrease_depth();
1539  // CkLoop does not support bitwise-OR reduction, so we just return the global error flag
1540  return cvm::get_error();
1541 }
1542 
1543 
1544 void calc_cv_biases_smp(int first, int last, void *result, int paramNum, void *param)
1545 {
1546  colvarproxy_namd *proxy = (colvarproxy_namd *) param;
1547  colvarmodule *cv = proxy->colvars;
1548 #if CMK_TRACE_ENABLED
1549  double before = CmiWallTimer();
1550 #endif
1551  cvm::increase_depth();
1552  for (int i = first; i <= last; i++) {
1553  colvarbias *b = (*(cv->biases_active()))[i];
1554  if (cvm::debug()) {
1555  cvm::log("["+cvm::to_str(proxy->smp_thread_id())+"/"+cvm::to_str(proxy->smp_num_threads())+
1556  "]: calc_cv_biases_smp(), first = "+cvm::to_str(first)+
1557  ", last = "+cvm::to_str(last)+", bias = "+
1558  b->name+"\n");
1559  }
1560  b->update();
1561  }
1562  cvm::decrease_depth();
1563 #if CMK_TRACE_ENABLED
1564  traceUserBracketEvent(GLOBAL_MASTER_CKLOOP_CALC_BIASES,before,CmiWallTimer());
1565 #endif
1566 }
1567 
1568 
1569 int colvarproxy_namd::smp_biases_loop()
1570 {
1571  colvarmodule *cv = this->colvars;
1572  const int numChunks = smp_num_threads() > cv->variables_active_smp()->size() ?
1573  cv->variables_active_smp()->size() :
1574  smp_num_threads();
1575  CkLoop_Parallelize(calc_cv_biases_smp, 1, this,
1576  numChunks, 0, cv->biases_active()->size()-1);
1577  return cvm::get_error();
1578 }
1579 
1580 
1581 void calc_cv_scripted_forces(int paramNum, void *param)
1582 {
1583  colvarproxy_namd *proxy = (colvarproxy_namd *) param;
1584  colvarmodule *cv = proxy->colvars;
1585 #if CMK_TRACE_ENABLED
1586  double before = CmiWallTimer();
1587 #endif
1588  if (cvm::debug()) {
1589  cvm::log("["+cvm::to_str(proxy->smp_thread_id())+"/"+cvm::to_str(proxy->smp_num_threads())+
1590  "]: calc_cv_scripted_forces()\n");
1591  }
1592  cv->calc_scripted_forces();
1593 #if CMK_TRACE_ENABLED
1594  traceUserBracketEvent(GLOBAL_MASTER_CKLOOP_CALC_SCRIPTED_BIASES,before,CmiWallTimer());
1595 #endif
1596 }
1597 
1598 
1599 int colvarproxy_namd::smp_biases_script_loop()
1600 {
1601  colvarmodule *cv = this->colvars;
1602  CkLoop_Parallelize(calc_cv_biases_smp, 1, this,
1603  cv->biases_active()->size(), 0, cv->biases_active()->size()-1,
1604  1, NULL, CKLOOP_NONE,
1605  calc_cv_scripted_forces, 1, this);
1606  return cvm::get_error();
1607 }
1608 
1609 #endif // #if CMK_SMP && USE_CKLOOP
1610 
1611 
1613 #if CMK_HAS_PARTITION
1614  return COLVARS_OK;
1615 #else
1616  return COLVARS_NOT_IMPLEMENTED;
1617 #endif
1618 }
1619 
1620 
1622  return CmiMyPartition();
1623 }
1624 
1625 
1627  return CmiNumPartitions();
1628 }
1629 
1630 
1632  replica_barrier();
1633 }
1634 
1635 
1636 int colvarproxy_namd::replica_comm_recv(char* msg_data, int buf_len,
1637  int src_rep) {
1638  DataMessage *recvMsg = NULL;
1639  replica_recv(&recvMsg, src_rep, CkMyPe());
1640  CmiAssert(recvMsg != NULL);
1641  int retval = recvMsg->size;
1642  if (buf_len >= retval) {
1643  memcpy(msg_data,recvMsg->data,retval);
1644  } else {
1645  retval = 0;
1646  }
1647  CmiFree(recvMsg);
1648  return retval;
1649 }
1650 
1651 
1652 int colvarproxy_namd::replica_comm_send(char* msg_data, int msg_len,
1653  int dest_rep) {
1654  replica_send(msg_data, msg_len, dest_rep, CkMyPe());
1655  return msg_len;
1656 }
1657 
1658 
1661  // This test is only valid for NAMD3
1662  if (freq % simparams->computeEnergies) {
1663  cvm::error("computeEnergies must be a divisor of lambda-dynamics period (" + cvm::to_str(freq) + ").\n");
1664  return COLVARS_INPUT_ERROR;
1665  }
1666  if (!simparams->alchOn) {
1667  cvm::error("alchOn must be enabled for lambda-dynamics.\n");
1668  return COLVARS_INPUT_ERROR;
1669  }
1670  if (!simparams->alchThermIntOn) {
1671  cvm::error("alchType must be set to TI for lambda-dynamics.\n");
1672  return COLVARS_INPUT_ERROR;
1673  }
1674  return COLVARS_OK;
1675 }
1676 
1677 
1679 int colvarproxy_namd::get_alch_lambda(cvm::real* lambda) {
1680  *lambda = simparams->getCurrentLambda(step);
1681  return COLVARS_OK;
1682 }
1683 
1684 
1687  if (simparams->alchLambdaFreq > 0) {
1688  cvm::error("Cannot set lambda from Colvars because alchLambdaFreq is enabled. "
1689  "Either remove biasing forces and extended Lagrangian dynamics on the alchemical coordinate, "
1690  "or disable alchLambdaFreq.\n");
1691  return COLVARS_INPUT_ERROR;
1692  } else {
1693  simparams->alchLambda = cached_alch_lambda;
1694  }
1695  return COLVARS_OK;
1696 }
1697 
1698 
1700 int colvarproxy_namd::get_dE_dlambda(cvm::real* dE_dlambda) {
1701  // Force data at step zero is garbage in NAMD3, zero in NAMD2
1702  if (cvm::step_relative() > 0) {
1703  *dE_dlambda = controller->getTIderivative();
1704  } else {
1705  *dE_dlambda = 0.0;
1706  }
1707  return COLVARS_OK;
1708 }
static Node * Object()
Definition: Node.h:86
Real atomcharge(int anum) const
Definition: Molecule.h:1124
int check_volmaps_available() override
std::ostream & output_stream(std::string const &output_name, std::string const description) override
const Controller & getController() const
Definition: NamdState.h:51
void clear_atom_group(int index) override
ForceList & modifyAppliedForces()
Definition: GlobalMaster.C:203
int get_volmap_id_from_name(char const *volmap_name) override
AtomIDList & modifyRequestedAtoms()
Definition: GlobalMaster.C:168
void clear_volmap(int index) override
GridforceGrid * get_gridfrc_grid(int gridnum) const
Definition: Molecule.h:1368
BigRealList & modifyGridObjForces()
Definition: GlobalMaster.C:220
Definition: PDB.h:36
int size(void) const
Definition: ResizeArray.h:131
SimParameters * simparams
Pointer to the NAMD simulation input object.
NAMD_HOST_DEVICE Vector c() const
Definition: Lattice.h:270
Random random
NAMD-style PRNG object.
void NAMD_err(const char *err_msg)
Definition: common.C:170
BigReal getTIderivative(void) const
Definition: Controller.h:121
int init_atom_group(std::vector< int > const &atoms_ids) override
NAMD_HOST_DEVICE int c_p() const
Definition: Lattice.h:291
int set_unit_system(std::string const &units_in, bool check_only) override
#define GLOBAL_MASTER_CKLOOP_CALC_ITEM
void error(std::string const &message) override
int replica_comm_send(char *msg_data, int msg_len, int dest_rep) override
int init_atom(int atom_number) override
AtomIDList::const_iterator getForceIdEnd()
Definition: GlobalMaster.C:301
AtomIDList::const_iterator getAtomIdBegin()
Definition: GlobalMaster.C:231
PositionList::const_iterator getGroupPositionEnd()
Definition: GlobalMaster.C:247
Communication between colvars and NAMD (implementation of colvarproxy)
Controller const * controller
Pointer to Controller object.
int init_volmap_by_id(int volmap_id) override
Definition: Vector.h:72
virtual void submit(void)=0
SimParameters * simParameters
Definition: Node.h:181
int check_volmap_by_name(char const *volmap_name) override
void clear_atom(int index) override
cvm::real amd_weight_factor
e_pdb_field pdb_field_str2enum(std::string const &pdb_field_str)
BigReal reassignTemp
BigReal & item(int i)
Definition: ReductionMgr.h:336
void replica_send(const char *sndbuf, int sendcount, int destPart, int destPE)
int get_alch_lambda(cvm::real *lambda)
Get value of alchemical lambda parameter from back-end.
std::ostream & endi(std::ostream &s)
Definition: InfoStream.C:54
BigReal z
Definition: Vector.h:74
ForceList::const_iterator getGroupTotalForceBegin()
Definition: GlobalMaster.C:251
cvm::step_number previous_NAMD_step
int numGridforceGrids
Definition: Molecule.h:626
SubmitReduction * willSubmit(int setID, int size=-1)
Definition: ReductionMgr.C:368
int update_target_temperature()
Get the target temperature from the NAMD thermostats supported so far.
static int debug
Definition: parm.C:36
static ReductionMgr * Object(void)
Definition: ReductionMgr.h:290
#define iout
Definition: InfoStream.h:51
NAMD_HOST_DEVICE int orthogonal() const
Definition: Lattice.h:273
NodeReduction * reduction
Definition: PatchData.h:139
int num_atoms(void)
Definition: PDB.C:323
void clear()
Definition: ResizeArray.h:91
int add(const Elem &elem)
Definition: ResizeArray.h:101
PositionList::const_iterator getAtomPositionBegin()
Definition: GlobalMaster.C:239
BigReal alchLambda
int compute_volmap(int flags, int volmap_id, cvm::atom_group *ag, cvm::real *value, cvm::real *atom_field) override
Molecule stores the structural information for the system.
Definition: Molecule.h:174
NAMD_HOST_DEVICE int b_p() const
Definition: Lattice.h:290
const ResizeArray< AtomIDList > & requestedGroups()
Definition: GlobalMaster.C:190
void request_total_force(bool yesno) override
char outputFilename[NAMD_FILENAME_BUFFER_SIZE]
int setup() override
AtomIDList::const_iterator getForceIdBegin()
Definition: GlobalMaster.C:296
void resize(int i)
Definition: ResizeArray.h:84
int init_volmap_by_name(const char *volmap_name) override
void setall(const Elem &elem)
Definition: ResizeArray.h:94
int check_atom_id(int atom_number) override
#define COLVARPROXY_VERSION
int run_colvar_gradient_callback(std::string const &name, std::vector< const colvarvalue *> const &cvcs, std::vector< cvm::matrix2d< cvm::real > > &gradient) override
#define GLOBAL_MASTER_CKLOOP_CALC_SCRIPTED_BIASES
std::vector< int > atoms_map
Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data...
void replica_recv(DataMessage **precvMsg, int srcPart, int srcPE)
int check_atom_name_selections_available() override
void getGridForceGridValue(int flags, T const *grid, cvm::atom_group *ag, cvm::real *value, cvm::real *atom_field)
Abstraction of the two types of NAMD volumetric maps.
Definition: Random.h:37
int index_for_key(const char *key)
void log(std::string const &message) override
e_pdb_field
int request_alch_energy_freq(int const freq)
Request energy computation every freq steps.
IntList::const_iterator getGridObjIndexBegin()
Definition: GlobalMaster.C:259
int get_atom_from_name(const char *segid, int resid, const char *aname) const
Definition: Molecule.C:126
SubmitReduction * reduction
Definition: GlobalMaster.h:198
bool accelMDOn
Used to submit restraint energy as MISC.
const AtomID * const_iterator
Definition: ResizeArray.h:38
int backup_file(char const *filename) override
PDBAtom * atom(int place)
Definition: PDB.C:393
BigReal rescaleTemp
void init_atoms_map()
Allocate an atoms map with the same size as the NAMD topology.
int flush_output_streams() override
NamdState * state
Definition: Node.h:184
IntList & modifyRequestedGridObjects()
Definition: GlobalMaster.C:214
BigReal accelMDdV
Definition: Controller.h:116
int reset() override
#define GLOBAL_MASTER_CKLOOP_CALC_BIASES
void replica_comm_barrier() override
BigReal langevinTemp
int send_alch_lambda(void)
Set value of alchemical lambda parameter in back-end.
BigReal x
Definition: Vector.h:74
int check_replicas_enabled() override
int close_output_stream(std::string const &output_name) override
void replica_barrier()
int replica_index() override
NAMD_HOST_DEVICE int a_p() const
Definition: Lattice.h:289
AtomIDList & modifyForcedAtoms()
Definition: GlobalMaster.C:198
int numAtoms
Definition: Molecule.h:586
ResizeArray< AtomIDList > & modifyRequestedGroups()
Definition: GlobalMaster.C:225
int update_group_properties(int index)
void NAMD_die(const char *err_msg)
Definition: common.C:147
int load_atoms_pdb(char const *filename, cvm::atom_group &atoms, std::string const &pdb_field, double pdb_field_value) override
int flush_output_stream(std::string const &output_name) override
ForceList & modifyGroupForces()
Definition: GlobalMaster.C:208
Real atommass(int anum) const
Definition: Molecule.h:1114
ConfigList * configList
Definition: Node.h:182
NAMD_HOST_DEVICE Vector b() const
Definition: Lattice.h:269
MGridforceParamsList mgridforcelist
void init_tcl_pointers() override
void update_atom_properties(int index)
void NAMD_backup_file(const char *filename, const char *extension)
Definition: common.C:235
virtual Vector get_scale(void) const =0
char restartFilename[NAMD_FILENAME_BUFFER_SIZE]
Bool set(const char *s)
Definition: Vector.h:261
BigRealList::const_iterator getGridObjValueBegin()
Definition: GlobalMaster.C:267
cvm::rvector position_distance(cvm::atom_pos const &pos1, cvm::atom_pos const &pos2) const
int check_volmap_by_id(int volmap_id) override
IntList::const_iterator getGridObjIndexEnd()
Definition: GlobalMaster.C:263
const Lattice * lattice
Definition: GlobalMaster.h:147
void add_energy(cvm::real energy) override
int get_dE_dlambda(cvm::real *dE_dlambda)
Get energy derivative with respect to lambda.
int close_output_streams() override
unsigned int randomSeed
void requestTotalForce(bool yesno=true)
Definition: GlobalMaster.h:135
char data[1]
Definition: DataExchanger.h:23
StringList * next
Definition: ConfigList.h:49
iterator begin(void)
Definition: ResizeArray.h:36
AtomIDList::const_iterator getAtomIdEnd()
Definition: GlobalMaster.C:235
BigReal getCurrentLambda(const int) const
char * data
Definition: ConfigList.h:48
iterator end(void)
Definition: ResizeArray.h:37
BigReal y
Definition: Vector.h:74
colvarmodule * colvars
Definition: Node.h:187
BigReal tCoupleTemp
int run_colvar_callback(std::string const &name, std::vector< const colvarvalue *> const &cvcs, colvarvalue &value) override
void addReductionEnergy(int reductionTag, BigReal energy)
Definition: GlobalMaster.C:143
GridforceGridType get_grid_type(void)
Definition: GridForceGrid.h:64
int num_replicas() override
const IntList & requestedGridObjs()
Definition: GlobalMaster.C:194
int run_force_callback() override
int replica_comm_recv(char *msg_data, int buf_len, int src_rep) override
StringList * find(const char *name) const
Definition: ConfigList.C:341
NAMD_HOST_DEVICE Vector a() const
Definition: Lattice.h:268
int update_atoms_map(AtomIDList::const_iterator begin, AtomIDList::const_iterator end)
PositionList::const_iterator getGroupPositionBegin()
Definition: GlobalMaster.C:243
Molecule * molecule
Definition: Node.h:179
ForceList::const_iterator getTotalForce()
Definition: GlobalMaster.C:306
BigReal stochRescaleTemp
int globalMasterFrequency
int load_coords_pdb(char const *filename, std::vector< cvm::atom_pos > &pos, const std::vector< int > &indices, std::string const &pdb_field, double const pdb_field_value) override
BigRealList::const_iterator getGridObjValueEnd()
Definition: GlobalMaster.C:271
NAMD_HOST_DEVICE Vector delta(const Position &pos1, const Position &pos2) const
Definition: Lattice.h:149
BigReal loweAndersenTemp
void GridForceGridLoop(T const *g, cvm::atom_group *ag, cvm::real *value, cvm::real *atom_field)
Implementation of inner loop; allows for atom list computation and use.
ForceList::const_iterator getGroupTotalForceEnd()
Definition: GlobalMaster.C:255