NAMD
Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
colvarproxy_namd Class Reference

Communication between colvars and NAMD (implementation of colvarproxy) More...

#include <colvarproxy_namd.h>

Inheritance diagram for colvarproxy_namd:
GlobalMaster

Public Member Functions

void init_tcl_pointers () override
 
 colvarproxy_namd ()
 
 ~colvarproxy_namd ()
 
int setup () override
 
int reset () override
 
int update_target_temperature ()
 Get the target temperature from the NAMD thermostats supported so far. More...
 
void init_atoms_map ()
 Allocate an atoms map with the same size as the NAMD topology. More...
 
int update_atoms_map (AtomIDList::const_iterator begin, AtomIDList::const_iterator end)
 
void calculate ()
 
void log (std::string const &message) override
 
void error (std::string const &message) override
 
int set_unit_system (std::string const &units_in, bool check_only) override
 
void add_energy (cvm::real energy) override
 
void request_total_force (bool yesno) override
 
bool total_forces_enabled () const override
 
int run_force_callback () override
 
int run_colvar_callback (std::string const &name, std::vector< const colvarvalue *> const &cvcs, colvarvalue &value) override
 
int run_colvar_gradient_callback (std::string const &name, std::vector< const colvarvalue *> const &cvcs, std::vector< cvm::matrix2d< cvm::real > > &gradient) override
 
cvm::real rand_gaussian () override
 
cvm::real get_accelMD_factor () const override
 
bool accelMD_enabled () const override
 
int check_replicas_enabled () override
 
int replica_index () override
 
int num_replicas () override
 
void replica_comm_barrier () override
 
int replica_comm_recv (char *msg_data, int buf_len, int src_rep) override
 
int replica_comm_send (char *msg_data, int msg_len, int dest_rep) override
 
int check_atom_name_selections_available () override
 
int init_atom (int atom_number) override
 
int check_atom_id (int atom_number) override
 
int init_atom (cvm::residue_id const &residue, std::string const &atom_name, std::string const &segment_id) override
 
int check_atom_id (cvm::residue_id const &residue, std::string const &atom_name, std::string const &segment_id) override
 
void clear_atom (int index) override
 
void update_atom_properties (int index)
 
cvm::rvector position_distance (cvm::atom_pos const &pos1, cvm::atom_pos const &pos2) const
 
int load_atoms_pdb (char const *filename, cvm::atom_group &atoms, std::string const &pdb_field, double const pdb_field_value) override
 
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
 
int scalable_group_coms () override
 
int init_atom_group (std::vector< int > const &atoms_ids) override
 
void clear_atom_group (int index) override
 
int update_group_properties (int index)
 
int check_volmaps_available () override
 
int init_volmap_by_id (int volmap_id) override
 
int init_volmap_by_name (const char *volmap_name) override
 
int check_volmap_by_id (int volmap_id) override
 
int check_volmap_by_name (char const *volmap_name) override
 
int get_volmap_id_from_name (char const *volmap_name) override
 
void clear_volmap (int index) override
 
int compute_volmap (int flags, int volmap_id, cvm::atom_iter atom_begin, cvm::atom_iter atom_end, cvm::real *value, cvm::real *atom_field) override
 
template<class T >
void getGridForceGridValue (int flags, T const *grid, cvm::atom_iter atom_begin, cvm::atom_iter atom_end, cvm::real *value, cvm::real *atom_field)
 Abstraction of the two types of NAMD volumetric maps. More...
 
template<class T , int flags>
void GridForceGridLoop (T const *g, cvm::atom_iter atom_begin, cvm::atom_iter atom_end, cvm::real *value, cvm::real *atom_field)
 Implementation of inner loop; allows for atom list computation and use. More...
 
std::ostream & output_stream (std::string const &output_name, std::string const description) override
 
int flush_output_stream (std::string const &output_name) override
 
int flush_output_streams () override
 
int close_output_stream (std::string const &output_name) override
 
int close_output_streams () override
 
int backup_file (char const *filename) override
 
- Public Member Functions inherited from GlobalMaster
void processData (AtomIDList::iterator a_i, AtomIDList::iterator a_e, PositionList::iterator p_i, PositionList::iterator g_i, PositionList::iterator g_e, BigRealList::iterator gm_i, BigRealList::iterator gm_e, ForceList::iterator gtf_i, ForceList::iterator gtf_e, IntList::iterator goi_i, IntList::iterator goi_e, BigRealList::iterator gov_i, BigRealList::iterator gov_e, AtomIDList::iterator last_atoms_forced_i, AtomIDList::iterator last_atoms_forced_e, ForceList::iterator last_forces_i, AtomIDList::iterator, AtomIDList::iterator, ForceList::iterator)
 
bool changedAtoms ()
 
const AtomIDListrequestedAtoms ()
 
bool changedForces ()
 
const AtomIDListforcedAtoms ()
 
const ForceListappliedForces ()
 
bool changedGroups ()
 
const ResizeArray< AtomIDList > & requestedGroups ()
 
const ForceListgroupForces ()
 
bool changedGridObjs ()
 
const IntListrequestedGridObjs ()
 
const BigRealListgridObjForces ()
 
bool requestedTotalForces ()
 
void clearChanged ()
 
virtual ~GlobalMaster ()
 
void check () const
 
void setLattice (const Lattice *lat)
 

Protected Member Functions

void update_accelMD_info ()
 
- Protected Member Functions inherited from GlobalMaster
 GlobalMaster ()
 
AtomIDListmodifyRequestedAtoms ()
 
AtomIDListmodifyForcedAtoms ()
 
ForceListmodifyAppliedForces ()
 
ResizeArray< AtomIDList > & modifyRequestedGroups ()
 
ForceListmodifyGroupForces ()
 
IntListmodifyRequestedGridObjects ()
 
BigRealListmodifyGridObjForces ()
 
AtomIDList::const_iterator getAtomIdBegin ()
 
AtomIDList::const_iterator getAtomIdEnd ()
 
PositionList::const_iterator getAtomPositionBegin ()
 
PositionList::const_iterator getGroupPositionBegin ()
 
PositionList::const_iterator getGroupPositionEnd ()
 
ForceList::const_iterator getGroupTotalForceBegin ()
 
ForceList::const_iterator getGroupTotalForceEnd ()
 
IntList::const_iterator getGridObjIndexBegin ()
 
IntList::const_iterator getGridObjIndexEnd ()
 
BigRealList::const_iterator getGridObjValueBegin ()
 
BigRealList::const_iterator getGridObjValueEnd ()
 
AtomIDList::const_iterator getLastAtomsForcedBegin ()
 
AtomIDList::const_iterator getLastAtomsForcedEnd ()
 
ForceList::const_iterator getLastForcesBegin ()
 
AtomIDList::const_iterator getForceIdBegin ()
 
AtomIDList::const_iterator getForceIdEnd ()
 
ForceList::const_iterator getTotalForce ()
 
void requestTotalForce (bool yesno=true)
 
BigRealList::const_iterator getGroupMassBegin ()
 
BigRealList::const_iterator getGroupMassEnd ()
 

Protected Attributes

std::vector< int > atoms_map
 Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data. More...
 
SimParameterssimparams
 Pointer to the NAMD simulation input object. More...
 
Random random
 NAMD-style PRNG object. More...
 
bool first_timestep
 
cvm::step_number previous_NAMD_step
 
SubmitReductionreduction
 Used to submit restraint energy as MISC. More...
 
bool accelMDOn
 Accelerated MD reweighting factor. More...
 
cvm::real amd_weight_factor
 
- Protected Attributes inherited from GlobalMaster
bool totalForceRequested
 
const Latticelattice
 
AtomIDList::iterator atomIdBegin
 
AtomIDList::iterator atomIdEnd
 
PositionList::iterator atomPositionBegin
 
PositionList::iterator groupPositionBegin
 
PositionList::iterator groupPositionEnd
 
BigRealList::iterator groupMassBegin
 
BigRealList::iterator groupMassEnd
 
ForceList::iterator groupTotalForceBegin
 
ForceList::iterator groupTotalForceEnd
 
IntList::iterator gridObjIndexBegin
 
IntList::iterator gridObjIndexEnd
 
BigRealList::iterator gridObjValueBegin
 
BigRealList::iterator gridObjValueEnd
 
AtomIDList::iterator lastAtomsForcedBegin
 
ForceList::iterator lastForcesBegin
 
AtomIDList::iterator lastAtomsForcedEnd
 
AtomIDList::iterator forceIdBegin
 
AtomIDList::iterator forceIdEnd
 
ForceList::iterator totalForceBegin
 
bool reqAtomsChanged
 
AtomIDList reqAtoms
 
bool appForcesChanged
 
AtomIDList fAtoms
 
ForceList appForces
 
bool reqGroupsChanged
 
ResizeArray< AtomIDListreqGroups
 
ForceList grpForces
 
bool reqGridObjsChanged
 
IntList reqGridObjs
 
BigRealList gridobjForces
 

Friends

class cvm::atom
 

Additional Inherited Members

- Public Attributes inherited from GlobalMaster
int step
 
int globalMasterStep
 
int old_num_groups_requested
 

Detailed Description

Communication between colvars and NAMD (implementation of colvarproxy)

Definition at line 38 of file colvarproxy_namd.h.

Constructor & Destructor Documentation

◆ colvarproxy_namd()

colvarproxy_namd::colvarproxy_namd ( )

Definition at line 47 of file colvarproxy_namd.C.

References accelMDOn, SimParameters::accelMDOn, amd_weight_factor, COLVARPROXY_VERSION, Node::colvars, Node::configList, StringList::data, debug, SimParameters::dt, endi(), error(), ConfigList::find(), first_timestep, SimParameters::firstTimestep, GLOBAL_MASTER_CKLOOP_CALC_BIASES, GLOBAL_MASTER_CKLOOP_CALC_ITEM, GLOBAL_MASTER_CKLOOP_CALC_SCRIPTED_BIASES, init_atoms_map(), iout, StringList::next, Node::Object(), ReductionMgr::Object(), SimParameters::outputFilename, random, SimParameters::randomSeed, reduction, PatchData::reduction, REDUCTIONS_BASIC, GlobalMaster::requestTotalForce(), SimParameters::restartFilename, SimParameters::restartFrequency, setup(), Node::simParameters, simparams, update_target_temperature(), and ReductionMgr::willSubmit().

