LCOV - code coverage report
Current view: top level - generic - FitToTemplate.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 105 108 97.2 %
Date: 2024-10-18 13:59:31 Functions: 5 7 71.4 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2014-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 "core/ActionWithValue.h"
      26             : #include "tools/Vector.h"
      27             : #include "tools/Matrix.h"
      28             : #include "tools/AtomNumber.h"
      29             : #include "tools/Tools.h"
      30             : #include "tools/RMSD.h"
      31             : #include "core/PlumedMain.h"
      32             : #include "core/ActionSet.h"
      33             : #include "core/GenericMolInfo.h"
      34             : #include "core/PbcAction.h"
      35             : #include "tools/PDB.h"
      36             : #include "tools/Pbc.h"
      37             : 
      38             : #include <vector>
      39             : #include <string>
      40             : #include <memory>
      41             : 
      42             : namespace PLMD {
      43             : namespace generic {
      44             : 
      45             : //+PLUMEDOC GENERIC FIT_TO_TEMPLATE
      46             : /*
      47             : This action is used to align a molecule to a template.
      48             : 
      49             : This can be used to move the coordinates stored in plumed
      50             : so as to be aligned with a provided template in PDB format. Pdb should contain
      51             : also weights for alignment (see the format of PDB files used e.g. for \ref RMSD).
      52             : Make sure your PDB file is correctly formatted as explained \ref pdbreader "in this page".
      53             : Weights for displacement are ignored, since no displacement is computed here.
      54             : Notice that all atoms (not only those in the template) are aligned.
      55             : To see what effect try
      56             : the \ref DUMPATOMS directive to output the atomic positions.
      57             : 
      58             : Also notice that PLUMED propagate forces correctly so that you can add a bias on a CV computed
      59             : after alignment. For many CVs this has no effect, but in some case the alignment can
      60             : change the result. Examples are:
      61             : - \ref POSITION CV since it is affected by a rigid shift of the system.
      62             : - \ref DISTANCE CV with COMPONENTS. Since the alignment could involve a rotation (with TYPE=OPTIMAL) the actual components could be different
      63             :   from the original ones.
      64             : - \ref CELL components for a similar reason.
      65             : - \ref DISTANCE from a \ref FIXEDATOM, provided the fixed atom is introduced _after_ the \ref FIT_TO_TEMPLATE action.
      66             : 
      67             : \attention
      68             : The implementation of TYPE=OPTIMAL is available but should be considered in testing phase. Please report any
      69             : strange behavior.
      70             : 
      71             : \attention
      72             : This directive modifies the stored position at the precise moment
      73             : it is executed. This means that only collective variables
      74             : which are below it in the input script will see the corrected positions.
      75             : As a general rule, put it at the top of the input file. Also, unless you
      76             : know exactly what you are doing, leave the default stride (1), so that
      77             : this action is performed at every MD step.
      78             : 
      79             : When running with periodic boundary conditions, the atoms should be
      80             : in the proper periodic image. This is done automatically since PLUMED 2.5,
      81             : by considering the ordered list of atoms and rebuilding the molecules using a procedure
      82             : that is equivalent to that done in \ref WHOLEMOLECULES . Notice that
      83             : rebuilding is local to this action. This is different from \ref WHOLEMOLECULES
      84             : which actually modifies the coordinates stored in PLUMED.
      85             : 
      86             : In case you want to recover the old behavior you should use the NOPBC flag.
      87             : In that case you need to take care that atoms are in the correct
      88             : periodic image.
      89             : 
      90             : \par Examples
      91             : 
      92             : Align the atomic position to a template then print them.
      93             : The following example is only translating the system so as
      94             : to align the center of mass of a molecule to the one in the reference
      95             : structure `ref.pdb`:
      96             : \plumedfile
      97             : # dump coordinates before fitting, to see the difference:
      98             : DUMPATOMS FILE=dump-before.xyz ATOMS=1-20
      99             : 
     100             : # fit coordinates to ref.pdb template
     101             : # this is a "TYPE=SIMPLE" fit, so that only translations are used.
     102             : FIT_TO_TEMPLATE STRIDE=1 REFERENCE=ref.pdb TYPE=SIMPLE
     103             : 
     104             : # dump coordinates after fitting, to see the difference:
     105             : DUMPATOMS FILE=dump-after.xyz ATOMS=1-20
     106             : \endplumedfile
     107             : 
     108             : The following example instead performs a rototranslational fit.
     109             : \plumedfile
     110             : # dump coordinates before fitting, to see the difference:
     111             : DUMPATOMS FILE=dump-before.xyz ATOMS=1-20
     112             : 
     113             : # fit coordinates to ref.pdb template
     114             : # this is a "TYPE=OPTIMAL" fit, so that rototranslations are used.
     115             : FIT_TO_TEMPLATE STRIDE=1 REFERENCE=ref.pdb TYPE=OPTIMAL
     116             : 
     117             : # dump coordinates after fitting, to see the difference:
     118             : DUMPATOMS FILE=dump-after.xyz ATOMS=1-20
     119             : \endplumedfile
     120             : 
     121             : In both these cases the reference structure should be provided in a reference pdb file such as the one below:
     122             : 
     123             : \auxfile{ref.pdb}
     124             : ATOM      8  HT3 ALA     2      -1.480  -1.560   1.212  1.00  1.00      DIA  H
     125             : ATOM      9  CAY ALA     2      -0.096   2.144  -0.669  1.00  1.00      DIA  C
     126             : ATOM     10  HY1 ALA     2       0.871   2.385  -0.588  1.00  1.00      DIA  H
     127             : ATOM     12  HY3 ALA     2      -0.520   2.679  -1.400  1.00  1.00      DIA  H
     128             : ATOM     14  OY  ALA     2      -1.139   0.931  -0.973  1.00  1.00      DIA  O
     129             : END
     130             : \endauxfile
     131             : 
     132             : In the following example you see two completely equivalent way
     133             : to restrain an atom close to a position that is defined in the reference
     134             : frame of an aligned molecule. You could for instance use this command to calculate the
     135             : position of the center of mass of a ligand after having aligned the atoms to the reference
     136             : frame of the protein that is determined by aligning the atoms in the protein to the coordinates
     137             : provided in the file ref.pdb
     138             : \plumedfile
     139             : # center of the ligand:
     140             : center: CENTER ATOMS=100-110
     141             : 
     142             : FIT_TO_TEMPLATE REFERENCE=ref.pdb TYPE=OPTIMAL
     143             : 
     144             : # place a fixed atom in the protein reference coordinates:
     145             : fix: FIXEDATOM AT=1.0,1.1,1.0
     146             : 
     147             : # take the distance between the fixed atom and the center of the ligand
     148             : d: DISTANCE ATOMS=center,fix
     149             : 
     150             : # apply a restraint
     151             : RESTRAINT ARG=d AT=0.0 KAPPA=100.0
     152             : \endplumedfile
     153             : 
     154             : Notice that you could have obtained an (almost) identical result by adding a fictitious
     155             : atom to `ref.pdb` with the serial number corresponding to the atom labelled `center` (there is no automatic way
     156             : to get it, but in this example it should be the number of atoms of the system plus one),
     157             : and properly setting the weights for alignment and displacement in \ref RMSD.
     158             : There are two differences to be expected:
     159             : (ab) \ref FIT_TO_TEMPLATE might be slower since it has to rototranslate all the available atoms and
     160             : (b) variables employing periodic boundary conditions (such as \ref DISTANCE without `NOPBC`, as in the example above)
     161             :   are allowed after \ref FIT_TO_TEMPLATE, whereas \ref RMSD expects the issues related to the periodic boundary conditions to be already solved.
     162             : The latter means that before the \ref RMSD statement one should use \ref WRAPAROUND or \ref WHOLEMOLECULES to properly place
     163             : the ligand.
     164             : 
     165             : 
     166             : */
     167             : //+ENDPLUMEDOC
     168             : 
     169             : 
     170             : class FitToTemplate:
     171             :   public ActionPilot,
     172             :   public ActionAtomistic,
     173             :   public ActionWithValue
     174             : {
     175             :   std::string type;
     176             :   bool nopbc;
     177             :   std::vector<double> weights;
     178             :   std::vector<std::pair<std::size_t,std::size_t> > p_aligned;
     179             :   Vector center;
     180             :   Vector shift;
     181             :   // optimal alignment related stuff
     182             :   std::unique_ptr<PLMD::RMSD> rmsd;
     183             :   Tensor rotation;
     184             :   Matrix< std::vector<Vector> > drotdpos;
     185             :   // not used anymore (see notes below at doNotRetrieve())
     186             :   // std::vector<Vector> positions;
     187             :   std::vector<Vector> DDistDRef;
     188             :   std::vector<Vector> ddistdpos;
     189             :   std::vector<Vector> centeredpositions;
     190             :   Vector center_positions;
     191             :   // Copy of the box value
     192             :   Value* boxValue;
     193             :   PbcAction* pbc_action;
     194             : public:
     195             :   explicit FitToTemplate(const ActionOptions&ao);
     196             :   static void registerKeywords( Keywords& keys );
     197          55 :   bool actionHasForces() override { return true; }
     198             :   void calculate() override;
     199             :   void apply() override;
     200           0 :   unsigned getNumberOfDerivatives() override {plumed_merror("You should not call this function");};
     201             : };
     202             : 
     203             : PLUMED_REGISTER_ACTION(FitToTemplate,"FIT_TO_TEMPLATE")
     204             : 
     205          11 : void FitToTemplate::registerKeywords( Keywords& keys ) {
     206          11 :   Action::registerKeywords( keys );
     207          11 :   ActionAtomistic::registerKeywords( keys );
     208          22 :   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!");
     209          22 :   keys.add("compulsory","REFERENCE","a file in pdb format containing the reference structure and the atoms involved in the CV.");
     210          22 :   keys.add("compulsory","TYPE","SIMPLE","the manner in which RMSD alignment is performed.  Should be OPTIMAL or SIMPLE.");
     211          22 :   keys.addFlag("NOPBC",false,"ignore the periodic boundary conditions when calculating distances");
     212          22 :   keys.setValueDescription("scalar","the RMSD distance");
     213          11 : }
     214             : 
     215           9 : FitToTemplate::FitToTemplate(const ActionOptions&ao):
     216             :   Action(ao),
     217             :   ActionPilot(ao),
     218             :   ActionAtomistic(ao),
     219             :   ActionWithValue(ao),
     220          18 :   nopbc(false)
     221             : {
     222             :   std::string reference;
     223           9 :   parse("REFERENCE",reference);
     224           9 :   type.assign("SIMPLE");
     225           9 :   parse("TYPE",type);
     226             : 
     227           9 :   parseFlag("NOPBC",nopbc);
     228             : // if(type!="SIMPLE") error("Only TYPE=SIMPLE is implemented in FIT_TO_TEMPLATE");
     229             : 
     230           9 :   checkRead();
     231             : 
     232           9 :   PDB pdb;
     233             : 
     234             :   // read everything in ang and transform to nm if we are not in natural units
     235           9 :   if( !pdb.read(reference,usingNaturalUnits(),0.1/getUnits().getLength()) )
     236           0 :     error("missing input file " + reference );
     237             : 
     238           9 :   requestAtoms(pdb.getAtomNumbers());
     239           9 :   log.printf("  found %zu atoms in input \n",pdb.getAtomNumbers().size());
     240           9 :   log.printf("  with indices : ");
     241          42 :   for(unsigned i=0; i<pdb.getAtomNumbers().size(); ++i) {
     242          33 :     if(i%25==0) log<<"\n";
     243          33 :     log.printf("%d ",pdb.getAtomNumbers()[i].serial());
     244             :   }
     245           9 :   log.printf("\n");
     246             : 
     247           9 :   std::vector<Vector> positions=pdb.getPositions();
     248           9 :   weights=pdb.getOccupancy();
     249           9 :   std::vector<AtomNumber> aligned=pdb.getAtomNumbers(); p_aligned.resize( aligned.size() );
     250          42 :   for(unsigned i=0; i<aligned.size(); ++i) p_aligned[i] = getValueIndices( aligned[i] );
     251             : 
     252             : 
     253             :   // normalize weights
     254          42 :   double n=0.0; for(unsigned i=0; i<weights.size(); ++i) n+=weights[i];
     255           9 :   if(n==0.0) {
     256           0 :     error("PDB file " + reference + " has zero weights. Please check the occupancy column.");
     257             :   }
     258           9 :   n=1.0/n;
     259          42 :   for(unsigned i=0; i<weights.size(); ++i) weights[i]*=n;
     260             : 
     261             :   // normalize weights for rmsd calculation
     262           9 :   std::vector<double> weights_measure=pdb.getBeta();
     263          42 :   n=0.0; for(unsigned i=0; i<weights_measure.size(); ++i) n+=weights_measure[i]; n=1.0/n;
     264          42 :   for(unsigned i=0; i<weights_measure.size(); ++i) weights_measure[i]*=n;
     265             : 
     266             :   // subtract the center
     267          42 :   for(unsigned i=0; i<weights.size(); ++i) center+=positions[i]*weights[i];
     268          42 :   for(unsigned i=0; i<weights.size(); ++i) positions[i]-=center;
     269             : 
     270          13 :   if(type=="OPTIMAL" or type=="OPTIMAL-FAST" ) {
     271           5 :     rmsd=Tools::make_unique<RMSD>();
     272           5 :     rmsd->set(weights,weights_measure,positions,type,false,false);// note: the reference is shifted now with center in the origin
     273          10 :     log<<"  Method chosen for fitting: "<<rmsd->getMethod()<<" \n";
     274             :   }
     275           9 :   if(nopbc) {
     276           1 :     log<<"  Ignoring PBCs when doing alignment, make sure your molecule is whole!<n";
     277             :   }
     278             :   // register the value of rmsd (might be useful sometimes)
     279          18 :   addValue(); setNotPeriodic();
     280             : 
     281             :   // I remove this optimization now in order to use makeWhole()
     282             :   // Notice that for FIT_TO_TEMPLATE TYPE=OPTIMAL a copy was made anyway
     283             :   // (due to the need to store position to propagate forces on rotational matrix later)
     284             :   // For FIT_TO_TEMPLATE TYPE=SIMPLE in principle we could use it and write an ad hoc
     285             :   // version of makeWhole that only computes the center. Too lazy to do it now.
     286             :   // In case we do it later, remember that uncommenting this line means that
     287             :   // getPositions will not work anymore! GB
     288             :   // doNotRetrieve();
     289             : 
     290             :   // this is required so as to allow modifyGlobalForce() to return correct
     291             :   // also for forces that are not owned (and thus not zeored) by all processors.
     292           9 :   pbc_action=plumed.getActionSet().selectWithLabel<PbcAction*>("Box");
     293           9 :   if( !pbc_action ) error("cannot align box has not been set");
     294           9 :   boxValue=pbc_action->copyOutput(0);
     295          18 : }
     296             : 
     297             : 
     298         108 : void FitToTemplate::calculate() {
     299             : 
     300         108 :   if(!nopbc) makeWhole();
     301             : 
     302         108 :   if (type=="SIMPLE") {
     303          48 :     Vector cc;
     304             : 
     305         144 :     for(unsigned i=0; i<p_aligned.size(); ++i) {
     306          96 :       cc+=weights[i]*getPosition(i);
     307             :     }
     308             : 
     309          48 :     shift=center-cc;
     310          48 :     setValue(shift.modulo());
     311          48 :     unsigned nat = getTotAtoms();
     312        6384 :     for(unsigned i=0; i<nat; i++) {
     313        6336 :       std::pair<std::size_t,std::size_t> a = getValueIndices( AtomNumber::index(i));
     314        6336 :       Vector ato=getGlobalPosition(a);
     315        6336 :       setGlobalPosition(a,ato+shift);
     316             :     }
     317             :   }
     318          60 :   else if( type=="OPTIMAL" or type=="OPTIMAL-FAST") {
     319             :     // specific stuff that provides all that is needed
     320          60 :     double r=rmsd->calc_FitElements( getPositions(), rotation,  drotdpos, centeredpositions, center_positions);
     321          60 :     setValue(r); unsigned nat = getTotAtoms();
     322        8004 :     for(unsigned i=0; i<nat; i++) {
     323        7944 :       std::pair<std::size_t,std::size_t> a = getValueIndices( AtomNumber::index(i));
     324        7944 :       Vector ato=getGlobalPosition(a);
     325        7944 :       setGlobalPosition(a,matmul(rotation,ato-center_positions)+center);
     326             :     }
     327             : // rotate box
     328          60 :     Pbc& pbc(pbc_action->getPbc());
     329          60 :     pbc.setBox(matmul(pbc_action->getPbc().getBox(),transpose(rotation)));
     330             :   }
     331         108 : }
     332             : 
     333         108 : void FitToTemplate::apply() {
     334         108 :   auto nat=getTotAtoms();
     335         108 :   if (type=="SIMPLE") {
     336          48 :     Vector totForce;
     337        6384 :     for(unsigned i=0; i<nat; i++) {
     338        6336 :       std::pair<std::size_t,std::size_t> a = getValueIndices( AtomNumber::index(i));
     339        6336 :       totForce+=getForce(a);
     340             :     }
     341          48 :     Tensor vv=Tensor(center,totForce);
     342         624 :     for(unsigned i=0; i<3; ++i) for(unsigned j=0; j<3; ++j) boxValue->addForce( 3*i+j, vv(i,j) );
     343         144 :     for(unsigned i=0; i<p_aligned.size(); ++i) { addForce( p_aligned[i], -totForce*weights[i]); }
     344          60 :   } else if ( type=="OPTIMAL" or type=="OPTIMAL-FAST") {
     345          60 :     Vector totForce;
     346        8004 :     for(unsigned i=0; i<nat; i++) {
     347        7944 :       std::pair<std::size_t,std::size_t> a = getValueIndices( AtomNumber::index(i));
     348        7944 :       Vector f=getForce(a);
     349             : // rotate back forces
     350        7944 :       Vector nf=matmul(transpose(rotation),f);
     351        7944 :       addForce(a, nf-f);
     352             : // accumulate rotated c.o.m. forces - this is already in the non rotated reference frame
     353        7944 :       totForce+=nf;
     354             :     }
     355          60 :     Tensor virial;
     356         780 :     for(unsigned i=0; i<3; ++i) for(unsigned j=0; j<3; ++j) virial[i][j] = boxValue->getForce( 3*i+j );
     357             : // notice that an extra Tensor(center,matmul(rotation,totForce)) is required to
     358             : // compute the derivatives of the rotation with respect to center
     359          60 :     Tensor ww=matmul(transpose(rotation),virial+Tensor(center,matmul(rotation,totForce)));
     360             : // rotate back virial
     361          60 :     virial=matmul(transpose(rotation),matmul(virial,rotation));
     362             : 
     363             : // now we compute the force due to alignment
     364         360 :     for(unsigned i=0; i<p_aligned.size(); i++) {
     365         300 :       Vector g;
     366        1200 :       for(unsigned k=0; k<3; k++) {
     367             : // this could be made faster computing only the diagonal of d
     368         900 :         Tensor d=matmul(ww,RMSD::getMatrixFromDRot(drotdpos,i,k));
     369         900 :         g[k]=(d(0,0)+d(1,1)+d(2,2));
     370             :       }
     371             : // here is the extra contribution
     372         300 :       addForce( p_aligned[i], -g-weights[i]*totForce );
     373             : // here it the contribution to the virial
     374             : // notice that here we can use absolute positions since, for the alignment to be defined,
     375             : // positions should be in one well defined periodic image
     376         300 :       virial+=extProduct(getPosition(i),g);
     377             :     }
     378             : // finally, correction to the virial
     379         120 :     boxValue->clearInputForce(); virial+=extProduct(matmul(transpose(rotation),center),totForce);
     380         780 :     for(unsigned i=0; i<3; ++i) for(unsigned j=0; j<3; ++j) boxValue->addForce( 3*i+j, virial(i,j) );
     381             :   }
     382         108 : }
     383             : 
     384             : }
     385             : }

Generated by: LCOV version 1.16