37 #include "colvarmodule.h"    39 #include "colvarbias.h"    40 #include "colvaratoms.h"    41 #include "colvarproxy.h"    44 #include "colvarscript.h"    49   engine_name_ = 
"NAMD";
    50 #if CMK_SMP && USE_CKLOOP    51   charm_lock_state = CmiCreateLock();
    56   if ( 0 == CkMyPe() ) {
    65   boltzmann_ = 0.001987191;
    73     iout << 
"Info: initializing the colvars proxy object.\n" << 
endi;
    80   colvarproxy_io::set_input_prefix(input_restart ? input_restart->
data : 
"");
    89   updated_masses_ = updated_charges_ = 
true;
   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");
   112   colvars = 
new colvarmodule(
this);
   114   cvm::log(
"Using NAMD interface, version "+
   116   colvars->cite_feature(
"NAMD engine");
   117   colvars->cite_feature(
"Colvars-NAMD interface");
   120   for ( ; config; config = config->
next ) {
   121     add_config(
"configfile", config->
data);
   125   colvarproxy::parse_module_config();
   127   colvars->update_engine_parameters();
   128   colvars->setup_input();
   129   colvars->setup_output();
   137   have_scripts = 
false;
   144 #if !defined (NAMD_UNIFIED_REDUCTION)   148   #if defined(NODEGROUP_FORCE_REGISTER) && !defined(NAMD_UNIFIED_REDUCTION)   149   CProxy_PatchData cpdata(CkpvAccess(BOCclass_group).patchData);
   150   PatchData *patchData = cpdata.ckLocalBranch();
   155     iout << 
"Info: done initializing the colvars proxy object.\n" << 
endi;
   161 #if CMK_SMP && USE_CKLOOP   162   CmiDestroyLock(charm_lock_state);
   164 #if !defined (NAMD_UNIFIED_REDUCTION)   172   int error_code = COLVARS_OK;
   186     error_code |= set_target_temperature(0.0);
   211     for (
size_t i = 0; i < atoms_ids.size(); i++) {
   212       if (atoms_ids[i] == *a_i) {
   221       int const index = add_atom_slot(*a_i);
   229     cvm::log(
"atoms_map = "+cvm::to_str(
atoms_map)+
".\n");
   238   int error_code = colvarproxy::setup();
   240   if (colvars->size() == 0) {
   245   log(
"Updating NAMD interface:\n");
   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");
   255   log(
"updating atomic data ("+cvm::to_str(atoms_ids.size())+
" atoms).\n");
   258   for (i = 0; i < atoms_ids.size(); i++) {
   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);
   267   size_t n_group_atoms = 0;
   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");
   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);
   287 #if NAMD_VERSION_NUMBER >= 34471681   288   log(
"updating grid object data ("+cvm::to_str(volmaps_ids.size())+
   289       " grid objects in total).\n");
   291     volmaps_new_colvar_forces[imap] = 0.0;
   295   size_t const new_features_hash =
   296     std::hash<std::string>{}(colvars->feature_report(0));
   297   if (new_features_hash != features_hash) {
   299     log(std::string(
"\n")+colvars->feature_report(0)+std::string(
"\n"));
   300     features_hash = new_features_hash;
   304   log(
"updating target temperature (T = "+
   305       cvm::to_str(target_temperature())+
" K).\n");
   318     cvm::log(
"colvarproxy_namd::reset()\n");
   321   int error_code = COLVARS_OK;
   331 #if NAMD_VERSION_NUMBER >= 34471681   341   error_code |= colvarproxy::reset();
   356     colvars->update_engine_parameters();
   357     colvars->setup_input();
   358     colvars->setup_output();
   368       colvars->it += time_step_factor();
   369       b_simulation_continuing = 
false;
   377       b_simulation_continuing = 
true;
   383       colvars->setup_output();
   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);
   401     boundaries_type = boundaries_non_periodic;
   405       boundaries_type = boundaries_pbc_ortho;
   407       boundaries_type = boundaries_pbc_triclinic;
   409     colvarproxy_system::update_pbc_lattice();
   411     boundaries_type = boundaries_unsupported;
   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");
   429   if ( (
atoms_map.size() != n_all_atoms) ||
   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);
   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);
   448 #if NAMD_VERSION_NUMBER >= 34471681   449   for (
int imap = 0; imap < volmaps_ids.size(); imap++) {
   450     volmaps_new_colvar_forces[imap] = 0.0;
   456       cvm::log(
"Updating positions arrays.\n");
   458     size_t n_positions = 0;
   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);
   474   if (total_force_requested && cvm::step_relative() > 0) {
   481         cvm::log(
"Updating total forces arrays.\n");
   483       size_t n_total_forces = 0;
   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);
   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);
   505         cvm::log(
"Updating group total forces arrays.\n");
   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);
   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);
   526       cvm::log(
"Updating group positions arrays.\n");
   534       atom_groups_coms[ig] = cvm::rvector(gp_i->
x, gp_i->
y, gp_i->
z);
   538 #if NAMD_VERSION_NUMBER >= 34471681   541       log(