48 {
49  engine_name_ = "NAMD";
50 
51  version_int = get_version_from_string(COLVARPROXY_VERSION);
52 #if CMK_TRACE_ENABLED
53  if ( 0 == CkMyPe() ) {
54  traceRegisterUserEvent("GM COLVAR item", GLOBAL_MASTER_CKLOOP_CALC_ITEM);
55  traceRegisterUserEvent("GM COLVAR bias", GLOBAL_MASTER_CKLOOP_CALC_BIASES );
56  traceRegisterUserEvent("GM COLVAR scripted bias", GLOBAL_MASTER_CKLOOP_CALC_SCRIPTED_BIASES );
57  }
58 #endif
59  first_timestep = true;
60  requestTotalForce(total_force_requested);
61 
62  boltzmann_ = 0.001987191;
63 
64  angstrom_value_ = 1.;
65 
66  // initialize pointers to NAMD configuration data
68 
69  if (cvm::debug())
70  iout << "Info: initializing the colvars proxy object.\n" << endi;
71 
72  // find the configuration file, if provided
73  StringList *config = Node::Object()->configList->find("colvarsConfig");
74 
75  // find the input state file
76  StringList *input_restart = Node::Object()->configList->find("colvarsInput");
77  colvarproxy_io::set_input_prefix(input_restart ? input_restart->data : "");
78 
80  set_integration_timestep(simparams->dt);
81 
83 
84  // both fields are taken from data structures already available
85  updated_masses_ = updated_charges_ = true;
86 
87  // Take the output prefixes from the NAMD input
88  colvarproxy_io::set_output_prefix(std::string(simparams->outputFilename));
89  colvarproxy_io::set_restart_output_prefix(std::string(simparams->restartFilename));
90  colvarproxy_io::set_default_restart_frequency(simparams->restartFrequency);
91 
92  if (simparams->accelMDOn) {
93  accelMDOn = true;
94  } else {
95  accelMDOn = false;
96  }
97  amd_weight_factor = 1.0;
98 
99  // check if it is possible to save output configuration
100  if ((!output_prefix_str.size()) && (!restart_output_prefix_str.size())) {
101  error("Error: neither the final output state file or "
102  "the output restart file could be defined, exiting.\n");
103  }
104 
105  init_atoms_map();
106 
107  // initialize module: this object will be the communication proxy
108  colvars = new colvarmodule(this);
109 
110  cvm::log("Using NAMD interface, version "+
111  cvm::to_str(COLVARPROXY_VERSION)+".\n");
112  colvars->cite_feature("NAMD engine");
113  colvars->cite_feature("Colvars-NAMD interface");
114 
115  errno = 0;
116  for ( ; config; config = config->next ) {
117  add_config("configfile", config->data);
118  }
119 
120  // Trigger immediate initialization of the module
121  colvarproxy::parse_module_config();
123  colvars->update_engine_parameters();
124  colvars->setup_input();
125  colvars->setup_output();
126 
127  // save to Node for Tcl script access
128  Node::Object()->colvars = colvars;
129 
130 #ifdef NAMD_TCL
131  have_scripts = true;
132 #else
133  have_scripts = false;
134 #endif
135 
136  if (simparams->firstTimestep != 0) {
137  colvars->set_initial_step(static_cast<cvm::step_number>(simparams->firstTimestep));
138  }
139 
141 
142  #ifdef NODEGROUP_FORCE_REGISTER
143  CProxy_PatchData cpdata(CkpvAccess(BOCclass_group).patchData);
144  PatchData *patchData = cpdata.ckLocalBranch();
145  nodeReduction = patchData->reduction;
146  #endif
147 
148  if (cvm::debug())
149  iout << "Info: done initializing the colvars proxy object.\n" << endi;
150 }
static Node * Object()
Definition: Node.h:86
SimParameters * simparams
Pointer to the NAMD simulation input object.
Random random
NAMD-style PRNG object.
#define GLOBAL_MASTER_CKLOOP_CALC_ITEM
void error(std::string const &message) override
SimParameters * simParameters
Definition: Node.h:181
cvm::real amd_weight_factor
std::ostream & endi(std::ostream &s)
Definition: InfoStream.C:54
SubmitReduction * willSubmit(int setID, int size=-1)
Definition: ReductionMgr.C:366
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:279
#define iout
Definition: InfoStream.h:51
NodeReduction * reduction
Definition: PatchData.h:133
char outputFilename[NAMD_FILENAME_BUFFER_SIZE]
int setup() override
#define COLVARPROXY_VERSION
#define GLOBAL_MASTER_CKLOOP_CALC_SCRIPTED_BIASES
Definition: Random.h:37
bool accelMDOn
Accelerated MD reweighting factor.
void init_atoms_map()
Allocate an atoms map with the same size as the NAMD topology.
#define GLOBAL_MASTER_CKLOOP_CALC_BIASES
ConfigList * configList
Definition: Node.h:182
char restartFilename[NAMD_FILENAME_BUFFER_SIZE]
unsigned int randomSeed
void requestTotalForce(bool yesno=true)
Definition: GlobalMaster.h:134
StringList * next
Definition: ConfigList.h:49
char * data
Definition: ConfigList.h:48
colvarmodule * colvars
Definition: Node.h:187
SubmitReduction * reduction
Used to submit restraint energy as MISC.
StringList * find(const char *name) const
Definition: ConfigList.C:341

◆ ~colvarproxy_namd()

colvarproxy_namd::~colvarproxy_namd ( )

Definition at line 153 of file colvarproxy_namd.C.

References reduction.

154 {
155  delete reduction;
156 }
SubmitReduction * reduction
Used to submit restraint energy as MISC.

Member Function Documentation

◆ accelMD_enabled()

bool colvarproxy_namd::accelMD_enabled ( ) const
inlineoverride

Definition at line 118 of file colvarproxy_namd.h.

References accelMDOn.

118  {
119  return accelMDOn;
120  }
bool accelMDOn
Accelerated MD reweighting factor.

◆ add_energy()

void colvarproxy_namd::add_energy ( cvm::real  energy)
override

Definition at line 654 of file colvarproxy_namd.C.

References SimParameters::CUDASOAintegrate, SubmitReduction::item(), reduction, REDUCTION_MISC_ENERGY, and simparams.

655 {
656  #ifdef NODEGROUP_FORCE_REGISTER
658  nodeReduction->item(REDUCTION_MISC_ENERGY) += energy;
659  } else {
661  }
662  #else
664  #endif
665 }
SimParameters * simparams
Pointer to the NAMD simulation input object.
BigReal & item(int i)
Definition: ReductionMgr.h:313
SubmitReduction * reduction
Used to submit restraint energy as MISC.

◆ backup_file()

int colvarproxy_namd::backup_file ( char const *  filename)
override

Definition at line 1175 of file colvarproxy_namd.C.

References NAMD_backup_file().

Referenced by output_stream().

1176 {
1177  if (std::string(filename).rfind(std::string(".colvars.state")) != std::string::npos) {
1178  NAMD_backup_file(filename, ".old");
1179  } else {
1180  NAMD_backup_file(filename, ".BAK");
1181  }
1182  return COLVARS_OK;
1183 }
void NAMD_backup_file(const char *filename, const char *extension)
Definition: common.C:235

◆ calculate()

void colvarproxy_namd::calculate ( )
virtual

Reimplemented from GlobalMaster.

Definition at line 336 of file colvarproxy_namd.C.

References Lattice::a(), Lattice::a_p(), ResizeArray< Elem >::add(), atoms_map, Lattice::b(), Lattice::b_p(), ResizeArray< Elem >::begin(), Lattice::c(), Lattice::c_p(), ResizeArray< Elem >::clear(), SimParameters::CUDASOAintegrate, debug, ResizeArray< Elem >::end(), first_timestep, GlobalMaster::getAtomIdBegin(), GlobalMaster::getAtomIdEnd(), GlobalMaster::getAtomPositionBegin(), GlobalMaster::getForceIdBegin(), GlobalMaster::getForceIdEnd(), GlobalMaster::getGridObjIndexBegin(), GlobalMaster::getGridObjIndexEnd(), GlobalMaster::getGridObjValueBegin(), GlobalMaster::getGridObjValueEnd(), GlobalMaster::getGroupPositionBegin(), GlobalMaster::getGroupPositionEnd(), GlobalMaster::getGroupTotalForceBegin(), GlobalMaster::getGroupTotalForceEnd(), GlobalMaster::getTotalForce(), GlobalMaster::lattice, log(), GlobalMaster::modifyAppliedForces(), GlobalMaster::modifyForcedAtoms(), GlobalMaster::modifyGridObjForces(), GlobalMaster::modifyGroupForces(), Node::molecule, SimParameters::N, Molecule::numAtoms, Node::Object(), Lattice::orthogonal(), SimParameters::outputFilename, previous_NAMD_step, reduction, GlobalMaster::requestedGridObjs(), GlobalMaster::requestedGroups(), ResizeArray< Elem >::resize(), SimParameters::restartFilename, SimParameters::restartFrequency, Vector::set(), ResizeArray< Elem >::setall(), setup(), simparams, GlobalMaster::step, SubmitReduction::submit(), update_accelMD_info(), update_atoms_map(), Vector::x, Vector::y, and Vector::z.

