LCOV - code coverage report
Current view: top level - generic - WholeMolecules.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 89 104 85.6 %
Date: 2025-04-08 21:11:17 Functions: 5 6 83.3 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2011-2023 The plumed team
       3             :    (see the PEOPLE file at the root of the distribution for a list of names)
       4             : 
       5             :    See http://www.plumed.org for more information.
       6             : 
       7             :    This file is part of plumed, version 2.
       8             : 
       9             :    plumed is free software: you can redistribute it and/or modify
      10             :    it under the terms of the GNU Lesser General Public License as published by
      11             :    the Free Software Foundation, either version 3 of the License, or
      12             :    (at your option) any later version.
      13             : 
      14             :    plumed is distributed in the hope that it will be useful,
      15             :    but WITHOUT ANY WARRANTY; without even the implied warranty of
      16             :    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
      17             :    GNU Lesser General Public License for more details.
      18             : 
      19             :    You should have received a copy of the GNU Lesser General Public License
      20             :    along with plumed.  If not, see <http://www.gnu.org/licenses/>.
      21             : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
      22             : #include "core/ActionAtomistic.h"
      23             : #include "core/ActionPilot.h"
      24             : #include "core/ActionRegister.h"
      25             : #include "tools/Vector.h"
      26             : #include "tools/AtomNumber.h"
      27             : #include "tools/Tools.h"
      28             : #include "core/PlumedMain.h"
      29             : #include "core/ActionSet.h"
      30             : #include "core/GenericMolInfo.h"
      31             : #include "tools/OpenMP.h"
      32             : #include "tools/Tree.h"
      33             : 
      34             : #include <vector>
      35             : #include <string>
      36             : 
      37             : namespace PLMD {
      38             : namespace generic {
      39             : 
      40             : //+PLUMEDOC GENERIC WHOLEMOLECULES
      41             : /*
      42             : This action is used to rebuild molecules that can become split by the periodic boundary conditions.
      43             : 
      44             : This command performs an operation that is similar what was done by the ALIGN_ATOMS keyword from plumed1.
      45             : This operation is needed as some MD dynamics code (e.g. GROMACS) can break molecules during the calculation.
      46             : Whenever we are able we try to ensure that molecules are reconstructed automatically.  You thus do not need
      47             : to use this action when you are using actions such as [COM](COM.md), [CENTER](CENTER.md), [GYRATION](GYRATION.md)
      48             : and so on as the reconstruction of molecules is done automatically in these actions.  It is, however, important to understand
      49             : molecule reconstruction as there are many cases (e.g. when you are calculating the end-to-end distance of a polymer)
      50             : where not using the WHOLEMOLCULES command can cause the artifacts discussed in the attached reference.
      51             : 
      52             : If you think that you need to use this command a good idea is to use the [DUMPATOMS](DUMPATOMS.md) directive
      53             : to output the atomic positions.  This will allow you to see the effect that including/not including WHOLEMOLECULES
      54             : has on the calculation.
      55             : 
      56             : > [!ATTENTION]
      57             : > This directive modifies the stored position at the precise moment
      58             : > it is executed. This means that only collective variables
      59             : > which are below it in the input script will see the corrected positions.
      60             : > As a general rule, put it at the top of the input file. Also, unless you
      61             : > know exactly what you are doing, leave the default stride (1), so that
      62             : > this action is performed at every MD step.
      63             : 
      64             : Notice that the behavior of WHOLEMOLECULES is affected by the last [MOLINFO](MOLINFO.md) action
      65             : that is present in the input file before the WHOLEMOLECULES command. Specifically, if the
      66             : [MOLINFO](MOLINFO.md) action does not have a `WHOLE` flag, then the behavior is the following:
      67             : 
      68             : - The first atom of the list is left in place
      69             : - Each atom of the list is shifted by a lattice vectors so that it becomes as close as possible
      70             :   to the previous one, iteratively.
      71             : 
      72             : In this way, if an entity consists of a list of atoms such that consecutive atoms in the
      73             : list are always closer than half a box side the entity will become whole.
      74             : This can be usually achieved selecting consecutive atoms (1-100), but it is also possible
      75             : to skip some atoms, provided consecutive chosen atoms are close enough.
      76             : 
      77             : If, by contrast, the [MOLINFO](MOLINFO.md) action does have a `WHOLE` flag, then a minimum spanning tree
      78             : is built based on the atoms passed to WHOLEMOLECULES using the coordinates in the PDB
      79             : passed to [MOLINFO](MOLINFO.md) as a reference, and this tree is used to reconstruct PBCs.
      80             : This approach is more robust when dealing with complexes of multiple molecules.
      81             : 
      82             : ## Examples
      83             : 
      84             : This command instructs plumed to reconstruct the molecule containing atoms 1-20
      85             : at every step of the calculation and dump them on a file.
      86             : 
      87             : ```plumed
      88             : # to see the effect, one could dump the atoms as they were before molecule reconstruction:
      89             : # DUMPATOMS FILE=dump-broken.xyz ATOMS=1-20
      90             : WHOLEMOLECULES ENTITY0=1-20
      91             : DUMPATOMS FILE=dump.xyz ATOMS=1-20
      92             : ```
      93             : 
      94             : This command instructs plumed to reconstruct two molecules containing atoms 1-20 and 30-40
      95             : 
      96             : ```plumed
      97             : WHOLEMOLECULES ENTITY0=1-20 ENTITY1=30-40
      98             : DUMPATOMS FILE=dump.xyz ATOMS=1-20,30-40
      99             : ```
     100             : 
     101             : This command instructs plumed to reconstruct the chain of backbone atoms in a
     102             : protein
     103             : 
     104             : ```plumed
     105             : #SETTINGS MOLFILE=regtest/basic/rt32/helix.pdb
     106             : MOLINFO STRUCTURE=regtest/basic/rt32/helix.pdb
     107             : WHOLEMOLECULES RESIDUES=all MOLTYPE=protein
     108             : ```
     109             : 
     110             : */
     111             : //+ENDPLUMEDOC
     112             : 
     113             : 
     114             : class WholeMolecules:
     115             :   public ActionPilot,
     116             :   public ActionAtomistic {
     117             :   std::vector<std::vector<std::pair<std::size_t,std::size_t> > > p_groups;
     118             :   std::vector<std::vector<std::pair<std::size_t,std::size_t> > > p_roots;
     119             :   std::vector<Vector> refs;
     120             :   bool doemst, addref;
     121             : public:
     122             :   explicit WholeMolecules(const ActionOptions&ao);
     123             :   static void registerKeywords( Keywords& keys );
     124         841 :   bool actionHasForces() override {
     125         841 :     return false;
     126             :   }
     127             :   void calculate() override;
     128       11946 :   void apply() override {}
     129             : };
     130             : 
     131             : PLUMED_REGISTER_ACTION(WholeMolecules,"WHOLEMOLECULES")
     132             : 
     133          93 : void WholeMolecules::registerKeywords( Keywords& keys ) {
     134          93 :   Action::registerKeywords( keys );
     135          93 :   ActionPilot::registerKeywords( keys );
     136          93 :   ActionAtomistic::registerKeywords( keys );
     137          93 :   keys.add("compulsory","STRIDE","1","the frequency with which molecules are reassembled.  Unless you are completely certain about what you are doing leave this set equal to 1!");
     138          93 :   keys.add("numbered","ENTITY","the atoms that make up a molecule that you wish to align. To specify multiple molecules use a list of ENTITY keywords: ENTITY0, ENTITY1,...");
     139         186 :   keys.reset_style("ENTITY","atoms");
     140          93 :   keys.add("residues","RESIDUES","this command specifies that the backbone atoms in a set of residues all must be aligned. It must be used in tandem with the \\ref MOLINFO "
     141             :            "action and the MOLTYPE keyword. If you wish to use all the residues from all the chains in your system you can do so by "
     142             :            "specifying all. Alternatively, if you wish to use a subset of the residues you can specify the particular residues "
     143             :            "you are interested in as a list of numbers");
     144          93 :   keys.add("optional","MOLTYPE","the type of molecule that is under study.  This is used to define the backbone atoms");
     145          93 :   keys.addFlag("EMST", false, "only for backward compatibility, as of PLUMED 2.11 this is the default when using MOLINFO with WHOLE");
     146          93 :   keys.addFlag("ADDREFERENCE", false, "Define the reference position of the first atom of each entity using a PDB file");
     147          93 :   keys.addDOI("10.1007/978-1-4939-9608-7_21");
     148          93 : }
     149             : 
     150          71 : WholeMolecules::WholeMolecules(const ActionOptions&ao):
     151             :   Action(ao),
     152             :   ActionPilot(ao),
     153             :   ActionAtomistic(ao),
     154          71 :   doemst(false), addref(false) {
     155             :   std::vector<std::vector<AtomNumber> > groups;
     156             :   std::vector<std::vector<AtomNumber> > roots;
     157             :   // parse optional flags
     158             :   bool doemst_tmp;
     159          71 :   parseFlag("EMST", doemst_tmp);
     160          71 :   if(doemst_tmp) {
     161           1 :     log << "EMST option is not needed any more as of PLUMED 2.11\n";
     162             :   }
     163          71 :   parseFlag("ADDREFERENCE", addref);
     164             : 
     165          71 :   auto* moldat=plumed.getActionSet().selectLatest<GenericMolInfo*>(this);
     166             : 
     167             :   // create groups from ENTITY
     168         422 :   for(int i=0;; i++) {
     169             :     std::vector<AtomNumber> group;
     170         986 :     parseAtomList("ENTITY",i,group);
     171         493 :     if( group.empty() ) {
     172             :       break;
     173             :     }
     174         422 :     groups.push_back(group);
     175         422 :   }
     176             : 
     177             :   // Read residues to align from MOLINFO
     178             :   std::vector<std::string> resstrings;
     179         142 :   parseVector("RESIDUES",resstrings);
     180          71 :   if( resstrings.size()>0 ) {
     181           0 :     if( resstrings.size()==1 ) {
     182           0 :       if( resstrings[0]=="all" ) {
     183             :         resstrings[0]="all-ter";  // Include terminal groups in alignment
     184             :       }
     185             :     }
     186             :     std::string moltype;
     187           0 :     parse("MOLTYPE",moltype);
     188           0 :     if(moltype.length()==0) {
     189           0 :       error("Found RESIDUES keyword without specification of the molecule - use MOLTYPE");
     190             :     }
     191           0 :     if( !moldat ) {
     192           0 :       error("MOLINFO is required to use RESIDUES");
     193             :     }
     194             :     std::vector< std::vector<AtomNumber> > backatoms;
     195           0 :     moldat->getBackbone( resstrings, moltype, backatoms );
     196           0 :     for(unsigned i=0; i<backatoms.size(); ++i) {
     197           0 :       groups.push_back( backatoms[i] );
     198             :     }
     199           0 :   }
     200             : 
     201             :   // check number of groups
     202          71 :   if(groups.size()==0) {
     203           0 :     error("no atoms found for WHOLEMOLECULES!");
     204             :   }
     205             : 
     206             :   // if using PDBs reorder atoms in groups based on proximity in PDB file
     207          71 :   if(moldat && moldat->isWhole()) {
     208           2 :     doemst=true;
     209             :   }
     210             : 
     211          71 :   if(doemst_tmp && ! doemst) {
     212           0 :     error("cannot enable EMST if MOLINFO is not WHOLE");
     213             :   }
     214             : 
     215          71 :   if(doemst) {
     216           2 :     if( !moldat ) {
     217           0 :       error("MOLINFO is required to use EMST");
     218             :     }
     219             :     // initialize tree
     220           2 :     Tree tree = Tree(moldat);
     221             :     // cycle on groups and reorder atoms
     222           4 :     for(unsigned i=0; i<groups.size(); ++i) {
     223           2 :       groups[i] = tree.getTree(groups[i]);
     224             :       // store root atoms
     225           2 :       roots.push_back(tree.getRoot());
     226             :     }
     227           2 :   } else {
     228             :     // fill root vector with previous atom in groups
     229         489 :     for(unsigned i=0; i<groups.size(); ++i) {
     230             :       std::vector<AtomNumber> root;
     231       11390 :       for(unsigned j=0; j<groups[i].size()-1; ++j) {
     232       10970 :         root.push_back(groups[i][j]);
     233             :       }
     234             :       // store root atoms
     235         420 :       roots.push_back(root);
     236             :     }
     237             :   }
     238             : 
     239             :   // adding reference if needed
     240          71 :   if(addref) {
     241           2 :     if( !moldat ) {
     242           0 :       error("MOLINFO is required to use ADDREFERENCE");
     243             :     }
     244           4 :     for(unsigned i=0; i<groups.size(); ++i) {
     245             :       // add reference position of first atom in entity
     246           4 :       refs.push_back(moldat->getPosition(groups[i][0]));
     247             :     }
     248             :   }
     249             : 
     250             :   // print out info
     251         493 :   for(unsigned i=0; i<groups.size(); ++i) {
     252         422 :     log.printf("  atoms in entity %d : ",i);
     253       12057 :     for(unsigned j=0; j<groups[i].size(); ++j) {
     254       11635 :       log.printf("%d ",groups[i][j].serial() );
     255             :     }
     256         422 :     log.printf("\n");
     257         422 :     if(addref) {
     258           2 :       log.printf("     with reference position : %lf %lf %lf\n",refs[i][0],refs[i][1],refs[i][2]);
     259             :     }
     260             :   }
     261             : 
     262             :   // collect all atoms
     263             :   std::vector<AtomNumber> merge;
     264         493 :   for(unsigned i=0; i<groups.size(); ++i) {
     265         422 :     merge.insert(merge.end(),groups[i].begin(),groups[i].end());
     266             :   }
     267             : 
     268             :   // Convert groups to p_groups
     269          71 :   p_groups.resize( groups.size() );
     270         493 :   for(unsigned i=0; i<groups.size(); ++i) {
     271         422 :     p_groups[i].resize( groups[i].size() );
     272       12057 :     for(unsigned j=0; j<groups[i].size(); ++j) {
     273       11635 :       p_groups[i][j] = getValueIndices( groups[i][j] );
     274             :     }
     275             :   }
     276             :   // Convert roots to p_roots
     277          71 :   p_roots.resize( roots.size() );
     278         493 :   for(unsigned i=0; i<roots.size(); ++i) {
     279         422 :     p_roots[i].resize( roots[i].size() );
     280       11635 :     for(unsigned j=0; j<roots[i].size(); ++j) {
     281       11213 :       p_roots[i][j] = getValueIndices( roots[i][j] );
     282             :     }
     283             :   }
     284             : 
     285             : 
     286          71 :   checkRead();
     287          71 :   Tools::removeDuplicates(merge);
     288          71 :   requestAtoms(merge);
     289             :   doNotRetrieve();
     290             :   doNotForce();
     291          71 : }
     292             : 
     293       11946 : void WholeMolecules::calculate() {
     294       24596 :   for(unsigned i=0; i<p_groups.size(); ++i) {
     295       12650 :     Vector first = getGlobalPosition(p_groups[i][0]);
     296       12650 :     if(addref) {
     297          12 :       first = refs[i]+pbcDistance(refs[i],first);
     298          12 :       setGlobalPosition( p_groups[i][0], first );
     299             :     }
     300       12650 :     if(!doemst) {
     301      662472 :       for(unsigned j=1; j<p_groups[i].size(); ++j) {
     302      649826 :         Vector second=getGlobalPosition(p_groups[i][j]);
     303      649826 :         first = first+pbcDistance(first,second);
     304      649826 :         setGlobalPosition(p_groups[i][j], first );
     305             :       }
     306             :     } else {
     307         490 :       for(unsigned j=1; j<p_groups[i].size(); ++j) {
     308         486 :         Vector first=getGlobalPosition(p_roots[i][j-1]);
     309         486 :         Vector second=getGlobalPosition(p_groups[i][j]);
     310         486 :         second=first+pbcDistance(first,second);
     311         486 :         setGlobalPosition(p_groups[i][j], second );
     312             :       }
     313             :     }
     314             :   }
     315       11946 : }
     316             : 
     317             : 
     318             : }
     319             : }

Generated by: LCOV version 1.16