"Updating grid objects.\n");
   548       for (
size_t imap = 0; imap < volmaps_ids.size(); imap++) {
   549         if (volmaps_ids[imap] == *goi_i) {
   550           volmaps_values[imap] = *gov_i;
   559     print_input_atomic_data();
   563   if (colvars->calc() != COLVARS_OK) {
   564     cvm::error(
"Error in the collective variables module.\n", COLVARS_ERROR);
   568     print_output_atomic_data();
   572   for (
size_t i = 0; i < atoms_ids.size(); i++) {
   573     cvm::rvector 
const &f = atoms_new_colvar_forces[i];
   578   if (atom_groups_new_colvar_forces.size() > 0) {
   582       cvm::rvector 
const &f = atom_groups_new_colvar_forces[ig];
   583       *gf_i = 
Vector(f.x, f.y, f.z);
   587 #if NAMD_VERSION_NUMBER >= 34471681   588   if (volmaps_new_colvar_forces.size() > 0) {
   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];
   605   #if defined(NODEGROUP_FORCE_REGISTER) && !defined(NAMD_UNIFIED_REDUCTION)   610   #if !defined(NAMD_UNIFIED_REDUCTION)   638   colvarproxy::init_tcl_pointers(); 
   644   return colvarproxy::tcl_run_force_callback();
   648                           std::string 
const &name,
   649                           std::vector<const colvarvalue *> 
const &cvc_values,
   652   return colvarproxy::tcl_run_colvar_callback(name, cvc_values, value);
   656                           std::string 
const &name,
   657                           std::vector<const colvarvalue *> 
const &cvc_values,
   658                           std::vector<cvm::matrix2d<cvm::real> > &gradient)
   660   return colvarproxy::tcl_run_colvar_gradient_callback(name, cvc_values,
   667   #if defined(NODEGROUP_FORCE_REGISTER) && !defined(NAMD_UNIFIED_REDUCTION)   674   #if !defined(NAMD_UNIFIED_REDUCTION)   685     cvm::log(
"colvarproxy_namd::request_total_force()\n");
   687   total_force_requested = yesno;
   690     cvm::log(
"colvarproxy_namd::request_total_force() end\n");
   697   std::istringstream is(message);
   699   while (std::getline(is, line))
   700     iout << 
"colvars: " << line << 
"\n";
   708   switch (cvm::get_error()) {
   709   case COLVARS_FILE_ERROR:
   711   case COLVARS_NOT_IMPLEMENTED:
   712     errno = ENOSYS; 
break;
   713   case COLVARS_MEMORY_ERROR:
   714     errno = ENOMEM; 
break;
   716   char const *msg = 
"Error in the collective variables module "   717     "(see above for details)";
   729   int const aid = (atom_number-1);
   732     cvm::log(
"Adding atom "+cvm::to_str(atom_number)+
   733         " for collective variables calculation.\n");
   736     cvm::error(
"Error: invalid atom number specified, "+
   737                cvm::to_str(atom_number)+
"\n", COLVARS_INPUT_ERROR);
   738     return COLVARS_INPUT_ERROR;
   755   int aid = (atom_number-1);
   757   for (
size_t i = 0; i < atoms_ids.size(); i++) {
   758     if (atoms_ids[i] == aid) {
   760       atoms_refcount[i] += 1;
   768     return COLVARS_INPUT_ERROR;
   771   int const index = add_atom_slot(aid);
   780                                     std::string 
const     &atom_name,
   781                                     std::string 
const     &segment_id)
   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+
"\"") :
   800                "\n", COLVARS_INPUT_ERROR);
   801     return COLVARS_INPUT_ERROR;
   813                                 std::string 
const     &atom_name,
   814                                 std::string 
const     &segment_id)
   816   int const aid = 
check_atom_id(residue, atom_name, segment_id);
   818   for (
size_t i = 0; i < atoms_ids.size(); i++) {
   819     if (atoms_ids[i] == aid) {
   821       atoms_refcount[i] += 1;
   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");
   833   int const index = add_atom_slot(aid);
   843   colvarproxy::clear_atom(index);
   852   double const mass = mol->
atommass(atoms_ids[index]);
   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");
   858   atoms_masses[index] = mass;
   860   atoms_charges[index] = mol->
atomcharge(atoms_ids[index]);
   865                                                  cvm::atom_pos 
const &pos2)
   868   Position const p1(pos1.x, pos1.y, pos1.z);
   869   Position const p2(pos2.x, pos2.y, pos2.z);
   872   return cvm::rvector(d.
x, d.
y, d.
z);
   892   if (colvarparse::to_lower_cppstr(pdb_field_str) ==
   893       colvarparse::to_lower_cppstr(
"O")) {
   897   if (colvarparse::to_lower_cppstr(pdb_field_str) ==
   898       colvarparse::to_lower_cppstr(
"B")) {
   902   if (colvarparse::to_lower_cppstr(pdb_field_str) ==
   903       colvarparse::to_lower_cppstr(
"X")) {
   907   if (colvarparse::to_lower_cppstr(pdb_field_str) ==
   908       colvarparse::to_lower_cppstr(
"Y")) {
   912   if (colvarparse::to_lower_cppstr(pdb_field_str) ==
   913       colvarparse::to_lower_cppstr(
"Z")) {
   918     cvm::error(
"Error: unsupported PDB field, \""+
   919                pdb_field_str+
"\".\n", COLVARS_INPUT_ERROR);
   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)
   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);
   938   bool const use_pdb_field = (pdb_field_str.size() > 0);
   944   std::vector<int>::const_iterator current_index = indices.begin();
   946   PDB *pdb = 
new PDB(pdb_filename);
   947   size_t const pdb_natoms = pdb->
num_atoms();
   949   if (pos.size() != pdb_natoms) {
   951     bool const pos_allocated = (pos.size() > 0);
   953     size_t ipos = 0, ipdb = 0;
   954     for ( ; ipdb < pdb_natoms; ipdb++) {
   958         double atom_pdb_field_value = 0.0;
   960         switch (pdb_field_index) {
   962           atom_pdb_field_value = (pdb->
atom(ipdb))->occupancy();
   965           atom_pdb_field_value = (pdb->
atom(ipdb))->temperaturefactor();
   968           atom_pdb_field_value = (pdb->
atom(ipdb))->xcoor();
   971           atom_pdb_field_value = (pdb->
atom(ipdb))->ycoor();
   974           atom_pdb_field_value = (pdb->
atom(ipdb))->zcoor();
   980         if ( (pdb_field_value) &&
   981              (atom_pdb_field_value != pdb_field_value) ) {
   983         } 
else if (atom_pdb_field_value == 0.0) {
   989         if (((
int) ipdb) != *current_index) {
   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);
  1006       pos[ipos] = cvm::atom_pos((pdb->
atom(ipdb))->xcoor(),
  1007                                 (pdb->
atom(ipdb))->ycoor(),
  1008                                 (pdb->
atom(ipdb))->zcoor());
  1010       if (!use_pdb_field && current_index == indices.end())
  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;
  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());
  1038                                      cvm::atom_group& atoms,
  1039                                      std::string 
const &pdb_field_str,
  1040                                      double const pdb_field_value)
  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);
  1046   PDB *pdb = 
new PDB(pdb_filename);
  1047   size_t const pdb_natoms = pdb->
num_atoms();
  1050   auto modify_atoms = atoms.get_atom_modifier();
  1051   for (
size_t ipdb = 0; ipdb < pdb_natoms; ipdb++) {
  1053     double atom_pdb_field_value = 0.0;
  1055     switch (pdb_field_index) {
  1057       atom_pdb_field_value = (pdb->
atom(ipdb))->occupancy();
  1060       atom_pdb_field_value = (pdb->
atom(ipdb))->temperaturefactor();
  1063       atom_pdb_field_value = (pdb->
atom(ipdb))->xcoor();
  1066       atom_pdb_field_value = (pdb->
atom(ipdb))->ycoor();
  1069       atom_pdb_field_value = (pdb->
atom(ipdb))->zcoor();
  1075     if ( (pdb_field_value) &&
  1076          (atom_pdb_field_value != pdb_field_value) ) {
  1078     } 
else if (atom_pdb_field_value == 0.0) {
  1082     if (atoms.is_enabled(colvardeps::f_ag_scalable)) {
  1083       modify_atoms.add_atom_id(ipdb);
  1085       modify_atoms.add_atom(cvm::atom_group::init_atom_from_proxy(
this, ipdb+1));
  1090   return (cvm::get_error() ? COLVARS_ERROR : COLVARS_OK);
  1095                                                std::string 
const description)
  1098     cvm::log(
"Using colvarproxy_namd::output_stream()\n");
  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_;
  1107   if (output_streams_.count(output_name) > 0) {
  1108     return *(output_streams_[output_name]);
  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);
  1119   return *(output_streams_[output_name]);
  1125   if (!io_available()) {
  1129   if (output_streams_.count(output_name) > 0) {
  1130     (
reinterpret_cast<ofstream_namd *
>(output_streams_[output_name]))->flush();
  1134   return cvm::error(
"Error: trying to flush an output file/channel "  1135                     "that wasn't open.\n", COLVARS_BUG_ERROR);
  1141   if (!io_available()) {
  1145   for (std::map<std::string, std::ostream *>::iterator osi = output_streams_.begin();
  1146        osi != output_streams_.end();
  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);
  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);
  1174   if (! io_available()) {
  1178   for (std::map<std::string, std::ostream *>::iterator osi = output_streams_.begin();
  1179        osi != output_streams_.end();
  1183   output_streams_.clear();
  1191   if (std::string(filename).rfind(std::string(
".colvars.state")) != std::string::npos) {
  1203     cvm::log(
"Requesting from NAMD a group of size "+cvm::to_str(atoms_ids.size())+
  1204         " for collective variables calculation.\n");
  1206   colvars->cite_feature(
"Scalable center-of-mass computation (NAMD)");
  1215     bool b_match = 
true;
  1217     if (namd_group.
size() != ((int) atoms_ids.size())) {
  1221       for (ia = 0; ia < namd_group.
size(); ia++) {
  1222         int const aid = atoms_ids[ia];
  1223         if (namd_group[ia] != aid) {
  1232         cvm::log(
"Group was already added.\n");
  1234       atom_groups_refcount[ig] += 1;
  1240   size_t const index = add_atom_group_slot(atom_groups_ids.size());
  1245   namd_group.
resize(atoms_ids.size());
  1247   for (
size_t ia = 0; ia < atoms_ids.size(); ia++) {
  1248     int const aid = atoms_ids[ia];
  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);
  1257     namd_group[ia] = aid;
  1263     cvm::log(
"Group has index "+cvm::to_str(index)+
"\n");
  1275   colvarproxy::clear_atom_group(index);
  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");
  1287   cvm::real total_mass = 0.0;
  1288   cvm::real total_charge = 0.0;
  1289   for (
int i = 0; i < namd_group.
size(); i++) {
  1293   atom_groups_masses[index] = total_mass;
  1294   atom_groups_charges[index] = total_charge;
  1297     cvm::log(
"total mass = "+cvm::to_str(total_mass)+
  1298              ", total charge = "+cvm::to_str(total_charge)+
"\n");
  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;
  1316 #if NAMD_VERSION_NUMBER >= 34471681  1327   for (
size_t i = 0; i < volmaps_ids.size(); i++) {
  1328     if (volmaps_ids[i] == volmap_id) {
  1330       volmaps_refcount[i] += 1;
  1337   if (error_code == COLVARS_OK) {
  1338     index = add_volmap_slot(volmap_id);
  1342   return (error_code == COLVARS_OK) ? index : -1;
  1348   if (volmap_name == NULL) {
  1349     return cvm::error(
"Error: no grid object name provided.", COLVARS_INPUT_ERROR);
  1352   int error_code = COLVARS_OK;
  1357   if (error_code == COLVARS_OK) {
  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);
  1372     for (
size_t i = 0; i < volmaps_ids.size(); i++) {
  1373       if (volmaps_ids[i] == volmap_id) {
  1375         volmaps_refcount[i] += 1;
  1380     index = add_volmap_slot(volmap_id);
  1384   return (error_code == COLVARS_OK) ? index : -1;
  1392     return cvm::error(
"Error: invalid numeric ID ("+cvm::to_str(volmap_id)+
  1393                       ") for map.\n", COLVARS_INPUT_ERROR);
  1395   colvars->cite_feature(
"GridForces volumetric map implementation for NAMD");
  1402   if (volmap_name == NULL) {
  1403     return cvm::error(
"Error: no grid object name provided.", COLVARS_INPUT_ERROR);
  1406   if (volmap_id < 0) {
  1407     return cvm::error(
"Error: invalid map name \""+std::string(volmap_name)+
  1408                       "\".\n", COLVARS_INPUT_ERROR);
  1410   colvars->cite_feature(
"GridForces volumetric map implementation for NAMD");
  1418   colvarproxy::clear_volmap(index);
  1424   int const volmap_id =
  1426   if (volmap_id < 0) {
  1434 template<
class T, 
int flags>
  1436                                          cvm::atom_group* ag,
  1438                                          cvm::real *atom_field)
  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)) {
  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;
  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;
  1472                                              cvm::atom_group* ag,
  1474                                              cvm::real *atom_field)
  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);
  1480     int const new_flags = volmap_flag_gradients;
  1481     GridForceGridLoop<T, new_flags>(g, ag, value, atom_field);
  1487                                      cvm::atom_group* ag,
  1489                                      cvm::real *atom_field)
  1496     getGridForceGridValue<GridforceFullMainGrid>(flags, g, ag,
  1500     getGridForceGridValue<GridforceLiteGrid>(flags, g, ag,
  1508 #if CMK_SMP && USE_CKLOOP // SMP only  1510 colvarproxy::smp_mode_t colvarproxy_namd::get_smp_mode()
 const {
  1514 int colvarproxy_namd::set_smp_mode(smp_mode_t mode) {
  1520 int colvarproxy_namd::smp_loop(
int n_items, std::function<
int (
int)> 
const &worker)
  1522   auto cmkWorker = [&](
int start, 
int end, 
void * ) {
  1523 #if CMK_TRACE_ENABLED  1524     double before = CmiWallTimer();
  1526     for (
int i = start; i <= end; i++) {
  1529 #if CMK_TRACE_ENABLED  1533   const int numChunks = smp_num_threads() > n_items ?
  1536   cvm::increase_depth();
  1537   CkLoop_Parallelize(numChunks, 0, n_items - 1, cmkWorker, 
nullptr, CKLOOP_NONE, 
nullptr);
  1538   cvm::decrease_depth();
  1540   return cvm::get_error();
  1544 void calc_cv_biases_smp(
int first, 
int last, 
void *result, 
int paramNum, 
void *param)
  1547   colvarmodule *cv = proxy->colvars;
  1548 #if CMK_TRACE_ENABLED  1549   double before = CmiWallTimer();
  1551   cvm::increase_depth();
  1552   for (
int i = first; i <= last; i++) {
  1553     colvarbias *b = (*(cv->biases_active()))[i];
  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 = "+
  1562   cvm::decrease_depth();
  1563 #if CMK_TRACE_ENABLED  1569 int colvarproxy_namd::smp_biases_loop()
  1571   colvarmodule *cv = this->colvars;
  1572   const int numChunks = smp_num_threads() > cv->variables_active_smp()->size() ?
  1573                         cv->variables_active_smp()->size() :
  1575   CkLoop_Parallelize(calc_cv_biases_smp, 1, 
this,
  1576                      numChunks, 0, cv->biases_active()->size()-1);
  1577   return cvm::get_error();
  1581 void calc_cv_scripted_forces(
int paramNum, 
void *param)
  1584   colvarmodule *cv = proxy->colvars;
  1585 #if CMK_TRACE_ENABLED  1586   double before = CmiWallTimer();
  1589     cvm::log(
"["+cvm::to_str(proxy->smp_thread_id())+
"/"+cvm::to_str(proxy->smp_num_threads())+
  1590              "]: calc_cv_scripted_forces()\n");
  1592   cv->calc_scripted_forces();
  1593 #if CMK_TRACE_ENABLED  1599 int colvarproxy_namd::smp_biases_script_loop()
  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();
  1609 #endif  // #if CMK_SMP && USE_CKLOOP  1613 #if CMK_HAS_PARTITION  1616   return COLVARS_NOT_IMPLEMENTED;
  1622   return CmiMyPartition();
  1627   return CmiNumPartitions();
  1640   CmiAssert(recvMsg != NULL);
  1641   int retval = recvMsg->
size;
  1642   if (buf_len >= retval) {
  1643     memcpy(msg_data,recvMsg->
data,retval);
  1663     cvm::error(
"computeEnergies must be a divisor of lambda-dynamics period (" + cvm::to_str(freq) + 
").\n");
  1664     return COLVARS_INPUT_ERROR;
  1667     cvm::error(
"alchOn must be enabled for lambda-dynamics.\n");
  1668     return COLVARS_INPUT_ERROR;
  1671     cvm::error(
"alchType must be set to TI for lambda-dynamics.\n");
  1672     return COLVARS_INPUT_ERROR;
  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;
  1702   if (cvm::step_relative() > 0) {
 
Real atomcharge(int anum) const
 
int check_volmaps_available() override
 
std::ostream & output_stream(std::string const &output_name, std::string const description) override
 
const Controller & getController() const
 
void clear_atom_group(int index) override
 
ForceList & modifyAppliedForces()
 
int get_volmap_id_from_name(char const *volmap_name) override
 
AtomIDList & modifyRequestedAtoms()
 
void clear_volmap(int index) override
 
GridforceGrid * get_gridfrc_grid(int gridnum) const
 
BigRealList & modifyGridObjForces()
 
SimParameters * simparams
Pointer to the NAMD simulation input object. 
 
NAMD_HOST_DEVICE Vector c() const
 
Random random
NAMD-style PRNG object. 
 
void NAMD_err(const char *err_msg)
 
BigReal getTIderivative(void) const
 
int init_atom_group(std::vector< int > const &atoms_ids) override
 
NAMD_HOST_DEVICE int c_p() const
 
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
 
void update_accelMD_info()
 
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()
 
AtomIDList::const_iterator getAtomIdBegin()
 
PositionList::const_iterator getGroupPositionEnd()
 
Communication between colvars and NAMD (implementation of colvarproxy) 
 
Controller const  * controller
Pointer to Controller object. 
 
int init_volmap_by_id(int volmap_id) override
 
virtual void submit(void)=0
 
SimParameters * simParameters
 
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)
 
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)
 
ForceList::const_iterator getGroupTotalForceBegin()
 
cvm::step_number previous_NAMD_step
 
SubmitReduction * willSubmit(int setID, int size=-1)
 
int update_target_temperature()
Get the target temperature from the NAMD thermostats supported so far. 
 
static ReductionMgr * Object(void)
 
NAMD_HOST_DEVICE int orthogonal() const
 
NodeReduction * reduction
 
int add(const Elem &elem)
 
PositionList::const_iterator getAtomPositionBegin()
 
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. 
 
NAMD_HOST_DEVICE int b_p() const
 
const ResizeArray< AtomIDList > & requestedGroups()
 
void request_total_force(bool yesno) override
 
char outputFilename[NAMD_FILENAME_BUFFER_SIZE]
 
AtomIDList::const_iterator getForceIdBegin()
 
int init_volmap_by_name(const char *volmap_name) override
 
void setall(const Elem &elem)
 
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. 
 
int index_for_key(const char *key)
 
void log(std::string const &message) override
 
int request_alch_energy_freq(int const freq)
Request energy computation every freq steps. 
 
IntList::const_iterator getGridObjIndexBegin()
 
int get_atom_from_name(const char *segid, int resid, const char *aname) const
 
SubmitReduction * reduction
 
bool accelMDOn
Used to submit restraint energy as MISC. 
 
const AtomID * const_iterator
 
int backup_file(char const *filename) override
 
PDBAtom * atom(int place)
 
void init_atoms_map()
Allocate an atoms map with the same size as the NAMD topology. 
 
int flush_output_streams() override
 
IntList & modifyRequestedGridObjects()
 
#define GLOBAL_MASTER_CKLOOP_CALC_BIASES
 
void replica_comm_barrier() override
 
int send_alch_lambda(void)
Set value of alchemical lambda parameter in back-end. 
 
int check_replicas_enabled() override
 
int close_output_stream(std::string const &output_name) override
 
int replica_index() override
 
NAMD_HOST_DEVICE int a_p() const
 
AtomIDList & modifyForcedAtoms()
 
ResizeArray< AtomIDList > & modifyRequestedGroups()
 
int update_group_properties(int index)
 
void NAMD_die(const char *err_msg)
 
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()
 
Real atommass(int anum) const
 
NAMD_HOST_DEVICE Vector b() const
 
MGridforceParamsList mgridforcelist
 
void init_tcl_pointers() override
 
void update_atom_properties(int index)
 
void NAMD_backup_file(const char *filename, const char *extension)
 
virtual Vector get_scale(void) const =0
 
char restartFilename[NAMD_FILENAME_BUFFER_SIZE]
 
BigRealList::const_iterator getGridObjValueBegin()
 
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()
 
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
 
void requestTotalForce(bool yesno=true)
 
AtomIDList::const_iterator getAtomIdEnd()
 
BigReal getCurrentLambda(const int) const
 
int run_colvar_callback(std::string const &name, std::vector< const colvarvalue *> const &cvcs, colvarvalue &value) override
 
void addReductionEnergy(int reductionTag, BigReal energy)
 
GridforceGridType get_grid_type(void)
 
int num_replicas() override
 
const IntList & requestedGridObjs()
 
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
 
NAMD_HOST_DEVICE Vector a() const
 
int update_atoms_map(AtomIDList::const_iterator begin, AtomIDList::const_iterator end)
 
PositionList::const_iterator getGroupPositionBegin()
 
ForceList::const_iterator getTotalForce()
 
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()
 
NAMD_HOST_DEVICE Vector delta(const Position &pos1, const Position &pos2) const
 
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()