337 {
338  errno = 0;
339 
340  if (first_timestep) {
341 
342  // First run after the proxy is constructed
343 
345  colvars->update_engine_parameters();
346  colvars->setup_input();
347  colvars->setup_output();
348 
349  first_timestep = false;
350 
351  } else {
352 
353  // Use the time step number inherited from GlobalMaster
354  if ( step - previous_NAMD_step == 1 ) {
355  colvars->it++;
356  b_simulation_continuing = false;
357  } else {
358 
359  // Cases covered by this condition:
360  // - run 0
361  // - beginning of a new run statement
362  // The internal counter is not incremented, and the objects are made
363  // aware of this via the following flag
364  b_simulation_continuing = true;
365 
366  // Update NAMD output and restart prefixes
367  colvarproxy_io::set_output_prefix(std::string(simparams->outputFilename));
368  colvarproxy_io::set_restart_output_prefix(std::string(simparams->restartFilename));
369  colvarproxy_io::set_default_restart_frequency(simparams->restartFrequency);
370  colvars->setup_output();
371 
372  }
373  }
374 
377 
378  {
379  Vector const a = lattice->a();
380  Vector const b = lattice->b();
381  Vector const c = lattice->c();
382  unit_cell_x.set(a.x, a.y, a.z);
383  unit_cell_y.set(b.x, b.y, c.z);
384  unit_cell_z.set(c.x, c.y, c.z);
385  }
386 
387  if (!lattice->a_p() && !lattice->b_p() && !lattice->c_p()) {
388  boundaries_type = boundaries_non_periodic;
389  reset_pbc_lattice();
390  } else if (lattice->a_p() && lattice->b_p() && lattice->c_p()) {
391  if (lattice->orthogonal()) {
392  boundaries_type = boundaries_pbc_ortho;
393  } else {
394  boundaries_type = boundaries_pbc_triclinic;
395  }
396  colvarproxy_system::update_pbc_lattice();
397  } else {
398  boundaries_type = boundaries_unsupported;
399  }
400 
401  if (cvm::debug()) {
402  cvm::log(std::string(cvm::line_marker)+
403  "colvarproxy_namd, step no. "+cvm::to_str(colvars->it)+"\n"+
404  "Updating atomic data arrays.\n");
405  }
406 
407  // must delete the forces applied at the previous step: we can do
408  // that because they have already been used and copied to other
409  // memory locations
412 
413  // If new atomic positions or forces have been requested by other
414  // GlobalMaster objects, add these to the atom map as well
415  size_t const n_all_atoms = Node::Object()->molecule->numAtoms;
416  if ( (atoms_map.size() != n_all_atoms) ||
417  (int(atoms_ids.size()) < (getAtomIdEnd() - getAtomIdBegin())) ||
418  (int(atoms_ids.size()) < (getForceIdEnd() - getForceIdBegin())) ) {
421  }
422 
423  // prepare local arrays
424  for (size_t i = 0; i < atoms_ids.size(); i++) {
425  atoms_positions[i] = cvm::rvector(0.0, 0.0, 0.0);
426  atoms_total_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
427  atoms_new_colvar_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
428  }
429 
430  for (size_t i = 0; i < atom_groups_ids.size(); i++) {
431  atom_groups_total_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
432  atom_groups_new_colvar_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
433  }
434 
435 #if NAMD_VERSION_NUMBER >= 34471681
436  for (int imap = 0; imap < volmaps_ids.size(); imap++) {
437  volmaps_new_colvar_forces[imap] = 0.0;
438  }
439 #endif
440 
441  {
442  if (cvm::debug()) {
443  cvm::log("Updating positions arrays.\n");
444  }
445  size_t n_positions = 0;
449 
450  for ( ; a_i != a_e; ++a_i, ++p_i ) {
451  atoms_positions[atoms_map[*a_i]] = cvm::rvector((*p_i).x, (*p_i).y, (*p_i).z);
452  n_positions++;
453  }
454 
455  // The following had to be relaxed because some atoms may be forced without their position being requested
456  // if (n_positions < atoms_ids.size()) {
457  // cvm::error("Error: did not receive the positions of all atoms.\n", COLVARS_BUG_ERROR);
458  // }
459  }
460 
461  if (total_force_requested && cvm::step_relative() > 0) {
462 
463  // sort the force arrays the previous step
464  // (but only do so if there *is* a previous step!)
465 
466  {
467  if (cvm::debug()) {
468  cvm::log("Updating total forces arrays.\n");
469  }
470  size_t n_total_forces = 0;
474 
475  for ( ; a_i != a_e; ++a_i, ++f_i ) {
476  atoms_total_forces[atoms_map[*a_i]] = cvm::rvector((*f_i).x, (*f_i).y, (*f_i).z);
477  n_total_forces++;
478  }
479 
480  if ( (! b_simulation_continuing) &&
481  (n_total_forces < atoms_ids.size()) ) {
482  cvm::error("Error: total forces were requested, but total forces "
483  "were not received for all atoms.\n"
484  "The most probable cause is combination of energy "
485  "minimization with a biasing method that requires MD (e.g. ABF).\n"
486  "Always run minimization and ABF separately.", COLVARS_INPUT_ERROR);
487  }
488  }
489 
490  {
491  if (cvm::debug()) {
492  cvm::log("Updating group total forces arrays.\n");
493  }
496  size_t i = 0;
497  if ( (! b_simulation_continuing) &&
498  ((f_e - f_i) != ((int) atom_groups_ids.size())) ) {
499  cvm::error("Error: total forces were requested for scalable groups, "
500  "but they are not in the same number from the number of groups.\n"
501  "The most probable cause is combination of energy "
502  "minimization with a biasing method that requires MD (e.g. ABF).\n"
503  "Always run minimization and ABF separately.", COLVARS_INPUT_ERROR);
504  }
505  for ( ; f_i != f_e; f_i++, i++) {
506  atom_groups_total_forces[i] = cvm::rvector((*f_i).x, (*f_i).y, (*f_i).z);
507  }
508  }
509  }
510 
511  {
512  if (cvm::debug()) {
513  cvm::log("Updating group positions arrays.\n");
514  }
515  // update group data (only coms available so far)
516  size_t ig;
517  // note: getGroupMassBegin() could be used here, but masses and charges
518  // have already been calculated from the last call to setup()
520  for (ig = 0; gp_i != getGroupPositionEnd(); gp_i++, ig++) {
521  atom_groups_coms[ig] = cvm::rvector(gp_i->x, gp_i->y, gp_i->z);
522  }
523  }
524 
525 #if NAMD_VERSION_NUMBER >= 34471681
526  {
527  if (cvm::debug()) {
528  log("Updating grid objects.\n");
529  }
530  // Using a simple nested loop: there probably won't be so many maps that
531  // this becomes performance-limiting
534  for ( ; gov_i != getGridObjValueEnd(); goi_i++, gov_i++) {
535  for (size_t imap = 0; imap < volmaps_ids.size(); imap++) {
536  if (volmaps_ids[imap] == *goi_i) {
537  volmaps_values[imap] = *gov_i;
538  break;
539  }
540  }
541  }
542  }
543 #endif
544 
545  if (cvm::debug()) {
546  print_input_atomic_data();
547  }
548 
549  // call the collective variable module
550  if (colvars->calc() != COLVARS_OK) {
551  cvm::error("Error in the collective variables module.\n", COLVARS_ERROR);
552  }
553 
554  if (cvm::debug()) {
555  print_output_atomic_data();
556  }
557 
558  // communicate all forces to the MD integrator
559  for (size_t i = 0; i < atoms_ids.size(); i++) {
560  cvm::rvector const &f = atoms_new_colvar_forces[i];
561  modifyForcedAtoms().add(atoms_ids[i]);
562  modifyAppliedForces().add(Vector(f.x, f.y, f.z));
563  }
564 
565  if (atom_groups_new_colvar_forces.size() > 0) {
568  for (int ig = 0; gf_i != modifyGroupForces().end(); gf_i++, ig++) {
569  cvm::rvector const &f = atom_groups_new_colvar_forces[ig];
570  *gf_i = Vector(f.x, f.y, f.z);
571  }
572  }
573 
574 #if NAMD_VERSION_NUMBER >= 34471681
575  if (volmaps_new_colvar_forces.size() > 0) {
580  for ( ; goi_i != getGridObjIndexEnd(); goi_i++, gof_i++) {
581  for (size_t imap = 0; imap < volmaps_ids.size(); imap++) {
582  if (volmaps_ids[imap] == *goi_i) {
583  *gof_i = volmaps_new_colvar_forces[imap];
584  break;
585  }
586  }
587  }
588  }
589 #endif
590 
591  // send MISC energy
592  #ifdef NODEGROUP_FORCE_REGISTER
594  reduction->submit();
595  }
596  #else
597  reduction->submit();
598  #endif
599 
600  // NAMD does not destruct GlobalMaster objects, so we must remember
601  // to write all output files at the end of a run
602  if (step == simparams->N) {
603  post_run();
604  }
605 }
static Node * Object()
Definition: Node.h:86
ForceList & modifyAppliedForces()
Definition: GlobalMaster.C:163
BigRealList & modifyGridObjForces()
Definition: GlobalMaster.C:180
SimParameters * simparams
Pointer to the NAMD simulation input object.
NAMD_HOST_DEVICE Vector c() const
Definition: Lattice.h:270
NAMD_HOST_DEVICE int c_p() const
Definition: Lattice.h:291
AtomIDList::const_iterator getForceIdEnd()
Definition: GlobalMaster.C:261
AtomIDList::const_iterator getAtomIdBegin()
Definition: GlobalMaster.C:191
PositionList::const_iterator getGroupPositionEnd()
Definition: GlobalMaster.C:207
Definition: Vector.h:72
BigReal z
Definition: Vector.h:74
ForceList::const_iterator getGroupTotalForceBegin()
Definition: GlobalMaster.C:211
cvm::step_number previous_NAMD_step
static int debug
Definition: parm.C:36
NAMD_HOST_DEVICE int orthogonal() const
Definition: Lattice.h:273
void clear()
Definition: ResizeArray.h:91
int add(const Elem &elem)
Definition: ResizeArray.h:101
PositionList::const_iterator getAtomPositionBegin()
Definition: GlobalMaster.C:199
NAMD_HOST_DEVICE int b_p() const
Definition: Lattice.h:290
const ResizeArray< AtomIDList > & requestedGroups()
Definition: GlobalMaster.C:150
char outputFilename[NAMD_FILENAME_BUFFER_SIZE]
int setup() override
AtomIDList::const_iterator getForceIdBegin()
Definition: GlobalMaster.C:256
void resize(int i)
Definition: ResizeArray.h:84
void setall(const Elem &elem)
Definition: ResizeArray.h:94
std::vector< int > atoms_map
Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data...
void log(std::string const &message) override
IntList::const_iterator getGridObjIndexBegin()
Definition: GlobalMaster.C:219
const Elem * const_iterator
Definition: ResizeArray.h:38
BigReal x
Definition: Vector.h:74
NAMD_HOST_DEVICE int a_p() const
Definition: Lattice.h:289
AtomIDList & modifyForcedAtoms()
Definition: GlobalMaster.C:158
int numAtoms
Definition: Molecule.h:585
ForceList & modifyGroupForces()
Definition: GlobalMaster.C:168
NAMD_HOST_DEVICE Vector b() const
Definition: Lattice.h:269
char restartFilename[NAMD_FILENAME_BUFFER_SIZE]
Bool set(const char *s)
Definition: Vector.h:261
BigRealList::const_iterator getGridObjValueBegin()
Definition: GlobalMaster.C:227
IntList::const_iterator getGridObjIndexEnd()
Definition: GlobalMaster.C:223
const Lattice * lattice
Definition: GlobalMaster.h:142
iterator begin(void)
Definition: ResizeArray.h:36
AtomIDList::const_iterator getAtomIdEnd()
Definition: GlobalMaster.C:195
iterator end(void)
Definition: ResizeArray.h:37
BigReal y
Definition: Vector.h:74
void submit(void)
Definition: ReductionMgr.h:324
const IntList & requestedGridObjs()
Definition: GlobalMaster.C:154
SubmitReduction * reduction
Used to submit restraint energy as MISC.
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:203
Molecule * molecule
Definition: Node.h:179
ForceList::const_iterator getTotalForce()
Definition: GlobalMaster.C:266
BigRealList::const_iterator getGridObjValueEnd()
Definition: GlobalMaster.C:231
ForceList::const_iterator getGroupTotalForceEnd()
Definition: GlobalMaster.C:215

◆ check_atom_id() [1/2]

int colvarproxy_namd::check_atom_id ( int  atom_number)
override

Definition at line 711 of file colvarproxy_namd.C.

References debug, Node::molecule, Molecule::numAtoms, and Node::Object().

Referenced by init_atom().

712 {
713  // NAMD's internal numbering starts from zero
714  int const aid = (atom_number-1);
715 
716  if (cvm::debug())
717  cvm::log("Adding atom "+cvm::to_str(atom_number)+
718  " for collective variables calculation.\n");
719 
720  if ( (aid < 0) || (aid >= Node::Object()->molecule->numAtoms) ) {
721  cvm::error("Error: invalid atom number specified, "+
722  cvm::to_str(atom_number)+"\n", COLVARS_INPUT_ERROR);
723  return COLVARS_INPUT_ERROR;
724  }
725 
726  return aid;
727 }
static Node * Object()
Definition: Node.h:86
static int debug
Definition: parm.C:36
int numAtoms
Definition: Molecule.h:585
Molecule * molecule
Definition: Node.h:179

◆ check_atom_id() [2/2]

int colvarproxy_namd::check_atom_id ( cvm::residue_id const &  residue,
std::string const &  atom_name,
std::string const &  segment_id 
)
override

Definition at line 764 of file colvarproxy_namd.C.

References Molecule::get_atom_from_name(), Node::molecule, and Node::Object().

767 {
768  int const aid =
769  (segment_id.size() ?
770  Node::Object()->molecule->get_atom_from_name(segment_id.c_str(),
771  residue,
772  atom_name.c_str()) :
774  residue,
775  atom_name.c_str()));
776 
777  if (aid < 0) {
778  // get_atom_from_name() has returned an error value
779  cvm::error("Error: could not find atom \""+
780  atom_name+"\" in residue "+
781  cvm::to_str(residue)+
782  ( (segment_id != "MAIN") ?
783  (", segment \""+segment_id+"\"") :
784  ("") )+
785  "\n", COLVARS_INPUT_ERROR);
786  return COLVARS_INPUT_ERROR;
787  }
788 
789  return aid;
790 }
static Node * Object()
Definition: Node.h:86
int get_atom_from_name(const char *segid, int resid, const char *aname) const
Definition: Molecule.C:126
Molecule * molecule
Definition: Node.h:179

◆ check_atom_name_selections_available()

int colvarproxy_namd::check_atom_name_selections_available ( )
override

Definition at line 730 of file colvarproxy_namd.C.

731 {
732  return COLVARS_OK;
733 }

◆ check_replicas_enabled()

int colvarproxy_namd::check_replicas_enabled ( )
override

Definition at line 1609 of file colvarproxy_namd.C.

1609  {
1610 #if CMK_HAS_PARTITION
1611  return COLVARS_OK;
1612 #else
1613  return COLVARS_NOT_IMPLEMENTED;
1614 #endif
1615 }

◆ check_volmap_by_id()

int colvarproxy_namd::check_volmap_by_id ( int  volmap_id)
override

Definition at line 1374 of file colvarproxy_namd.C.

References Node::molecule, Molecule::numGridforceGrids, and Node::Object().

Referenced by init_volmap_by_id().

1375 {
1376  Molecule *mol = Node::Object()->molecule;
1377  if ((volmap_id < 0) || (volmap_id >= mol->numGridforceGrids)) {
1378  return cvm::error("Error: invalid numeric ID ("+cvm::to_str(volmap_id)+
1379  ") for map.\n", COLVARS_INPUT_ERROR);
1380  }
1381  colvars->cite_feature("GridForces volumetric map implementation for NAMD");
1382  return COLVARS_OK;
1383 }
static Node * Object()
Definition: Node.h:86
int numGridforceGrids
Definition: Molecule.h:624
Molecule stores the structural information for the system.
Definition: Molecule.h:175
Molecule * molecule
Definition: Node.h:179

◆ check_volmap_by_name()

int colvarproxy_namd::check_volmap_by_name ( char const *  volmap_name)
override

Definition at line 1386 of file colvarproxy_namd.C.

References MGridforceParamsList::index_for_key(), SimParameters::mgridforcelist, and simparams.

Referenced by get_volmap_id_from_name(), and init_volmap_by_name().

1387 {
1388  if (volmap_name == NULL) {
1389  return cvm::error("Error: no grid object name provided.", COLVARS_INPUT_ERROR);
1390  }
1391  int volmap_id = simparams->mgridforcelist.index_for_key(volmap_name);
1392  if (volmap_id < 0) {
1393  return cvm::error("Error: invalid map name \""+std::string(volmap_name)+
1394  "\".\n", COLVARS_INPUT_ERROR);
1395  }
1396  colvars->cite_feature("GridForces volumetric map implementation for NAMD");
1397  return COLVARS_OK;
1398 }
SimParameters * simparams
Pointer to the NAMD simulation input object.
int index_for_key(const char *key)
MGridforceParamsList mgridforcelist

◆ check_volmaps_available()

int colvarproxy_namd::check_volmaps_available ( )
override

Definition at line 1305 of file colvarproxy_namd.C.

1306 {
1307  return COLVARS_OK;
1308 }

◆ clear_atom()

void colvarproxy_namd::clear_atom ( int  index)
override

Definition at line 826 of file colvarproxy_namd.C.

827 {
828  colvarproxy::clear_atom(index);
829  // TODO remove it from GlobalMaster arrays?
830 }

◆ clear_atom_group()

void colvarproxy_namd::clear_atom_group ( int  index)
override

Definition at line 1258 of file colvarproxy_namd.C.

1259 {
1260  // do nothing, keep the NAMD arrays in sync with the colvarproxy ones
1261  colvarproxy::clear_atom_group(index);
1262 }

◆ clear_volmap()

void colvarproxy_namd::clear_volmap ( int  index)
override

Definition at line 1401 of file colvarproxy_namd.C.

1402 {
1403  // TODO remove from GlobalMaster
1404  colvarproxy::clear_volmap(index);
1405 }

◆ close_output_stream()

int colvarproxy_namd::close_output_stream ( std::string const &  output_name)
override

Definition at line 1141 of file colvarproxy_namd.C.

1142 {
1143  if (!io_available()) {
1144  return cvm::error("Error: trying to access an output file/channel "
1145  "from the wrong thread.\n", COLVARS_BUG_ERROR);
1146  }
1147 
1148  if (output_streams_.count(output_name) > 0) {
1149  (reinterpret_cast<ofstream_namd *>(output_streams_[output_name]))->close();
1150  delete output_streams_[output_name];
1151  output_streams_.erase(output_name);
1152  }
1153 
1154  return COLVARS_OK;
1155 }

◆ close_output_streams()

int colvarproxy_namd::close_output_streams ( )
override

Definition at line 1158 of file colvarproxy_namd.C.

1159 {
1160  if (! io_available()) {
1161  return COLVARS_OK;
1162  }
1163 
1164  for (std::map<std::string, std::ostream *>::iterator osi = output_streams_.begin();
1165  osi != output_streams_.end();
1166  osi++) {
1167  (reinterpret_cast<ofstream_namd *>(osi->second))->close();
1168  }
1169  output_streams_.clear();
1170 
1171  return COLVARS_OK;
1172 }

◆ compute_volmap()

int colvarproxy_namd::compute_volmap ( int  flags,
int  volmap_id,
cvm::atom_iter  atom_begin,
cvm::atom_iter  atom_end,
cvm::real *  value,
cvm::real *  atom_field 
)
override

Definition at line 1473 of file colvarproxy_namd.C.

References GridforceGrid::get_grid_type(), Molecule::get_gridfrc_grid(), GridforceGrid::GridforceGridTypeFull, GridforceGrid::GridforceGridTypeLite, Node::molecule, and Node::Object().

1479 {
1480  Molecule *mol = Node::Object()->molecule;
1481  GridforceGrid *grid = mol->get_gridfrc_grid(volmap_id);
1482  // Inheritance is not possible with GridForceGrid's design
1484  GridforceFullMainGrid *g = dynamic_cast<GridforceFullMainGrid *>(grid);
1485  getGridForceGridValue<GridforceFullMainGrid>(flags, g, atom_begin, atom_end,
1486  value, atom_field);
1487  } else if (grid->get_grid_type() == GridforceGrid::GridforceGridTypeLite) {
1488  GridforceLiteGrid *g = dynamic_cast<GridforceLiteGrid *>(grid);
1489  getGridForceGridValue<GridforceLiteGrid>(flags, g, atom_begin, atom_end,
1490  value, atom_field);
1491  }
1492  return COLVARS_OK;
1493 }
static Node * Object()
Definition: Node.h:86
GridforceGrid * get_gridfrc_grid(int gridnum) const
Definition: Molecule.h:1363
Molecule stores the structural information for the system.
Definition: Molecule.h:175
GridforceGridType get_grid_type(void)
Definition: GridForceGrid.h:64
Molecule * molecule
Definition: Node.h:179

◆ error()

void colvarproxy_namd::error ( std::string const &  message)
override

Definition at line 690 of file colvarproxy_namd.C.

References log(), NAMD_die(), and NAMD_err().

Referenced by colvarproxy_namd().

691 {
692  log(message);
693  switch (cvm::get_error()) {
694  case COLVARS_FILE_ERROR:
695  errno = EIO; break;
696  case COLVARS_NOT_IMPLEMENTED:
697  errno = ENOSYS; break;
698  case COLVARS_MEMORY_ERROR:
699  errno = ENOMEM; break;
700  }
701  char const *msg = "Error in the collective variables module "
702  "(see above for details)";
703  if (errno) {
704  NAMD_err(msg);
705  } else {
706  NAMD_die(msg);
707  }
708 }
void NAMD_err(const char *err_msg)
Definition: common.C:170
void log(std::string const &message) override
void NAMD_die(const char *err_msg)
Definition: common.C:147

◆ flush_output_stream()

int colvarproxy_namd::flush_output_stream ( std::string const &  output_name)
override

Definition at line 1109 of file colvarproxy_namd.C.

1110 {
1111  if (!io_available()) {
1112  return COLVARS_OK;
1113  }
1114 
1115  if (output_streams_.count(output_name) > 0) {
1116  (reinterpret_cast<ofstream_namd *>(output_streams_[output_name]))->flush();
1117  return COLVARS_OK;
1118  }
1119 
1120  return cvm::error("Error: trying to flush an output file/channel "
1121  "that wasn't open.\n", COLVARS_BUG_ERROR);
1122 }

◆ flush_output_streams()

int colvarproxy_namd::flush_output_streams ( )
override

Definition at line 1125 of file colvarproxy_namd.C.

1126 {
1127  if (!io_available()) {
1128  return COLVARS_OK;
1129  }
1130 
1131  for (std::map<std::string, std::ostream *>::iterator osi = output_streams_.begin();
1132  osi != output_streams_.end();
1133  osi++) {
1134  (reinterpret_cast<ofstream_namd *>(osi->second))->flush();
1135  }
1136 
1137  return COLVARS_OK;
1138 }

◆ get_accelMD_factor()

cvm::real colvarproxy_namd::get_accelMD_factor ( ) const
inlineoverride

Definition at line 114 of file colvarproxy_namd.h.

References amd_weight_factor.

114  {
115  return amd_weight_factor;
116  }
cvm::real amd_weight_factor

◆ get_volmap_id_from_name()

int colvarproxy_namd::get_volmap_id_from_name ( char const *  volmap_name)
override

Definition at line 1408 of file colvarproxy_namd.C.

References check_volmap_by_name(), MGridforceParamsList::index_for_key(), SimParameters::mgridforcelist, and simparams.

1409 {
1410  int const volmap_id =
1411  simparams->mgridforcelist.index_for_key(volmap_name);
1412  if (volmap_id < 0) {
1413  // Print error
1414  check_volmap_by_name(volmap_name);
1415  }
1416  return volmap_id;
1417 }
SimParameters * simparams
Pointer to the NAMD simulation input object.
int check_volmap_by_name(char const *volmap_name) override
int index_for_key(const char *key)
MGridforceParamsList mgridforcelist

◆ getGridForceGridValue()

template<class T >
void colvarproxy_namd::getGridForceGridValue ( int  flags,
T const *  grid,
cvm::atom_iter  atom_begin,
cvm::atom_iter  atom_end,
cvm::real *  value,
cvm::real *  atom_field 
)

Abstraction of the two types of NAMD volumetric maps.

Definition at line 1454 of file colvarproxy_namd.C.

1460 {
1461  if (flags & volmap_flag_use_atom_field) {
1462  int const new_flags = volmap_flag_use_atom_field | volmap_flag_gradients;
1463  GridForceGridLoop<T, new_flags>(g, atom_begin, atom_end,
1464  value, atom_field);
1465  } else {
1466  int const new_flags = volmap_flag_gradients;
1467  GridForceGridLoop<T, new_flags>(g, atom_begin, atom_end,
1468  value, atom_field);
1469  }
1470 }

◆ GridForceGridLoop()

template<class T , int flags>
void colvarproxy_namd::GridForceGridLoop ( T const *  g,
cvm::atom_iter  atom_begin,
cvm::atom_iter  atom_end,
cvm::real *  value,
cvm::real *  atom_field 
)

Implementation of inner loop; allows for atom list computation and use.

Definition at line 1421 of file colvarproxy_namd.C.

References Vector::x, Vector::y, and Vector::z.

1426 {
1427  float V = 0.0f;
1428  Vector dV(0.0);
1429  int i = 0;
1430  cvm::atom_iter ai = atom_begin;
1431  for ( ; ai != atom_end; ai++, i++) {
1432  if (g->compute_VdV(Position(ai->pos.x, ai->pos.y, ai->pos.z), V, dV)) {
1433  // out-of-bounds atom
1434  V = 0.0f;
1435  dV = 0.0;
1436  } else {
1437  if (flags & volmap_flag_use_atom_field) {
1438  *value += V * atom_field[i];
1439  if (flags & volmap_flag_gradients) {
1440  ai->grad += atom_field[i] * cvm::rvector(dV.x, dV.y, dV.z);
1441  }
1442  } else {
1443  *value += V;
1444  if (flags & volmap_flag_gradients) {
1445  ai->grad += cvm::rvector(dV.x, dV.y, dV.z);
1446  }
1447  }
1448  }
1449  }
1450 }
Definition: Vector.h:72

◆ init_atom() [1/2]

int colvarproxy_namd::init_atom ( int  atom_number)
override

Definition at line 736 of file colvarproxy_namd.C.

References ResizeArray< Elem >::add(), atoms_map, check_atom_id(), GlobalMaster::modifyRequestedAtoms(), and update_atom_properties().

737 {
738  // save time by checking first whether this atom has been requested before
739  // (this is more common than a non-valid atom number)
740  int aid = (atom_number-1);
741 
742  for (size_t i = 0; i < atoms_ids.size(); i++) {
743  if (atoms_ids[i] == aid) {
744  // this atom id was already recorded
745  atoms_refcount[i] += 1;
746  return i;
747  }
748  }
749 
750  aid = check_atom_id(atom_number);
751 
752  if (aid < 0) {
753  return COLVARS_INPUT_ERROR;
754  }
755 
756  int const index = add_atom_slot(aid);
757  atoms_map[aid] = index;
758  modifyRequestedAtoms().add(aid);
759  update_atom_properties(index);
760  return index;
761 }
AtomIDList & modifyRequestedAtoms()
Definition: GlobalMaster.C:128
int add(const Elem &elem)
Definition: ResizeArray.h:101
int check_atom_id(int atom_number) override
std::vector< int > atoms_map
Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data...
void update_atom_properties(int index)

◆ init_atom() [2/2]

int colvarproxy_namd::init_atom ( cvm::residue_id const &  residue,
std::string const &  atom_name,
std::string const &  segment_id 
)
override

For AMBER topologies, the segment id is automatically set to "MAIN" (the segment id assigned by NAMD's AMBER topology parser), and is therefore optional when an AMBER topology is used

Definition at line 797 of file colvarproxy_namd.C.

References ResizeArray< Elem >::add(), atoms_map, check_atom_id(), debug, GlobalMaster::modifyRequestedAtoms(), and update_atom_properties().

800 {
801  int const aid = check_atom_id(residue, atom_name, segment_id);
802 
803  for (size_t i = 0; i < atoms_ids.size(); i++) {
804  if (atoms_ids[i] == aid) {
805  // this atom id was already recorded
806  atoms_refcount[i] += 1;
807  return i;
808  }
809  }
810 
811  if (cvm::debug())
812  cvm::log("Adding atom \""+
813  atom_name+"\" in residue "+
814  cvm::to_str(residue)+
815  " (index "+cvm::to_str(aid)+
816  ") for collective variables calculation.\n");
817 
818  int const index = add_atom_slot(aid);
819  atoms_map[aid] = index;
820  modifyRequestedAtoms().add(aid);
821  update_atom_properties(index);
822  return index;
823 }
AtomIDList & modifyRequestedAtoms()
Definition: GlobalMaster.C:128
static int debug
Definition: parm.C:36
int add(const Elem &elem)
Definition: ResizeArray.h:101
int check_atom_id(int atom_number) override
std::vector< int > atoms_map
Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data...
void update_atom_properties(int index)

◆ init_atom_group()

int colvarproxy_namd::init_atom_group ( std::vector< int > const &  atoms_ids)
override

Definition at line 1186 of file colvarproxy_namd.C.

References debug, GlobalMaster::modifyGroupForces(), GlobalMaster::modifyRequestedGroups(), Node::molecule, Molecule::numAtoms, Node::Object(), ResizeArray< Elem >::resize(), ResizeArray< Elem >::size(), and update_group_properties().

1187 {
1188  if (cvm::debug())
1189  cvm::log("Requesting from NAMD a group of size "+cvm::to_str(atoms_ids.size())+
1190  " for collective variables calculation.\n");
1191 
1192  colvars->cite_feature("Scalable center-of-mass computation (NAMD)");
1193 
1194  // Note: modifyRequestedGroups is supposed to be in sync with the colvarproxy arrays,
1195  // and to stay that way during a simulation
1196 
1197  // compare this new group to those already allocated inside GlobalMaster
1198  int ig;
1199  for (ig = 0; ig < modifyRequestedGroups().size(); ig++) {
1200  AtomIDList const &namd_group = modifyRequestedGroups()[ig];
1201  bool b_match = true;
1202 
1203  if (namd_group.size() != ((int) atoms_ids.size())) {
1204  b_match = false;
1205  } else {
1206  int ia;
1207  for (ia = 0; ia < namd_group.size(); ia++) {
1208  int const aid = atoms_ids[ia];
1209  if (namd_group[ia] != aid) {
1210  b_match = false;
1211  break;
1212  }
1213  }
1214  }
1215 
1216  if (b_match) {
1217  if (cvm::debug())
1218  cvm::log("Group was already added.\n");
1219  // this group already exists
1220  atom_groups_refcount[ig] += 1;
1221  return ig;
1222  }
1223  }
1224 
1225  // add this group (note: the argument of add_atom_group_slot() is redundant for NAMD, and provided only for consistency)
1226  size_t const index = add_atom_group_slot(atom_groups_ids.size());
1227  modifyRequestedGroups().resize(atom_groups_ids.size());
1228  // the following is done in calculate()
1229  // modifyGroupForces().resize(atom_groups_ids.size());
1230  AtomIDList &namd_group = modifyRequestedGroups()[index];
1231  namd_group.resize(atoms_ids.size());
1232  int const n_all_atoms = Node::Object()->molecule->numAtoms;
1233  for (size_t ia = 0; ia < atoms_ids.size(); ia++) {
1234  int const aid = atoms_ids[ia];
1235  if (cvm::debug())
1236  cvm::log("Adding atom "+cvm::to_str(aid+1)+
1237  " for collective variables calculation.\n");
1238  if ( (aid < 0) || (aid >= n_all_atoms) ) {
1239  cvm::error("Error: invalid atom number specified, "+
1240  cvm::to_str(aid+1)+"\n", COLVARS_INPUT_ERROR);
1241  return -1;
1242  }
1243  namd_group[ia] = aid;
1244  }
1245 
1246  update_group_properties(index);
1247 
1248  if (cvm::debug()) {
1249  cvm::log("Group has index "+cvm::to_str(index)+"\n");
1250  cvm::log("modifyRequestedGroups length = "+cvm::to_str(modifyRequestedGroups().size())+
1251  ", modifyGroupForces length = "+cvm::to_str(modifyGroupForces().size())+"\n");
1252  }
1253 
1254  return index;
1255 }
static Node * Object()
Definition: Node.h:86
int size(void) const
Definition: ResizeArray.h:131
static int debug
Definition: parm.C:36
void resize(int i)
Definition: ResizeArray.h:84
int numAtoms
Definition: Molecule.h:585
ResizeArray< AtomIDList > & modifyRequestedGroups()
Definition: GlobalMaster.C:185
int update_group_properties(int index)
ForceList & modifyGroupForces()
Definition: GlobalMaster.C:168
Molecule * molecule
Definition: Node.h:179

◆ init_atoms_map()

void colvarproxy_namd::init_atoms_map ( )

Allocate an atoms map with the same size as the NAMD topology.

Definition at line 182 of file colvarproxy_namd.C.

References atoms_map, Node::molecule, Molecule::numAtoms, and Node::Object().

Referenced by colvarproxy_namd(), and update_atoms_map().

183 {
184  size_t const n_all_atoms = Node::Object()->molecule->numAtoms;
185  atoms_map.assign(n_all_atoms, -1);
186 }
static Node * Object()
Definition: Node.h:86
std::vector< int > atoms_map
Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data...
int numAtoms
Definition: Molecule.h:585
Molecule * molecule
Definition: Node.h:179

◆ init_tcl_pointers()

void colvarproxy_namd::init_tcl_pointers ( )
override

Definition at line 621 of file colvarproxy_namd.C.

References Node::Object().

622 {
623 #ifdef NAMD_TCL
624  // Store pointer to NAMD's Tcl interpreter
625  set_tcl_interp(Node::Object()->getScript()->interp);
626 #else
627  colvarproxy::init_tcl_pointers(); // Create dedicated interpreter
628 #endif
629 }
static Node * Object()
Definition: Node.h:86

◆ init_volmap_by_id()

int colvarproxy_namd::init_volmap_by_id ( int  volmap_id)
override

Definition at line 1311 of file colvarproxy_namd.C.

References ResizeArray< Elem >::add(), check_volmap_by_id(), and GlobalMaster::modifyRequestedGridObjects().

1312 {
1313  for (size_t i = 0; i < volmaps_ids.size(); i++) {
1314  if (volmaps_ids[i] == volmap_id) {
1315  // this map has already been requested
1316  volmaps_refcount[i] += 1;
1317  return i;
1318  }
1319  }
1320 
1321  int error_code = check_volmap_by_id(volmap_id);
1322  int index = -1;
1323  if (error_code == COLVARS_OK) {
1324  index = add_volmap_slot(volmap_id);
1325  modifyRequestedGridObjects().add(volmap_id);
1326  }
1327 
1328  return (error_code == COLVARS_OK) ? index : -1;
1329 }
int add(const Elem &elem)
Definition: ResizeArray.h:101
IntList & modifyRequestedGridObjects()
Definition: GlobalMaster.C:174
int check_volmap_by_id(int volmap_id) override

◆ init_volmap_by_name()

int colvarproxy_namd::init_volmap_by_name ( const char *  volmap_name)
override

Definition at line 1332 of file colvarproxy_namd.C.

References ResizeArray< Elem >::add(), check_volmap_by_name(), Molecule::get_gridfrc_grid(), GridforceGrid::get_scale(), MGridforceParamsList::index_for_key(), SimParameters::mgridforcelist, GlobalMaster::modifyRequestedGridObjects(), Node::molecule, Node::Object(), simparams, Vector::x, Vector::y, and Vector::z.

1333 {
1334  if (volmap_name == NULL) {
1335  return cvm::error("Error: no grid object name provided.", COLVARS_INPUT_ERROR);
1336  }
1337 
1338  int error_code = COLVARS_OK;
1339 
1340  error_code |= check_volmap_by_name(volmap_name);
1341 
1342  int index = -1;
1343  if (error_code == COLVARS_OK) {
1344 
1345  int volmap_id = simparams->mgridforcelist.index_for_key(volmap_name);
1346 
1347  // Check that the scale factor is correctly set to zero
1348  Molecule *mol = Node::Object()->molecule;
1349  GridforceGrid const *grid = mol->get_gridfrc_grid(volmap_id);
1350  Vector const gfScale = grid->get_scale();
1351  if ((gfScale.x != 0.0) || (gfScale.y != 0.0) || (gfScale.z != 0.0)) {
1352  error_code |= cvm::error("Error: GridForce map \""+
1353  std::string(volmap_name)+
1354  "\" has non-zero scale factors.\n",
1355  COLVARS_INPUT_ERROR);
1356  }
1357 
1358  for (size_t i = 0; i < volmaps_ids.size(); i++) {
1359  if (volmaps_ids[i] == volmap_id) {
1360  // this map has already been requested
1361  volmaps_refcount[i] += 1;
1362  return i;
1363  }
1364  }
1365 
1366  index = add_volmap_slot(volmap_id);
1367  modifyRequestedGridObjects().add(volmap_id);
1368  }
1369 
1370  return (error_code == COLVARS_OK) ? index : -1;
1371 }
static Node * Object()
Definition: Node.h:86
GridforceGrid * get_gridfrc_grid(int gridnum) const
Definition: Molecule.h:1363
SimParameters * simparams
Pointer to the NAMD simulation input object.
Definition: Vector.h:72
int check_volmap_by_name(char const *volmap_name) override
BigReal z
Definition: Vector.h:74
int add(const Elem &elem)
Definition: ResizeArray.h:101
Molecule stores the structural information for the system.
Definition: Molecule.h:175
int index_for_key(const char *key)
IntList & modifyRequestedGridObjects()
Definition: GlobalMaster.C:174
BigReal x
Definition: Vector.h:74
MGridforceParamsList mgridforcelist
virtual Vector get_scale(void) const =0
BigReal y
Definition: Vector.h:74
Molecule * molecule
Definition: Node.h:179

◆ load_atoms_pdb()

int colvarproxy_namd::load_atoms_pdb ( char const *  filename,
cvm::atom_group &  atoms,
std::string const &  pdb_field,
double const  pdb_field_value 
)
override

Definition at line 1023 of file colvarproxy_namd.C.

References PDB::atom(), cvm::atom, e_pdb_beta, e_pdb_occ, e_pdb_x, e_pdb_y, e_pdb_z, PDB::num_atoms(), and pdb_field_str2enum().

1027 {
1028  if (pdb_field_str.size() == 0)
1029  cvm::error("Error: must define which PDB field to use "
1030  "in order to define atoms from a PDB file.\n", COLVARS_INPUT_ERROR);
1031 
1032  PDB *pdb = new PDB(pdb_filename);
1033  size_t const pdb_natoms = pdb->num_atoms();
1034 
1035  e_pdb_field pdb_field_index = pdb_field_str2enum(pdb_field_str);
1036 
1037  for (size_t ipdb = 0; ipdb < pdb_natoms; ipdb++) {
1038 
1039  double atom_pdb_field_value = 0.0;
1040 
1041  switch (pdb_field_index) {
1042  case e_pdb_occ:
1043  atom_pdb_field_value = (pdb->atom(ipdb))->occupancy();
1044  break;
1045  case e_pdb_beta:
1046  atom_pdb_field_value = (pdb->atom(ipdb))->temperaturefactor();
1047  break;
1048  case e_pdb_x:
1049  atom_pdb_field_value = (pdb->atom(ipdb))->xcoor();
1050  break;
1051  case e_pdb_y:
1052  atom_pdb_field_value = (pdb->atom(ipdb))->ycoor();
1053  break;
1054  case e_pdb_z:
1055  atom_pdb_field_value = (pdb->atom(ipdb))->zcoor();
1056  break;
1057  default:
1058  break;
1059  }
1060 
1061  if ( (pdb_field_value) &&
1062  (atom_pdb_field_value != pdb_field_value) ) {
1063  continue;
1064  } else if (atom_pdb_field_value == 0.0) {
1065  continue;
1066  }
1067 
1068  if (atoms.is_enabled(colvardeps::f_ag_scalable)) {
1069  atoms.add_atom_id(ipdb);
1070  } else {
1071  atoms.add_atom(cvm::atom(ipdb+1));
1072  }
1073  }
1074 
1075  delete pdb;
1076  return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
1077 }
Definition: PDB.h:36
e_pdb_field pdb_field_str2enum(std::string const &pdb_field_str)
int num_atoms(void)
Definition: PDB.C:323
e_pdb_field
PDBAtom * atom(int place)
Definition: PDB.C:393

◆ load_coords_pdb()

int colvarproxy_namd::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

Definition at line 911 of file colvarproxy_namd.C.

References PDB::atom(), e_pdb_beta, e_pdb_occ, e_pdb_x, e_pdb_y, e_pdb_z, PDB::num_atoms(), and pdb_field_str2enum().

916 {
917  if (pdb_field_str.size() == 0 && indices.size() == 0) {
918  cvm::error("Bug alert: either PDB field should be defined or list of "
919  "atom IDs should be available when loading atom coordinates!\n", COLVARS_BUG_ERROR);
920  }
921 
922  e_pdb_field pdb_field_index;
923  bool const use_pdb_field = (pdb_field_str.size() > 0);
924  if (use_pdb_field) {
925  pdb_field_index = pdb_field_str2enum(pdb_field_str);
926  }
927 
928  // next index to be looked up in PDB file (if list is supplied)
929  std::vector<int>::const_iterator current_index = indices.begin();
930 
931  PDB *pdb = new PDB(pdb_filename);
932  size_t const pdb_natoms = pdb->num_atoms();
933 
934  if (pos.size() != pdb_natoms) {
935 
936  bool const pos_allocated = (pos.size() > 0);
937 
938  size_t ipos = 0, ipdb = 0;
939  for ( ; ipdb < pdb_natoms; ipdb++) {
940 
941  if (use_pdb_field) {
942  // PDB field mode: skip atoms with wrong value in PDB field
943  double atom_pdb_field_value = 0.0;
944 
945  switch (pdb_field_index) {
946  case e_pdb_occ:
947  atom_pdb_field_value = (pdb->atom(ipdb))->occupancy();
948  break;
949  case e_pdb_beta:
950  atom_pdb_field_value = (pdb->atom(ipdb))->temperaturefactor();
951  break;
952  case e_pdb_x:
953  atom_pdb_field_value = (pdb->atom(ipdb))->xcoor();
954  break;
955  case e_pdb_y:
956  atom_pdb_field_value = (pdb->atom(ipdb))->ycoor();
957  break;
958  case e_pdb_z:
959  atom_pdb_field_value = (pdb->atom(ipdb))->zcoor();
960  break;
961  default:
962  break;
963  }
964 
965  if ( (pdb_field_value) &&
966  (atom_pdb_field_value != pdb_field_value) ) {
967  continue;
968  } else if (atom_pdb_field_value == 0.0) {
969  continue;
970  }
971 
972  } else {
973  // Atom ID mode: use predefined atom IDs from the atom group
974  if (((int) ipdb) != *current_index) {
975  // Skip atoms not in the list
976  continue;
977  } else {
978  current_index++;
979  }
980  }
981 
982  if (!pos_allocated) {
983  pos.push_back(cvm::atom_pos(0.0, 0.0, 0.0));
984  } else if (ipos >= pos.size()) {
985  cvm::error("Error: the PDB file \""+
986  std::string(pdb_filename)+
987  "\" contains coordinates for "
988  "more atoms than needed.\n", COLVARS_BUG_ERROR);
989  }
990 
991  pos[ipos] = cvm::atom_pos((pdb->atom(ipdb))->xcoor(),
992  (pdb->atom(ipdb))->ycoor(),
993  (pdb->atom(ipdb))->zcoor());
994  ipos++;
995  if (!use_pdb_field && current_index == indices.end())
996  break;
997  }
998 
999  if (ipos < pos.size() || (!use_pdb_field && current_index != indices.end())) {
1000  size_t n_requested = use_pdb_field ? pos.size() : indices.size();
1001  cvm::error("Error: number of matching records in the PDB file \""+
1002  std::string(pdb_filename)+"\" ("+cvm::to_str(ipos)+
1003  ") does not match the number of requested coordinates ("+
1004  cvm::to_str(n_requested)+").\n", COLVARS_INPUT_ERROR);
1005  return COLVARS_ERROR;
1006  }
1007  } else {
1008 
1009  // when the PDB contains exactly the number of atoms of the array,
1010  // ignore the fields and just read coordinates
1011  for (size_t ia = 0; ia < pos.size(); ia++) {
1012  pos[ia] = cvm::atom_pos((pdb->atom(ia))->xcoor(),
1013  (pdb->atom(ia))->ycoor(),
1014  (pdb->atom(ia))->zcoor());
1015  }
1016  }
1017 
1018  delete pdb;
1019  return COLVARS_OK;
1020 }
Definition: PDB.h:36
e_pdb_field pdb_field_str2enum(std::string const &pdb_field_str)
int num_atoms(void)
Definition: PDB.C:323
e_pdb_field
PDBAtom * atom(int place)
Definition: PDB.C:393

◆ log()

void colvarproxy_namd::log ( std::string const &  message)
override

Definition at line 680 of file colvarproxy_namd.C.

References endi(), and iout.

Referenced by calculate(), error(), setup(), and update_atom_properties().

681 {
682  std::istringstream is(message);
683  std::string line;
684  while (std::getline(is, line))
685  iout << "colvars: " << line << "\n";
686  iout << endi;
687 }
std::ostream & endi(std::ostream &s)
Definition: InfoStream.C:54
#define iout
Definition: InfoStream.h:51

◆ num_replicas()

int colvarproxy_namd::num_replicas ( )
override

Definition at line 1623 of file colvarproxy_namd.C.

1623  {
1624  return CmiNumPartitions();
1625 }

◆ output_stream()

std::ostream & colvarproxy_namd::output_stream ( std::string const &  output_name,
std::string const  description 
)
override

Definition at line 1080 of file colvarproxy_namd.C.

References backup_file(), and debug.

1082 {
1083  if (cvm::debug()) {
1084  cvm::log("Using colvarproxy_namd::output_stream()\n");
1085  }
1086 
1087  if (!io_available()) {
1088  cvm::error("Error: trying to access an output file/channel "
1089  "from the wrong thread.\n", COLVARS_BUG_ERROR);
1090  return *output_stream_error_;
1091  }
1092 
1093  if (output_streams_.count(output_name) > 0) {
1094  return *(output_streams_[output_name]);
1095  }
1096 
1097  backup_file(output_name.c_str());
1098 
1099  output_streams_[output_name] = new ofstream_namd(output_name.c_str(), std::ios::binary);
1100  if (! output_streams_[output_name]->good()) {
1101  cvm::error("Error: cannot write to "+description+" \""+output_name+"\".\n",
1102  COLVARS_FILE_ERROR);
1103  }
1104 
1105  return *(output_streams_[output_name]);
1106 }
static int debug
Definition: parm.C:36
int backup_file(char const *filename) override

◆ position_distance()

cvm::rvector colvarproxy_namd::position_distance ( cvm::atom_pos const &  pos1,
cvm::atom_pos const &  pos2 
) const

Definition at line 849 of file colvarproxy_namd.C.

References Lattice::delta(), GlobalMaster::lattice, Vector::x, Vector::y, and Vector::z.

852 {
853  Position const p1(pos1.x, pos1.y, pos1.z);
854  Position const p2(pos2.x, pos2.y, pos2.z);
855  // return p2 - p1
856  Vector const d = this->lattice->delta(p2, p1);
857  return cvm::rvector(d.x, d.y, d.z);
858 }
Definition: Vector.h:72
BigReal z
Definition: Vector.h:74
BigReal x
Definition: Vector.h:74
const Lattice * lattice
Definition: GlobalMaster.h:142
BigReal y
Definition: Vector.h:74
NAMD_HOST_DEVICE Vector delta(const Position &pos1, const Position &pos2) const
Definition: Lattice.h:149

◆ rand_gaussian()

cvm::real colvarproxy_namd::rand_gaussian ( )
inlineoverride

Definition at line 109 of file colvarproxy_namd.h.

References Random::gaussian(), and random.

110  {
111  return random.gaussian();
112  }
Random random
NAMD-style PRNG object.
BigReal gaussian(void)
Definition: Random.h:116

◆ replica_comm_barrier()

void colvarproxy_namd::replica_comm_barrier ( )
override

Definition at line 1628 of file colvarproxy_namd.C.

References replica_barrier().

1628  {
1629  replica_barrier();
1630 }
void replica_barrier()

◆ replica_comm_recv()

int colvarproxy_namd::replica_comm_recv ( char *  msg_data,
int  buf_len,
int  src_rep 
)
override

Definition at line 1633 of file colvarproxy_namd.C.

References DataMessage::data, replica_recv(), and DataMessage::size.

1634  {
1635  DataMessage *recvMsg = NULL;
1636  replica_recv(&recvMsg, src_rep, CkMyPe());
1637  CmiAssert(recvMsg != NULL);
1638  int retval = recvMsg->size;
1639  if (buf_len >= retval) {
1640  memcpy(msg_data,recvMsg->data,retval);
1641  } else {
1642  retval = 0;
1643  }
1644  CmiFree(recvMsg);
1645  return retval;
1646 }
void replica_recv(DataMessage **precvMsg, int srcPart, int srcPE)
char data[1]
Definition: DataExchanger.h:23

◆ replica_comm_send()

int colvarproxy_namd::replica_comm_send ( char *  msg_data,
int  msg_len,
int  dest_rep 
)
override

Definition at line 1649 of file colvarproxy_namd.C.

References replica_send().

1650  {
1651  replica_send(msg_data, msg_len, dest_rep, CkMyPe());
1652  return msg_len;
1653 }
void replica_send(const char *sndbuf, int sendcount, int destPart, int destPE)

◆ replica_index()

int colvarproxy_namd::replica_index ( )
override

Definition at line 1618 of file colvarproxy_namd.C.

1618  {
1619  return CmiMyPartition();
1620 }

◆ request_total_force()

void colvarproxy_namd::request_total_force ( bool  yesno)
override

Definition at line 667 of file colvarproxy_namd.C.

References debug, and GlobalMaster::requestTotalForce().

668 {
669  if (cvm::debug()) {
670  cvm::log("colvarproxy_namd::request_total_force()\n");
671  }
672  total_force_requested = yesno;
673  requestTotalForce(total_force_requested);
674  if (cvm::debug()) {
675  cvm::log("colvarproxy_namd::request_total_force() end\n");
676  }
677 }
static int debug
Definition: parm.C:36
void requestTotalForce(bool yesno=true)
Definition: GlobalMaster.h:134

◆ reset()

int colvarproxy_namd::reset ( )
override

Definition at line 304 of file colvarproxy_namd.C.

References atoms_map, ResizeArray< Elem >::clear(), debug, GlobalMaster::modifyAppliedForces(), GlobalMaster::modifyForcedAtoms(), GlobalMaster::modifyGridObjForces(), GlobalMaster::modifyGroupForces(), GlobalMaster::modifyRequestedAtoms(), GlobalMaster::modifyRequestedGridObjects(), GlobalMaster::modifyRequestedGroups(), and GlobalMaster::requestTotalForce().

305 {
306  if (cvm::debug()) {
307  cvm::log("colvarproxy_namd::reset()\n");
308  }
309 
310  int error_code = COLVARS_OK;
311 
312  // Unrequest all positions, total forces, etc from NAMD
316 
317  modifyRequestedGroups().clear();
319 
320 #if NAMD_VERSION_NUMBER >= 34471681
323 #endif
324 
325  requestTotalForce(false);
326 
327  atoms_map.clear();
328 
329  // Clear internal atomic data
330  error_code |= colvarproxy::reset();
331 
332  return error_code;
333 }
ForceList & modifyAppliedForces()
Definition: GlobalMaster.C:163
AtomIDList & modifyRequestedAtoms()
Definition: GlobalMaster.C:128
BigRealList & modifyGridObjForces()
Definition: GlobalMaster.C:180
static int debug
Definition: parm.C:36
void clear()
Definition: ResizeArray.h:91
std::vector< int > atoms_map
Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data...
IntList & modifyRequestedGridObjects()
Definition: GlobalMaster.C:174
AtomIDList & modifyForcedAtoms()
Definition: GlobalMaster.C:158
ResizeArray< AtomIDList > & modifyRequestedGroups()
Definition: GlobalMaster.C:185
ForceList & modifyGroupForces()
Definition: GlobalMaster.C:168
void requestTotalForce(bool yesno=true)
Definition: GlobalMaster.h:134

◆ run_colvar_callback()

int colvarproxy_namd::run_colvar_callback ( std::string const &  name,
std::vector< const colvarvalue *> const &  cvcs,
colvarvalue &  value 
)
override

Definition at line 636 of file colvarproxy_namd.C.

640 {
641  return colvarproxy::tcl_run_colvar_callback(name, cvc_values, value);
642 }

◆ run_colvar_gradient_callback()

int colvarproxy_namd::run_colvar_gradient_callback ( std::string const &  name,
std::vector< const colvarvalue *> const &  cvcs,
std::vector< cvm::matrix2d< cvm::real > > &  gradient 
)
override

Definition at line 644 of file colvarproxy_namd.C.

648 {
649  return colvarproxy::tcl_run_colvar_gradient_callback(name, cvc_values,
650  gradient);
651 }

◆ run_force_callback()

int colvarproxy_namd::run_force_callback ( )
override

Definition at line 631 of file colvarproxy_namd.C.

632 {
633  return colvarproxy::tcl_run_force_callback();
634 }

◆ scalable_group_coms()

int colvarproxy_namd::scalable_group_coms ( )
inlineoverride

Definition at line 205 of file colvarproxy_namd.h.

206  {
207  return COLVARS_OK;
208  }

◆ set_unit_system()

int colvarproxy_namd::set_unit_system ( std::string const &  units_in,
bool  check_only 
)
override

Definition at line 1292 of file colvarproxy_namd.C.

1293 {
1294  if (units_in != "real") {
1295  cvm::error("Error: Specified unit system \"" + units_in + "\" is unsupported in NAMD. Supported units are \"real\" (A, kcal/mol).\n");
1296  return COLVARS_ERROR;
1297  }
1298  return COLVARS_OK;
1299 }

◆ setup()

int colvarproxy_namd::setup ( )
override

Definition at line 225 of file colvarproxy_namd.C.

References SimParameters::dt, log(), GlobalMaster::modifyGridObjForces(), GlobalMaster::modifyRequestedGroups(), simparams, ResizeArray< Elem >::size(), update_atom_properties(), update_group_properties(), update_target_temperature(), and SimParameters::wrapAll.

Referenced by calculate(), and colvarproxy_namd().

226 {
227  int error_code = colvarproxy::setup();
228 
229  if (colvars->size() == 0) {
230  // Module is empty, nothing to do
231  return COLVARS_OK;
232  }
233 
234  log("Updating NAMD interface:\n");
235 
236  errno = 0;
237 
238  if (simparams->wrapAll) {
239  log("Warning: enabling wrapAll can lead to inconsistent results "
240  "for Colvars calculations: please disable wrapAll, "
241  "as is the default option in NAMD.\n");
242  }
243 
244  log("updating atomic data ("+cvm::to_str(atoms_ids.size())+" atoms).\n");
245 
246  size_t i;
247  for (i = 0; i < atoms_ids.size(); i++) {
249 
250  // zero out mutable arrays
251  atoms_positions[i] = cvm::rvector(0.0, 0.0, 0.0);
252  atoms_total_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
253  atoms_new_colvar_forces[i] = cvm::rvector(0.0, 0.0, 0.0);
254  }
255 
256  size_t n_group_atoms = 0;
257  for (int ig = 0; ig < modifyRequestedGroups().size(); ig++) {
258  n_group_atoms += modifyRequestedGroups()[ig].size();
259  }
260 
261  log("updating group data ("+cvm::to_str(atom_groups_ids.size())+
262  " scalable groups, "+
263  cvm::to_str(n_group_atoms)+" atoms in total).\n");
264 
265  // Note: groupMassBegin, groupMassEnd may be used here, but they won't work for charges
266  for (int ig = 0; ig < modifyRequestedGroups().size(); ig++) {
267 
268  // update mass and charge
270 
271  atom_groups_coms[ig] = cvm::rvector(0.0, 0.0, 0.0);
272  atom_groups_total_forces[ig] = cvm::rvector(0.0, 0.0, 0.0);
273  atom_groups_new_colvar_forces[ig] = cvm::rvector(0.0, 0.0, 0.0);
274  }
275 
276 #if NAMD_VERSION_NUMBER >= 34471681
277  log("updating grid object data ("+cvm::to_str(volmaps_ids.size())+
278  " grid objects in total).\n");
279  for (int imap = 0; imap < modifyGridObjForces().size(); imap++) {
280  volmaps_new_colvar_forces[imap] = 0.0;
281  }
282 #endif
283 
284  size_t const new_features_hash =
285  std::hash<std::string>{}(colvars->feature_report(0));
286  if (new_features_hash != features_hash) {
287  // Nag only once, there may be many run commands
288  log(std::string("\n")+colvars->feature_report(0)+std::string("\n"));
289  features_hash = new_features_hash;
290  }
291 
293  log("updating target temperature (T = "+
294  cvm::to_str(target_temperature())+" K).\n");
295 
296  // Note: not needed currently, but may be in the future if NAMD allows
297  // redefining the timestep
298  set_integration_timestep(simparams->dt);
299 
300  return error_code;
301 }
BigRealList & modifyGridObjForces()
Definition: GlobalMaster.C:180
int size(void) const
Definition: ResizeArray.h:131
SimParameters * simparams
Pointer to the NAMD simulation input object.
int update_target_temperature()
Get the target temperature from the NAMD thermostats supported so far.
void log(std::string const &message) override
ResizeArray< AtomIDList > & modifyRequestedGroups()
Definition: GlobalMaster.C:185
int update_group_properties(int index)
void update_atom_properties(int index)

◆ total_forces_enabled()

bool colvarproxy_namd::total_forces_enabled ( ) const
inlineoverride

Definition at line 96 of file colvarproxy_namd.h.

97  {
98  return total_force_requested;
99  }

◆ update_accelMD_info()

void colvarproxy_namd::update_accelMD_info ( )
protected

Definition at line 607 of file colvarproxy_namd.C.

References Controller::accelMDdV, accelMDOn, amd_weight_factor, NamdState::getController(), Node::Object(), and Node::state.

Referenced by calculate().

607  {
608  if (accelMDOn == false) {
609  return;
610  }
611  const Controller& c = Node::Object()->state->getController();
612  // This aMD factor is from previous step!
613  amd_weight_factor = std::exp(c.accelMDdV /
614  (target_temperature() * boltzmann()));
615 // std::cout << "Step: " << cvm::to_str(colvars->it) << " accelMD dV in colvars: " << c.accelMDdV << std::endl;
616 }
static Node * Object()
Definition: Node.h:86
const Controller & getController() const
Definition: NamdState.h:51
cvm::real amd_weight_factor
bool accelMDOn
Accelerated MD reweighting factor.
NamdState * state
Definition: Node.h:184
BigReal accelMDdV
Definition: Controller.h:117

◆ update_atom_properties()

void colvarproxy_namd::update_atom_properties ( int  index)

Definition at line 833 of file colvarproxy_namd.C.

References Molecule::atomcharge(), Molecule::atommass(), log(), Node::molecule, and Node::Object().

Referenced by init_atom(), setup(), and update_atoms_map().

834 {
835  Molecule *mol = Node::Object()->molecule;
836  // update mass
837  double const mass = mol->atommass(atoms_ids[index]);
838  if (mass <= 0.001) {
839  this->log("Warning: near-zero mass for atom "+
840  cvm::to_str(atoms_ids[index]+1)+
841  "; expect unstable dynamics if you apply forces to it.\n");
842  }
843  atoms_masses[index] = mass;
844  // update charge
845  atoms_charges[index] = mol->atomcharge(atoms_ids[index]);
846 }
static Node * Object()
Definition: Node.h:86
Real atomcharge(int anum) const
Definition: Molecule.h:1117
Molecule stores the structural information for the system.
Definition: Molecule.h:175
void log(std::string const &message) override
Real atommass(int anum) const
Definition: Molecule.h:1107
Molecule * molecule
Definition: Node.h:179

◆ update_atoms_map()

int colvarproxy_namd::update_atoms_map ( AtomIDList::const_iterator  begin,
AtomIDList::const_iterator  end 
)

Definition at line 189 of file colvarproxy_namd.C.

References ResizeArray< Elem >::add(), atoms_map, debug, init_atoms_map(), GlobalMaster::modifyRequestedAtoms(), Node::molecule, Molecule::numAtoms, Node::Object(), and update_atom_properties().

Referenced by calculate().

191 {
192  if (atoms_map.size() != Node::Object()->molecule->numAtoms) {
193  init_atoms_map();
194  }
195 
196  for (AtomIDList::const_iterator a_i = begin; a_i != end; a_i++) {
197 
198  if (atoms_map[*a_i] >= 0) continue;
199 
200  for (size_t i = 0; i < atoms_ids.size(); i++) {
201  if (atoms_ids[i] == *a_i) {
202  atoms_map[*a_i] = i;
203  break;
204  }
205  }
206 
207  if (atoms_map[*a_i] < 0) {
208  // this atom is probably managed by another GlobalMaster:
209  // add it here anyway to avoid having to test for array boundaries at each step
210  int const index = add_atom_slot(*a_i);
211  atoms_map[*a_i] = index;
212  modifyRequestedAtoms().add(*a_i);
213  update_atom_properties(index);
214  }
215  }
216 
217  if (cvm::debug()) {
218  cvm::log("atoms_map = "+cvm::to_str(atoms_map)+".\n");
219  }
220 
221  return COLVARS_OK;
222 }
static Node * Object()
Definition: Node.h:86
AtomIDList & modifyRequestedAtoms()
Definition: GlobalMaster.C:128
static int debug
Definition: parm.C:36
int add(const Elem &elem)
Definition: ResizeArray.h:101
std::vector< int > atoms_map
Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data...
const Elem * const_iterator
Definition: ResizeArray.h:38
void init_atoms_map()
Allocate an atoms map with the same size as the NAMD topology.
int numAtoms
Definition: Molecule.h:585
void update_atom_properties(int index)
Molecule * molecule
Definition: Node.h:179

◆ update_group_properties()

int colvarproxy_namd::update_group_properties ( int  index)

Definition at line 1265 of file colvarproxy_namd.C.

References Molecule::atomcharge(), Molecule::atommass(), debug, GlobalMaster::modifyRequestedGroups(), Node::molecule, Node::Object(), and ResizeArray< Elem >::size().

Referenced by init_atom_group(), and setup().

1266 {
1267  AtomIDList const &namd_group = modifyRequestedGroups()[index];
1268  if (cvm::debug()) {
1269  cvm::log("Re-calculating total mass and charge for scalable group no. "+cvm::to_str(index+1)+" ("+
1270  cvm::to_str(namd_group.size())+" atoms).\n");
1271  }
1272 
1273  cvm::real total_mass = 0.0;
1274  cvm::real total_charge = 0.0;
1275  for (int i = 0; i < namd_group.size(); i++) {
1276  total_mass += Node::Object()->molecule->atommass(namd_group[i]);
1277  total_charge += Node::Object()->molecule->atomcharge(namd_group[i]);
1278  }
1279  atom_groups_masses[index] = total_mass;
1280  atom_groups_charges[index] = total_charge;
1281 
1282  if (cvm::debug()) {
1283  cvm::log("total mass = "+cvm::to_str(total_mass)+
1284  ", total charge = "+cvm::to_str(total_charge)+"\n");
1285  }
1286 
1287  return COLVARS_OK;
1288 }
static Node * Object()
Definition: Node.h:86
Real atomcharge(int anum) const
Definition: Molecule.h:1117
int size(void) const
Definition: ResizeArray.h:131
static int debug
Definition: parm.C:36
ResizeArray< AtomIDList > & modifyRequestedGroups()
Definition: GlobalMaster.C:185
Real atommass(int anum) const
Definition: Molecule.h:1107
Molecule * molecule
Definition: Node.h:179

◆ update_target_temperature()

int colvarproxy_namd::update_target_temperature ( )

Get the target temperature from the NAMD thermostats supported so far.

Definition at line 159 of file colvarproxy_namd.C.

References SimParameters::langevinOn, SimParameters::langevinTemp, SimParameters::loweAndersenOn, SimParameters::loweAndersenTemp, SimParameters::reassignFreq, SimParameters::reassignTemp, SimParameters::rescaleFreq, SimParameters::rescaleTemp, simparams, SimParameters::stochRescaleOn, SimParameters::stochRescaleTemp, SimParameters::tCoupleOn, and SimParameters::tCoupleTemp.

Referenced by colvarproxy_namd(), and setup().

160 {
161  int error_code = COLVARS_OK;
162  if (simparams->rescaleFreq > 0) {
163  error_code |= set_target_temperature(simparams->rescaleTemp);
164  } else if (simparams->reassignFreq > 0) {
165  error_code |= set_target_temperature(simparams->reassignTemp);
166  } else if (simparams->langevinOn) {
167  error_code |= set_target_temperature(simparams->langevinTemp);
168  } else if (simparams->tCoupleOn) {
169  error_code |= set_target_temperature(simparams->tCoupleTemp);
170  } else if (simparams->loweAndersenOn) {
171  error_code |= set_target_temperature(simparams->loweAndersenTemp);
172  } else if (simparams->stochRescaleOn) {
173  error_code |= set_target_temperature(simparams->stochRescaleTemp);
174  } else {
175  error_code |= set_target_temperature(0.0);
176  }
177  return error_code;
178 }
SimParameters * simparams
Pointer to the NAMD simulation input object.
BigReal reassignTemp
BigReal rescaleTemp
BigReal langevinTemp
BigReal tCoupleTemp
BigReal stochRescaleTemp
BigReal loweAndersenTemp

Friends And Related Function Documentation

◆ cvm::atom

friend class cvm::atom
friend

Definition at line 70 of file colvarproxy_namd.h.

Referenced by load_atoms_pdb().

Member Data Documentation

◆ accelMDOn

bool colvarproxy_namd::accelMDOn
protected

Accelerated MD reweighting factor.

Definition at line 62 of file colvarproxy_namd.h.

Referenced by accelMD_enabled(), colvarproxy_namd(), and update_accelMD_info().

◆ amd_weight_factor

cvm::real colvarproxy_namd::amd_weight_factor
protected

Definition at line 63 of file colvarproxy_namd.h.

Referenced by colvarproxy_namd(), get_accelMD_factor(), and update_accelMD_info().

◆ atoms_map

std::vector<int> colvarproxy_namd::atoms_map
protected

Array of atom indices (relative to the colvarproxy arrays), usedfor faster copy of atomic data.

Definition at line 44 of file colvarproxy_namd.h.

Referenced by calculate(), init_atom(), init_atoms_map(), reset(), and update_atoms_map().

◆ first_timestep

bool colvarproxy_namd::first_timestep
protected

Definition at line 52 of file colvarproxy_namd.h.

Referenced by calculate(), and colvarproxy_namd().

◆ previous_NAMD_step

cvm::step_number colvarproxy_namd::previous_NAMD_step
protected

Definition at line 53 of file colvarproxy_namd.h.

Referenced by calculate().

◆ random

Random colvarproxy_namd::random
protected

NAMD-style PRNG object.

Definition at line 50 of file colvarproxy_namd.h.

Referenced by colvarproxy_namd(), and rand_gaussian().

◆ reduction

SubmitReduction* colvarproxy_namd::reduction
protected

Used to submit restraint energy as MISC.

Definition at line 56 of file colvarproxy_namd.h.

Referenced by add_energy(), calculate(), colvarproxy_namd(), and ~colvarproxy_namd().

◆ simparams

SimParameters* colvarproxy_namd::simparams
protected

The documentation for this class was generated from the following files: